3.6-git
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, DiscretizationMethods::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 static constexpr DiscretizationMethod discMethod{};
142
144 template< class Problem, class ElementVolumeVariables, class ElementFluxVarsCache >
145 static Scalar flux(const Problem& problem,
146 const Element& element,
147 const FVElementGeometry& fvGeometry,
148 const ElementVolumeVariables& elemVolVars,
149 const SubControlVolumeFace& scvf,
150 const ElementFluxVarsCache& elemFluxVarsCache)
151 {
152 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
153 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
154
155 // get inside/outside volume variables
156 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
157 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
158 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
159 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
160
161 // the inside and outside temperatures
162 const Scalar tInside = insideVolVars.temperature();
163 const Scalar tFacet = facetVolVars.temperature();
164 const Scalar tOutside = outsideVolVars.temperature();
165
166 Scalar flux = fluxVarsCache.heatConductionTijInside()*tInside
167 + fluxVarsCache.heatConductionTijFacet()*tFacet;
168
169 if (!scvf.boundary())
170 flux += fluxVarsCache.heatConductionTijOutside()*tOutside;
171 return flux;
172 }
173
174 // The flux variables cache has to be bound to an element prior to flux calculations
175 // During the binding, the transmissibility will be computed and stored using the method below.
176 template< class Problem, class ElementVolumeVariables >
177 static typename Cache::HeatConductionTransmissibilityContainer
178 calculateTransmissibility(const Problem& problem,
179 const Element& element,
180 const FVElementGeometry& fvGeometry,
181 const ElementVolumeVariables& elemVolVars,
182 const SubControlVolumeFace& scvf)
183 {
184 typename Cache::HeatConductionTransmissibilityContainer tij;
185 if (!problem.couplingManager().isCoupled(element, scvf))
186 {
188 tij[Cache::insideTijIdx] = ParentType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
189 return tij;
190 }
191
193 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(), "FacetCoupling.Xi", 1.0);
194
195 const auto insideScvIdx = scvf.insideScvIdx();
196 const auto& insideScv = fvGeometry.scv(insideScvIdx);
197 const auto& insideVolVars = elemVolVars[insideScvIdx];
198 const auto wIn = Extrusion::area(fvGeometry, scvf)
199 *computeTpfaTransmissibility(fvGeometry, scvf, insideScv,
200 insideVolVars.effectiveThermalConductivity(),
201 insideVolVars.extrusionFactor());
202
203 // proceed depending on the interior BC types used
204 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
205
206 // neumann-coupling
207 if (iBcTypes.hasOnlyNeumann())
208 {
209 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
210 const auto wFacet = 2.0*Extrusion::area(fvGeometry, scvf)*insideVolVars.extrusionFactor()
211 /facetVolVars.extrusionFactor()
212 *vtmv(scvf.unitOuterNormal(),
213 facetVolVars.effectiveThermalConductivity(),
214 scvf.unitOuterNormal());
215
216 // The fluxes across this face and the outside face can be expressed in matrix form:
217 // \f$\mathbf{C} \bar{\mathbf{u}} + \mathbf{D} \mathbf{u} + \mathbf{E} \mathbf{u}_\gamma\f$,
218 // where \f$\gamma$\f denotes the domain living on the facets and \f$\bar{\mathbf{u}}$\f are
219 // intermediate face unknowns in the matrix domain. Equivalently, flux continuity reads:
220 // \f$\mathbf{A} \bar{\mathbf{u}} = \mathbf{B} \mathbf{u} + \mathbf{M} \mathbf{u}_\gamma\f$.
221 // Combining the two, we can eliminate the intermediate unknowns and compute the transmissibilities
222 // that allow the description of the fluxes as functions of the cell and Dirichlet temperatures only.
223 if (!scvf.boundary())
224 {
225 const auto outsideScvIdx = scvf.outsideScvIdx();
226 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
227 const auto wOut = -1.0*Extrusion::area(fvGeometry, scvf)
228 *computeTpfaTransmissibility(fvGeometry, scvf, fvGeometry.scv(outsideScvIdx),
229 outsideVolVars.effectiveThermalConductivity(),
230 outsideVolVars.extrusionFactor());
231
232 if ( !Dune::FloatCmp::eq(xi, 1.0, 1e-6) )
233 {
234 // optimized implementation: factorization obtained using sympy
235 // see CCTpfaFacetCouplingDarcysLaw for more details
236 const Scalar factor = wIn * wFacet / ( wIn * wOut * ( 2.0 * xi - 1.0 ) + wFacet * ( xi * ( wIn + wOut ) + wFacet ) );
237 tij[Cache::insideTijIdx] = factor * ( wOut * xi + wFacet );
238 tij[Cache::outsideTijIdx] = factor * ( wOut * ( 1.0 - xi ) );
239 tij[Cache::facetTijIdx] = factor * ( - wOut - wFacet );
240 }
241 else
242 {
243 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
244 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
245 tij[Cache::outsideTijIdx] = 0.0;
246 }
247 }
248 else
249 {
250 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
251 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
252 tij[Cache::outsideTijIdx] = 0.0;
253 }
254 }
255 else if (iBcTypes.hasOnlyDirichlet())
256 {
257 tij[Cache::insideTijIdx] = wIn;
258 tij[Cache::outsideTijIdx] = 0.0;
259 tij[Cache::facetTijIdx] = -wIn;
260 }
261 else
262 DUNE_THROW(Dune::NotImplemented, "Interior boundary types other than pure Dirichlet or Neumann");
263
264 return tij;
265 }
266};
267
272template<class TypeTag>
273class CCTpfaFacetCouplingFouriersLawImpl<TypeTag, /*isNetwork*/true>
274: public FouriersLawImplementation<TypeTag, DiscretizationMethods::CCTpfa>
275{
278
281 using FVElementGeometry = typename GridGeometry::LocalView;
282 using SubControlVolume = typename GridGeometry::SubControlVolume;
283 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
284 using Extrusion = Extrusion_t<GridGeometry>;
285
286 using GridView = typename GridGeometry::GridView;
287 using Element = typename GridView::template Codim<0>::Entity;
288
292 class FacetCouplingFouriersLawCache
293 {
294 public:
296 using Filler = typename ParentType::Cache::Filler;
297
301 static constexpr int insideTijIdx = 0;
302 static constexpr int facetTijIdx = 1;
303
305 using HeatConductionTransmissibilityContainer = std::array<Scalar, 2>;
306
308 template< class Problem, class ElementVolumeVariables >
309 void updateHeatConduction(const Problem& problem,
310 const Element& element,
311 const FVElementGeometry& fvGeometry,
312 const ElementVolumeVariables& elemVolVars,
313 const SubControlVolumeFace &scvf)
314 {
315 tij_ = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
316 }
317
322 Scalar heatConductionTij() const
323 { return tij_[insideTijIdx]; }
324
326 Scalar heatConductionTijInside() const
327 { return tij_[insideTijIdx]; }
328
330 Scalar heatConductionTijFacet() const
331 { return tij_[facetTijIdx]; }
332
333 private:
334 HeatConductionTransmissibilityContainer tij_;
335 };
336
337public:
339 using Cache = FacetCouplingFouriersLawCache;
340
343 static constexpr DiscretizationMethod discMethod{};
344
346 template< class Problem, class ElementVolumeVariables, class ElementFluxVarsCache >
347 static Scalar flux(const Problem& problem,
348 const Element& element,
349 const FVElementGeometry& fvGeometry,
350 const ElementVolumeVariables& elemVolVars,
351 const SubControlVolumeFace& scvf,
352 const ElementFluxVarsCache& elemFluxVarsCache)
353 {
354 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
355 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
356
357 // get inside/facet volume variables
358 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
359 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
360 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
361
362 return fluxVarsCache.heatConductionTijInside()*insideVolVars.temperature()
363 + fluxVarsCache.heatConductionTijFacet()*facetVolVars.temperature();
364 }
365
366 // The flux variables cache has to be bound to an element prior to flux calculations
367 // During the binding, the transmissibility will be computed and stored using the method below.
368 template< class Problem, class ElementVolumeVariables >
369 static typename Cache::HeatConductionTransmissibilityContainer
370 calculateTransmissibility(const Problem& problem,
371 const Element& element,
372 const FVElementGeometry& fvGeometry,
373 const ElementVolumeVariables& elemVolVars,
374 const SubControlVolumeFace& scvf)
375 {
376 typename Cache::HeatConductionTransmissibilityContainer tij;
377 if (!problem.couplingManager().isCoupled(element, scvf))
378 {
380 tij[Cache::insideTijIdx] = ParentType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
381 return tij;
382 }
383
385 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(), "FacetCoupling.Xi", 1.0);
386
387 // On surface grids only xi = 1.0 can be used, as the coupling condition
388 // for xi != 1.0 does not generalize for surface grids where the normal
389 // vectors of the inside/outside elements have different orientations.
390 if (Dune::FloatCmp::ne(xi, 1.0, 1e-6))
391 DUNE_THROW(Dune::InvalidStateException, "Xi != 1.0 cannot be used on surface grids");
392
393 const auto insideScvIdx = scvf.insideScvIdx();
394 const auto& insideScv = fvGeometry.scv(insideScvIdx);
395 const auto& insideVolVars = elemVolVars[insideScvIdx];
396 const auto wIn = Extrusion::area(fvGeometry, scvf)
397 *computeTpfaTransmissibility(fvGeometry, scvf, insideScv,
398 insideVolVars.effectiveThermalConductivity(),
399 insideVolVars.extrusionFactor());
400
401 // proceed depending on the interior BC types used
402 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
403
404 // neumann-coupling
405 if (iBcTypes.hasOnlyNeumann())
406 {
407 // Here we use the square root of the facet extrusion factor
408 // as an approximate average distance from scvf ip to facet center
409 using std::sqrt;
410 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
411 const auto wFacet = 2.0*Extrusion::area(fvGeometry, scvf)*insideVolVars.extrusionFactor()
412 /sqrt(facetVolVars.extrusionFactor())
413 *vtmv(scvf.unitOuterNormal(),
414 facetVolVars.effectiveThermalConductivity(),
415 scvf.unitOuterNormal());
416
417 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
418 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
419 }
420 else if (iBcTypes.hasOnlyDirichlet())
421 {
422 tij[Cache::insideTijIdx] = wIn;
423 tij[Cache::facetTijIdx] = -wIn;
424 }
425 else
426 DUNE_THROW(Dune::NotImplemented, "Interior boundary types other than pure Dirichlet or Neumann");
427
428 return tij;
429 }
430};
431
432} // end namespace Dumux
433
434#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Define some often used mathematical functions.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
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:48
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:863
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
Definition: method.hh:37
forward declaration of the method-specific implementation
Definition: flux/box/fourierslaw.hh:38
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:112
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:178
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:145
Specialization of CCTpfaFacetCouplingFouriersLawImpl for dim<dimWorld.
Definition: multidomain/facet/cellcentered/tpfa/fourierslaw.hh:275
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:347
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:370
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.