version 3.8
multidomain/facet/cellcentered/tpfa/fourierslaw.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_DISCRETIZATION_CC_TPFA_FACET_COUPLING_FOURIERS_LAW_HH
13#define DUMUX_DISCRETIZATION_CC_TPFA_FACET_COUPLING_FOURIERS_LAW_HH
14
15#include <array>
16#include <cmath>
17
18#include <dune/common/float_cmp.hh>
19#include <dune/common/fvector.hh>
20
21#include <dumux/common/math.hh>
24
28
30
31namespace Dumux {
32
34template<class TypeTag,
35 bool isNetwork>
37
44template<class TypeTag>
48
53template<class TypeTag>
54class CCTpfaFacetCouplingFouriersLawImpl<TypeTag, /*isNetwork*/false>
55: public FouriersLawImplementation<TypeTag, DiscretizationMethods::CCTpfa>
56{
59
62 using FVElementGeometry = typename GridGeometry::LocalView;
63 using SubControlVolume = typename GridGeometry::SubControlVolume;
64 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
65 using Extrusion = Extrusion_t<GridGeometry>;
66
67 using GridView = typename GridGeometry::GridView;
68 using Element = typename GridView::template Codim<0>::Entity;
69
73 class FacetCouplingFouriersLawCache
74 {
75 public:
77 using Filler = typename ParentType::Cache::Filler;
78
82 static constexpr int insideTijIdx = 0;
83 static constexpr int outsideTijIdx = 1;
84 static constexpr int facetTijIdx = 2;
85
87 using HeatConductionTransmissibilityContainer = std::array<Scalar, 3>;
88
90 template< class Problem, class ElementVolumeVariables >
91 void updateHeatConduction(const Problem& problem,
92 const Element& element,
93 const FVElementGeometry& fvGeometry,
94 const ElementVolumeVariables& elemVolVars,
95 const SubControlVolumeFace &scvf)
96 {
97 tij_ = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
98 }
99
104 Scalar heatConductionTij() const
105 { return tij_[insideTijIdx]; }
106
108 Scalar heatConductionTijInside() const
109 { return tij_[insideTijIdx]; }
110
112 Scalar heatConductionTijOutside() const
113 { return tij_[outsideTijIdx]; }
114
116 Scalar heatConductionTijFacet() const
117 { return tij_[facetTijIdx]; }
118
119 private:
120 HeatConductionTransmissibilityContainer tij_;
121 };
122
123public:
125 using Cache = FacetCouplingFouriersLawCache;
126
129 static constexpr DiscretizationMethod discMethod{};
130
132 template< class Problem, class ElementVolumeVariables, class ElementFluxVarsCache >
133 static Scalar flux(const Problem& problem,
134 const Element& element,
135 const FVElementGeometry& fvGeometry,
136 const ElementVolumeVariables& elemVolVars,
137 const SubControlVolumeFace& scvf,
138 const ElementFluxVarsCache& elemFluxVarsCache)
139 {
140 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
141 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
142
143 // get inside/outside volume variables
144 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
145 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
146 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
147 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
148
149 // the inside and outside temperatures
150 const Scalar tInside = insideVolVars.temperature();
151 const Scalar tFacet = facetVolVars.temperature();
152 const Scalar tOutside = outsideVolVars.temperature();
153
154 Scalar flux = fluxVarsCache.heatConductionTijInside()*tInside
155 + fluxVarsCache.heatConductionTijFacet()*tFacet;
156
157 if (!scvf.boundary())
158 flux += fluxVarsCache.heatConductionTijOutside()*tOutside;
159 return flux;
160 }
161
162 // The flux variables cache has to be bound to an element prior to flux calculations
163 // During the binding, the transmissibility will be computed and stored using the method below.
164 template< class Problem, class ElementVolumeVariables >
165 static typename Cache::HeatConductionTransmissibilityContainer
166 calculateTransmissibility(const Problem& problem,
167 const Element& element,
168 const FVElementGeometry& fvGeometry,
169 const ElementVolumeVariables& elemVolVars,
170 const SubControlVolumeFace& scvf)
171 {
172 typename Cache::HeatConductionTransmissibilityContainer tij;
173 if (!problem.couplingManager().isCoupled(element, scvf))
174 {
176 tij[Cache::insideTijIdx] = ParentType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
177 return tij;
178 }
179
181 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(), "FacetCoupling.Xi", 1.0);
182
183 const auto insideScvIdx = scvf.insideScvIdx();
184 const auto& insideScv = fvGeometry.scv(insideScvIdx);
185 const auto& insideVolVars = elemVolVars[insideScvIdx];
186 const auto wIn = Extrusion::area(fvGeometry, scvf)
187 *computeTpfaTransmissibility(fvGeometry, scvf, insideScv,
188 insideVolVars.effectiveThermalConductivity(),
189 insideVolVars.extrusionFactor());
190
191 // proceed depending on the interior BC types used
192 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
193
194 // neumann-coupling
195 if (iBcTypes.hasOnlyNeumann())
196 {
197 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
198 const auto wFacet = 2.0*Extrusion::area(fvGeometry, scvf)*insideVolVars.extrusionFactor()
199 /facetVolVars.extrusionFactor()
200 *vtmv(scvf.unitOuterNormal(),
201 facetVolVars.effectiveThermalConductivity(),
202 scvf.unitOuterNormal());
203
204 // The fluxes across this face and the outside face can be expressed in matrix form:
205 // \f$\mathbf{C} \bar{\mathbf{u}} + \mathbf{D} \mathbf{u} + \mathbf{E} \mathbf{u}_\gamma\f$,
206 // where \f$\gamma$\f denotes the domain living on the facets and \f$\bar{\mathbf{u}}$\f are
207 // intermediate face unknowns in the matrix domain. Equivalently, flux continuity reads:
208 // \f$\mathbf{A} \bar{\mathbf{u}} = \mathbf{B} \mathbf{u} + \mathbf{M} \mathbf{u}_\gamma\f$.
209 // Combining the two, we can eliminate the intermediate unknowns and compute the transmissibilities
210 // that allow the description of the fluxes as functions of the cell and Dirichlet temperatures only.
211 if (!scvf.boundary())
212 {
213 const auto outsideScvIdx = scvf.outsideScvIdx();
214 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
215 const auto wOut = -1.0*Extrusion::area(fvGeometry, scvf)
216 *computeTpfaTransmissibility(fvGeometry, scvf, fvGeometry.scv(outsideScvIdx),
217 outsideVolVars.effectiveThermalConductivity(),
218 outsideVolVars.extrusionFactor());
219
220 if ( !Dune::FloatCmp::eq(xi, 1.0, 1e-6) )
221 {
222 // optimized implementation: factorization obtained using sympy
223 // see CCTpfaFacetCouplingDarcysLaw for more details
224 const Scalar factor = wIn * wFacet / ( wIn * wOut * ( 2.0 * xi - 1.0 ) + wFacet * ( xi * ( wIn + wOut ) + wFacet ) );
225 tij[Cache::insideTijIdx] = factor * ( wOut * xi + wFacet );
226 tij[Cache::outsideTijIdx] = factor * ( wOut * ( 1.0 - xi ) );
227 tij[Cache::facetTijIdx] = factor * ( - wOut - wFacet );
228 }
229 else
230 {
231 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
232 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
233 tij[Cache::outsideTijIdx] = 0.0;
234 }
235 }
236 else
237 {
238 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
239 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
240 tij[Cache::outsideTijIdx] = 0.0;
241 }
242 }
243 else if (iBcTypes.hasOnlyDirichlet())
244 {
245 tij[Cache::insideTijIdx] = wIn;
246 tij[Cache::outsideTijIdx] = 0.0;
247 tij[Cache::facetTijIdx] = -wIn;
248 }
249 else
250 DUNE_THROW(Dune::NotImplemented, "Interior boundary types other than pure Dirichlet or Neumann");
251
252 return tij;
253 }
254};
255
260template<class TypeTag>
261class CCTpfaFacetCouplingFouriersLawImpl<TypeTag, /*isNetwork*/true>
262: public FouriersLawImplementation<TypeTag, DiscretizationMethods::CCTpfa>
263{
266
269 using FVElementGeometry = typename GridGeometry::LocalView;
270 using SubControlVolume = typename GridGeometry::SubControlVolume;
271 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
272 using Extrusion = Extrusion_t<GridGeometry>;
273
274 using GridView = typename GridGeometry::GridView;
275 using Element = typename GridView::template Codim<0>::Entity;
276
280 class FacetCouplingFouriersLawCache
281 {
282 public:
284 using Filler = typename ParentType::Cache::Filler;
285
289 static constexpr int insideTijIdx = 0;
290 static constexpr int facetTijIdx = 1;
291
293 using HeatConductionTransmissibilityContainer = std::array<Scalar, 2>;
294
296 template< class Problem, class ElementVolumeVariables >
297 void updateHeatConduction(const Problem& problem,
298 const Element& element,
299 const FVElementGeometry& fvGeometry,
300 const ElementVolumeVariables& elemVolVars,
301 const SubControlVolumeFace &scvf)
302 {
303 tij_ = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
304 }
305
310 Scalar heatConductionTij() const
311 { return tij_[insideTijIdx]; }
312
314 Scalar heatConductionTijInside() const
315 { return tij_[insideTijIdx]; }
316
318 Scalar heatConductionTijFacet() const
319 { return tij_[facetTijIdx]; }
320
321 private:
322 HeatConductionTransmissibilityContainer tij_;
323 };
324
325public:
327 using Cache = FacetCouplingFouriersLawCache;
328
331 static constexpr DiscretizationMethod discMethod{};
332
334 template< class Problem, class ElementVolumeVariables, class ElementFluxVarsCache >
335 static Scalar flux(const Problem& problem,
336 const Element& element,
337 const FVElementGeometry& fvGeometry,
338 const ElementVolumeVariables& elemVolVars,
339 const SubControlVolumeFace& scvf,
340 const ElementFluxVarsCache& elemFluxVarsCache)
341 {
342 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
343 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
344
345 // get inside/facet volume variables
346 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
347 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
348 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
349
350 return fluxVarsCache.heatConductionTijInside()*insideVolVars.temperature()
351 + fluxVarsCache.heatConductionTijFacet()*facetVolVars.temperature();
352 }
353
354 // The flux variables cache has to be bound to an element prior to flux calculations
355 // During the binding, the transmissibility will be computed and stored using the method below.
356 template< class Problem, class ElementVolumeVariables >
357 static typename Cache::HeatConductionTransmissibilityContainer
358 calculateTransmissibility(const Problem& problem,
359 const Element& element,
360 const FVElementGeometry& fvGeometry,
361 const ElementVolumeVariables& elemVolVars,
362 const SubControlVolumeFace& scvf)
363 {
364 typename Cache::HeatConductionTransmissibilityContainer tij;
365 if (!problem.couplingManager().isCoupled(element, scvf))
366 {
368 tij[Cache::insideTijIdx] = ParentType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
369 return tij;
370 }
371
373 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(), "FacetCoupling.Xi", 1.0);
374
375 // On surface grids only xi = 1.0 can be used, as the coupling condition
376 // for xi != 1.0 does not generalize for surface grids where the normal
377 // vectors of the inside/outside elements have different orientations.
378 if (Dune::FloatCmp::ne(xi, 1.0, 1e-6))
379 DUNE_THROW(Dune::InvalidStateException, "Xi != 1.0 cannot be used on surface grids");
380
381 const auto insideScvIdx = scvf.insideScvIdx();
382 const auto& insideScv = fvGeometry.scv(insideScvIdx);
383 const auto& insideVolVars = elemVolVars[insideScvIdx];
384 const auto wIn = Extrusion::area(fvGeometry, scvf)
385 *computeTpfaTransmissibility(fvGeometry, scvf, insideScv,
386 insideVolVars.effectiveThermalConductivity(),
387 insideVolVars.extrusionFactor());
388
389 // proceed depending on the interior BC types used
390 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
391
392 // neumann-coupling
393 if (iBcTypes.hasOnlyNeumann())
394 {
395 // Here we use the square root of the facet extrusion factor
396 // as an approximate average distance from scvf ip to facet center
397 using std::sqrt;
398 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
399 const auto wFacet = 2.0*Extrusion::area(fvGeometry, scvf)*insideVolVars.extrusionFactor()
400 /sqrt(facetVolVars.extrusionFactor())
401 *vtmv(scvf.unitOuterNormal(),
402 facetVolVars.effectiveThermalConductivity(),
403 scvf.unitOuterNormal());
404
405 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
406 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
407 }
408 else if (iBcTypes.hasOnlyDirichlet())
409 {
410 tij[Cache::insideTijIdx] = wIn;
411 tij[Cache::facetTijIdx] = -wIn;
412 }
413 else
414 DUNE_THROW(Dune::NotImplemented, "Interior boundary types other than pure Dirichlet or Neumann");
415
416 return tij;
417 }
418};
419
420} // end namespace Dumux
421
422#endif
Specialization of CCTpfaFacetCouplingFouriersLawImpl for dim=dimWorld.
Definition: multidomain/facet/cellcentered/tpfa/fourierslaw.hh:56
static Cache::HeatConductionTransmissibilityContainer calculateTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Definition: multidomain/facet/cellcentered/tpfa/fourierslaw.hh:166
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVarsCache &elemFluxVarsCache)
Compute the conductive heat flux.
Definition: multidomain/facet/cellcentered/tpfa/fourierslaw.hh:133
Specialization of CCTpfaFacetCouplingFouriersLawImpl for dim<dimWorld.
Definition: multidomain/facet/cellcentered/tpfa/fourierslaw.hh:263
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVarsCache &elemFluxVarsCache)
Compute the conductive heat flux.
Definition: multidomain/facet/cellcentered/tpfa/fourierslaw.hh:335
static Cache::HeatConductionTransmissibilityContainer calculateTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Definition: multidomain/facet/cellcentered/tpfa/fourierslaw.hh:358
Forward declaration of the implementation.
Definition: multidomain/facet/cellcentered/tpfa/fourierslaw.hh:36
Fourier's law for cell-centered finite volume schemes with two-point flux approximation.
Definition: flux/cctpfa/fourierslaw.hh:34
TpfaFouriersLawCache Cache
export the type for the corresponding cache
Definition: flux/cctpfa/fourierslaw.hh:100
forward declaration of the method-specific implementation
Definition: flux/box/fourierslaw.hh:26
Defines all properties used in Dumux.
Helper classes to compute the integration elements.
Fourier's law for cell-centered finite volume schemes with two-point flux approximation.
Tensor::field_type computeTpfaTransmissibility(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const typename FVElementGeometry::SubControlVolume &scv, const Tensor &T, typename FVElementGeometry::SubControlVolume::Traits::Scalar extrusionFactor)
Free function to evaluate the Tpfa transmissibility associated with the flux (in the form of flux = T...
Definition: tpfa/computetransmissibility.hh:36
Dune::DenseMatrix< MAT >::value_type vtmv(const Dune::DenseVector< V1 > &v1, const Dune::DenseMatrix< MAT > &M, const Dune::DenseVector< V2 > &v2)
Evaluates the scalar product of a vector v2, projected by a matrix M, with a vector v1.
Definition: math.hh:851
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Define some often used mathematical functions.
The available discretization methods in Dumux.
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Definition: method.hh:25
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...