3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_DISCRETIZATION_CC_TPFA_FACET_COUPLING_FOURIERS_LAW_HH
25#define DUMUX_DISCRETIZATION_CC_TPFA_FACET_COUPLING_FOURIERS_LAW_HH
26
27#include <array>
28#include <cmath>
29
30#include <dune/common/float_cmp.hh>
31#include <dune/common/fvector.hh>
32
33#include <dumux/common/math.hh>
36
40
42
43namespace Dumux {
44
46template<class TypeTag,
47 bool isNetwork>
49
56template<class TypeTag>
60
65template<class TypeTag>
66class CCTpfaFacetCouplingFouriersLawImpl<TypeTag, /*isNetwork*/false>
67: public FouriersLawImplementation<TypeTag, DiscretizationMethod::cctpfa>
68{
71
74 using FVElementGeometry = typename GridGeometry::LocalView;
75 using SubControlVolume = typename GridGeometry::SubControlVolume;
76 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
77 using Extrusion = Extrusion_t<GridGeometry>;
78
79 using GridView = typename GridGeometry::GridView;
80 using Element = typename GridView::template Codim<0>::Entity;
81
85 class FacetCouplingFouriersLawCache
86 {
87 public:
89 using Filler = typename ParentType::Cache::Filler;
90
94 static constexpr int insideTijIdx = 0;
95 static constexpr int outsideTijIdx = 1;
96 static constexpr int facetTijIdx = 2;
97
99 using HeatConductionTransmissibilityContainer = std::array<Scalar, 3>;
100
102 template< class Problem, class ElementVolumeVariables >
103 void updateHeatConduction(const Problem& problem,
104 const Element& element,
105 const FVElementGeometry& fvGeometry,
106 const ElementVolumeVariables& elemVolVars,
107 const SubControlVolumeFace &scvf)
108 {
109 tij_ = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
110 }
111
116 Scalar heatConductionTij() const
117 { return tij_[insideTijIdx]; }
118
120 Scalar heatConductionTijInside() const
121 { return tij_[insideTijIdx]; }
122
124 Scalar heatConductionTijOutside() const
125 { return tij_[outsideTijIdx]; }
126
128 Scalar heatConductionTijFacet() const
129 { return tij_[facetTijIdx]; }
130
131 private:
132 HeatConductionTransmissibilityContainer tij_;
133 };
134
135public:
137 using Cache = FacetCouplingFouriersLawCache;
138
141
143 template< class Problem, class ElementVolumeVariables, class ElementFluxVarsCache >
144 static Scalar flux(const Problem& problem,
145 const Element& element,
146 const FVElementGeometry& fvGeometry,
147 const ElementVolumeVariables& elemVolVars,
148 const SubControlVolumeFace& scvf,
149 const ElementFluxVarsCache& elemFluxVarsCache)
150 {
151 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
152 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
153
154 // get inside/outside volume variables
155 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
156 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
157 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
158 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
159
160 // the inside and outside temperatures
161 const Scalar tInside = insideVolVars.temperature();
162 const Scalar tFacet = facetVolVars.temperature();
163 const Scalar tOutside = outsideVolVars.temperature();
164
165 Scalar flux = fluxVarsCache.heatConductionTijInside()*tInside
166 + fluxVarsCache.heatConductionTijFacet()*tFacet;
167
168 if (!scvf.boundary())
169 flux += fluxVarsCache.heatConductionTijOutside()*tOutside;
170 return flux;
171 }
172
173 // The flux variables cache has to be bound to an element prior to flux calculations
174 // During the binding, the transmissibility will be computed and stored using the method below.
175 template< class Problem, class ElementVolumeVariables >
176 static typename Cache::HeatConductionTransmissibilityContainer
177 calculateTransmissibility(const Problem& problem,
178 const Element& element,
179 const FVElementGeometry& fvGeometry,
180 const ElementVolumeVariables& elemVolVars,
181 const SubControlVolumeFace& scvf)
182 {
183 typename Cache::HeatConductionTransmissibilityContainer tij;
184 if (!problem.couplingManager().isCoupled(element, scvf))
185 {
187 tij[Cache::insideTijIdx] = ParentType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
188 return tij;
189 }
190
192 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(), "FacetCoupling.Xi", 1.0);
193
194 const auto insideScvIdx = scvf.insideScvIdx();
195 const auto& insideScv = fvGeometry.scv(insideScvIdx);
196 const auto& insideVolVars = elemVolVars[insideScvIdx];
197 const auto wIn = Extrusion::area(scvf)
198 *computeTpfaTransmissibility(scvf, insideScv,
199 insideVolVars.effectiveThermalConductivity(),
200 insideVolVars.extrusionFactor());
201
202 // proceed depending on the interior BC types used
203 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
204
205 // neumann-coupling
206 if (iBcTypes.hasOnlyNeumann())
207 {
208 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
209 const auto wFacet = 2.0*Extrusion::area(scvf)*insideVolVars.extrusionFactor()
210 /facetVolVars.extrusionFactor()
211 *vtmv(scvf.unitOuterNormal(),
212 facetVolVars.effectiveThermalConductivity(),
213 scvf.unitOuterNormal());
214
215 // The fluxes across this face and the outside face can be expressed in matrix form:
216 // \f$\mathbf{C} \bar{\mathbf{u}} + \mathbf{D} \mathbf{u} + \mathbf{E} \mathbf{u}_\gamma\f$,
217 // where \f$\gamma$\f denotes the domain living on the facets and \f$\bar{\mathbf{u}}$\f are
218 // intermediate face unknowns in the matrix domain. Equivalently, flux continuity reads:
219 // \f$\mathbf{A} \bar{\mathbf{u}} = \mathbf{B} \mathbf{u} + \mathbf{M} \mathbf{u}_\gamma\f$.
220 // Combining the two, we can eliminate the intermediate unknowns and compute the transmissibilities
221 // that allow the description of the fluxes as functions of the cell and Dirichlet temperatures only.
222 if (!scvf.boundary())
223 {
224 const auto outsideScvIdx = scvf.outsideScvIdx();
225 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
226 const auto wOut = -1.0*Extrusion::area(scvf)
227 *computeTpfaTransmissibility(scvf, fvGeometry.scv(outsideScvIdx),
228 outsideVolVars.effectiveThermalConductivity(),
229 outsideVolVars.extrusionFactor());
230
231 if ( !Dune::FloatCmp::eq(xi, 1.0, 1e-6) )
232 {
233 // optimized implementation: factorization obtained using sympy
234 // see CCTpfaFacetCouplingDarcysLaw for more details
235 const Scalar factor = wIn * wFacet / ( wIn * wOut * ( 2.0 * xi - 1.0 ) + wFacet * ( xi * ( wIn + wOut ) + wFacet ) );
236 tij[Cache::insideTijIdx] = factor * ( wOut * xi + wFacet );
237 tij[Cache::outsideTijIdx] = factor * ( wOut * ( 1.0 - xi ) );
238 tij[Cache::facetTijIdx] = factor * ( - wOut - wFacet );
239 }
240 else
241 {
242 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
243 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
244 tij[Cache::outsideTijIdx] = 0.0;
245 }
246 }
247 else
248 {
249 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
250 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
251 tij[Cache::outsideTijIdx] = 0.0;
252 }
253 }
254 else if (iBcTypes.hasOnlyDirichlet())
255 {
256 tij[Cache::insideTijIdx] = wIn;
257 tij[Cache::outsideTijIdx] = 0.0;
258 tij[Cache::facetTijIdx] = -wIn;
259 }
260 else
261 DUNE_THROW(Dune::NotImplemented, "Interior boundary types other than pure Dirichlet or Neumann");
262
263 return tij;
264 }
265};
266
271template<class TypeTag>
272class CCTpfaFacetCouplingFouriersLawImpl<TypeTag, /*isNetwork*/true>
273: public FouriersLawImplementation<TypeTag, DiscretizationMethod::cctpfa>
274{
277
280 using FVElementGeometry = typename GridGeometry::LocalView;
281 using SubControlVolume = typename GridGeometry::SubControlVolume;
282 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
283 using Extrusion = Extrusion_t<GridGeometry>;
284
285 using GridView = typename GridGeometry::GridView;
286 using Element = typename GridView::template Codim<0>::Entity;
287
291 class FacetCouplingFouriersLawCache
292 {
293 public:
295 using Filler = typename ParentType::Cache::Filler;
296
300 static constexpr int insideTijIdx = 0;
301 static constexpr int facetTijIdx = 1;
302
304 using HeatConductionTransmissibilityContainer = std::array<Scalar, 2>;
305
307 template< class Problem, class ElementVolumeVariables >
308 void updateHeatConduction(const Problem& problem,
309 const Element& element,
310 const FVElementGeometry& fvGeometry,
311 const ElementVolumeVariables& elemVolVars,
312 const SubControlVolumeFace &scvf)
313 {
314 tij_ = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
315 }
316
321 Scalar heatConductionTij() const
322 { return tij_[insideTijIdx]; }
323
325 Scalar heatConductionTijInside() const
326 { return tij_[insideTijIdx]; }
327
329 Scalar heatConductionTijFacet() const
330 { return tij_[facetTijIdx]; }
331
332 private:
333 HeatConductionTransmissibilityContainer tij_;
334 };
335
336public:
338 using Cache = FacetCouplingFouriersLawCache;
339
342
344 template< class Problem, class ElementVolumeVariables, class ElementFluxVarsCache >
345 static Scalar flux(const Problem& problem,
346 const Element& element,
347 const FVElementGeometry& fvGeometry,
348 const ElementVolumeVariables& elemVolVars,
349 const SubControlVolumeFace& scvf,
350 const ElementFluxVarsCache& elemFluxVarsCache)
351 {
352 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
353 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
354
355 // get inside/facet volume variables
356 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
357 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
358 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
359
360 return fluxVarsCache.heatConductionTijInside()*insideVolVars.temperature()
361 + fluxVarsCache.heatConductionTijFacet()*facetVolVars.temperature();
362 }
363
364 // The flux variables cache has to be bound to an element prior to flux calculations
365 // During the binding, the transmissibility will be computed and stored using the method below.
366 template< class Problem, class ElementVolumeVariables >
367 static typename Cache::HeatConductionTransmissibilityContainer
368 calculateTransmissibility(const Problem& problem,
369 const Element& element,
370 const FVElementGeometry& fvGeometry,
371 const ElementVolumeVariables& elemVolVars,
372 const SubControlVolumeFace& scvf)
373 {
374 typename Cache::HeatConductionTransmissibilityContainer tij;
375 if (!problem.couplingManager().isCoupled(element, scvf))
376 {
378 tij[Cache::insideTijIdx] = ParentType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
379 return tij;
380 }
381
383 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(), "FacetCoupling.Xi", 1.0);
384
385 // On surface grids only xi = 1.0 can be used, as the coupling condition
386 // for xi != 1.0 does not generalize for surface grids where the normal
387 // vectors of the inside/outside elements have different orientations.
388 if (Dune::FloatCmp::ne(xi, 1.0, 1e-6))
389 DUNE_THROW(Dune::InvalidStateException, "Xi != 1.0 cannot be used on surface grids");
390
391 const auto insideScvIdx = scvf.insideScvIdx();
392 const auto& insideScv = fvGeometry.scv(insideScvIdx);
393 const auto& insideVolVars = elemVolVars[insideScvIdx];
394 const auto wIn = Extrusion::area(scvf)
395 *computeTpfaTransmissibility(scvf, insideScv,
396 insideVolVars.effectiveThermalConductivity(),
397 insideVolVars.extrusionFactor());
398
399 // proceed depending on the interior BC types used
400 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
401
402 // neumann-coupling
403 if (iBcTypes.hasOnlyNeumann())
404 {
405 // Here we use the square root of the facet extrusion factor
406 // as an approximate average distance from scvf ip to facet center
407 using std::sqrt;
408 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
409 const auto wFacet = 2.0*Extrusion::area(scvf)*insideVolVars.extrusionFactor()
410 /sqrt(facetVolVars.extrusionFactor())
411 *vtmv(scvf.unitOuterNormal(),
412 facetVolVars.effectiveThermalConductivity(),
413 scvf.unitOuterNormal());
414
415 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
416 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
417 }
418 else if (iBcTypes.hasOnlyDirichlet())
419 {
420 tij[Cache::insideTijIdx] = wIn;
421 tij[Cache::facetTijIdx] = -wIn;
422 }
423 else
424 DUNE_THROW(Dune::NotImplemented, "Interior boundary types other than pure Dirichlet or Neumann");
425
426 return tij;
427 }
428};
429
430} // end namespace Dumux
431
432#endif
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
Helper classes to compute the integration elements.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
Tensor::field_type computeTpfaTransmissibility(const SubControlVolumeFace &scvf, const SubControlVolume &scv, const Tensor &T, typename 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:47
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:849
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
forward declaration of the method-specific implementation
Definition: flux/fourierslaw.hh:37
Fourier's law for cell-centered finite volume schemes with two-point flux approximation.
Definition: flux/cctpfa/fourierslaw.hh:46
TpfaFouriersLawCache Cache
export the type for the corresponding cache
Definition: flux/cctpfa/fourierslaw.hh:109
Forward declaration of the implementation.
Definition: multidomain/facet/cellcentered/tpfa/fourierslaw.hh:48
Specialization of CCTpfaFacetCouplingFouriersLawImpl for dim=dimWorld.
Definition: multidomain/facet/cellcentered/tpfa/fourierslaw.hh:68
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:177
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:144
Specialization of CCTpfaFacetCouplingFouriersLawImpl for dim<dimWorld.
Definition: multidomain/facet/cellcentered/tpfa/fourierslaw.hh:274
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:345
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:368
Declares all properties used in Dumux.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...
Fourier's law for cell-centered finite volume schemes with two-point flux approximation.