3.2-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
39
41
42namespace Dumux {
43
45template<class TypeTag,
46 bool isNetwork>
48
55template<class TypeTag>
59
64template<class TypeTag>
65class CCTpfaFacetCouplingFouriersLawImpl<TypeTag, /*isNetwork*/false>
66: public FouriersLawImplementation<TypeTag, DiscretizationMethod::cctpfa>
67{
70
73 using FVElementGeometry = typename GridGeometry::LocalView;
74 using SubControlVolume = typename GridGeometry::SubControlVolume;
75 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
76
77 using GridView = typename GridGeometry::GridView;
78 using Element = typename GridView::template Codim<0>::Entity;
79
83 class FacetCouplingFouriersLawCache
84 {
85 public:
87 using Filler = typename ParentType::Cache::Filler;
88
92 static constexpr int insideTijIdx = 0;
93 static constexpr int outsideTijIdx = 1;
94 static constexpr int facetTijIdx = 2;
95
97 using HeatConductionTransmissibilityContainer = std::array<Scalar, 3>;
98
100 template< class Problem, class ElementVolumeVariables >
101 void updateHeatConduction(const Problem& problem,
102 const Element& element,
103 const FVElementGeometry& fvGeometry,
104 const ElementVolumeVariables& elemVolVars,
105 const SubControlVolumeFace &scvf)
106 {
107 tij_ = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
108 }
109
114 Scalar heatConductionTij() const
115 { return tij_[insideTijIdx]; }
116
118 Scalar heatConductionTijInside() const
119 { return tij_[insideTijIdx]; }
120
122 Scalar heatConductionTijOutside() const
123 { return tij_[outsideTijIdx]; }
124
126 Scalar heatConductionTijFacet() const
127 { return tij_[facetTijIdx]; }
128
129 private:
130 HeatConductionTransmissibilityContainer tij_;
131 };
132
133public:
135 using Cache = FacetCouplingFouriersLawCache;
136
139
141 template< class Problem, class ElementVolumeVariables, class ElementFluxVarsCache >
142 static Scalar flux(const Problem& problem,
143 const Element& element,
144 const FVElementGeometry& fvGeometry,
145 const ElementVolumeVariables& elemVolVars,
146 const SubControlVolumeFace& scvf,
147 const ElementFluxVarsCache& elemFluxVarsCache)
148 {
149 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
150 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
151
152 // get inside/outside volume variables
153 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
154 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
155 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
156 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
157
158 // the inside and outside temperatures
159 const Scalar tInside = insideVolVars.temperature();
160 const Scalar tFacet = facetVolVars.temperature();
161 const Scalar tOutside = outsideVolVars.temperature();
162
163 Scalar flux = fluxVarsCache.heatConductionTijInside()*tInside
164 + fluxVarsCache.heatConductionTijFacet()*tFacet;
165
166 if (!scvf.boundary())
167 flux += fluxVarsCache.heatConductionTijOutside()*tOutside;
168 return flux;
169 }
170
171 // The flux variables cache has to be bound to an element prior to flux calculations
172 // During the binding, the transmissibility will be computed and stored using the method below.
173 template< class Problem, class ElementVolumeVariables >
174 static typename Cache::HeatConductionTransmissibilityContainer
175 calculateTransmissibility(const Problem& problem,
176 const Element& element,
177 const FVElementGeometry& fvGeometry,
178 const ElementVolumeVariables& elemVolVars,
179 const SubControlVolumeFace& scvf)
180 {
181 typename Cache::HeatConductionTransmissibilityContainer tij;
182 if (!problem.couplingManager().isCoupled(element, scvf))
183 {
185 tij[Cache::insideTijIdx] = ParentType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
186 return tij;
187 }
188
190 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(), "FacetCoupling.Xi", 1.0);
191
192 const auto insideScvIdx = scvf.insideScvIdx();
193 const auto& insideScv = fvGeometry.scv(insideScvIdx);
194 const auto& insideVolVars = elemVolVars[insideScvIdx];
195 const auto wIn = scvf.area()*computeTpfaTransmissibility(scvf,
196 insideScv,
197 insideVolVars.effectiveThermalConductivity(),
198 insideVolVars.extrusionFactor());
199
200 // proceed depending on the interior BC types used
201 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
202
203 // neumann-coupling
204 if (iBcTypes.hasOnlyNeumann())
205 {
206 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
207 const auto wFacet = 2.0*scvf.area()*insideVolVars.extrusionFactor()
208 /facetVolVars.extrusionFactor()
209 *vtmv(scvf.unitOuterNormal(),
210 facetVolVars.effectiveThermalConductivity(),
211 scvf.unitOuterNormal());
212
213 // The fluxes across this face and the outside face can be expressed in matrix form:
214 // \f$\mathbf{C} \bar{\mathbf{u}} + \mathbf{D} \mathbf{u} + \mathbf{E} \mathbf{u}_\gamma\f$,
215 // where \f$\gamma$\f denotes the domain living on the facets and \f$\bar{\mathbf{u}}$\f are
216 // intermediate face unknowns in the matrix domain. Equivalently, flux continuity reads:
217 // \f$\mathbf{A} \bar{\mathbf{u}} = \mathbf{B} \mathbf{u} + \mathbf{M} \mathbf{u}_\gamma\f$.
218 // Combining the two, we can eliminate the intermediate unknowns and compute the transmissibilities
219 // that allow the description of the fluxes as functions of the cell and Dirichlet temperatures only.
220 if (!scvf.boundary())
221 {
222 const auto outsideScvIdx = scvf.outsideScvIdx();
223 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
224 const auto wOut = -1.0*scvf.area()*computeTpfaTransmissibility(scvf,
225 fvGeometry.scv(outsideScvIdx),
226 outsideVolVars.effectiveThermalConductivity(),
227 outsideVolVars.extrusionFactor());
228
229 if ( !Dune::FloatCmp::eq(xi, 1.0, 1e-6) )
230 {
231 // optimized implementation: factorization obtained using sympy
232 // see CCTpfaFacetCouplingDarcysLaw for more details
233 const Scalar factor = wIn * wFacet / ( wIn * wOut * ( 2.0 * xi - 1.0 ) + wFacet * ( xi * ( wIn + wOut ) + wFacet ) );
234 tij[Cache::insideTijIdx] = factor * ( wOut * xi + wFacet );
235 tij[Cache::outsideTijIdx] = factor * ( wOut * ( 1.0 - xi ) );
236 tij[Cache::facetTijIdx] = factor * ( - wOut - wFacet );
237 }
238 else
239 {
240 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
241 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
242 tij[Cache::outsideTijIdx] = 0.0;
243 }
244 }
245 else
246 {
247 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
248 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
249 tij[Cache::outsideTijIdx] = 0.0;
250 }
251 }
252 else if (iBcTypes.hasOnlyDirichlet())
253 {
254 tij[Cache::insideTijIdx] = wIn;
255 tij[Cache::outsideTijIdx] = 0.0;
256 tij[Cache::facetTijIdx] = -wIn;
257 }
258 else
259 DUNE_THROW(Dune::NotImplemented, "Interior boundary types other than pure Dirichlet or Neumann");
260
261 return tij;
262 }
263};
264
269template<class TypeTag>
270class CCTpfaFacetCouplingFouriersLawImpl<TypeTag, /*isNetwork*/true>
271: public FouriersLawImplementation<TypeTag, DiscretizationMethod::cctpfa>
272{
275
278 using FVElementGeometry = typename GridGeometry::LocalView;
279 using SubControlVolume = typename GridGeometry::SubControlVolume;
280 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
281
282 using GridView = typename GridGeometry::GridView;
283 using Element = typename GridView::template Codim<0>::Entity;
284
288 class FacetCouplingFouriersLawCache
289 {
290 public:
292 using Filler = typename ParentType::Cache::Filler;
293
297 static constexpr int insideTijIdx = 0;
298 static constexpr int facetTijIdx = 1;
299
301 using HeatConductionTransmissibilityContainer = std::array<Scalar, 2>;
302
304 template< class Problem, class ElementVolumeVariables >
305 void updateHeatConduction(const Problem& problem,
306 const Element& element,
307 const FVElementGeometry& fvGeometry,
308 const ElementVolumeVariables& elemVolVars,
309 const SubControlVolumeFace &scvf)
310 {
311 tij_ = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
312 }
313
318 Scalar heatConductionTij() const
319 { return tij_[insideTijIdx]; }
320
322 Scalar heatConductionTijInside() const
323 { return tij_[insideTijIdx]; }
324
326 Scalar heatConductionTijFacet() const
327 { return tij_[facetTijIdx]; }
328
329 private:
330 HeatConductionTransmissibilityContainer tij_;
331 };
332
333public:
335 using Cache = FacetCouplingFouriersLawCache;
336
339
341 template< class Problem, class ElementVolumeVariables, class ElementFluxVarsCache >
342 static Scalar flux(const Problem& problem,
343 const Element& element,
344 const FVElementGeometry& fvGeometry,
345 const ElementVolumeVariables& elemVolVars,
346 const SubControlVolumeFace& scvf,
347 const ElementFluxVarsCache& elemFluxVarsCache)
348 {
349 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
350 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
351
352 // get inside/facet volume variables
353 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
354 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
355 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
356
357 return fluxVarsCache.heatConductionTijInside()*insideVolVars.temperature()
358 + fluxVarsCache.heatConductionTijFacet()*facetVolVars.temperature();
359 }
360
361 // The flux variables cache has to be bound to an element prior to flux calculations
362 // During the binding, the transmissibility will be computed and stored using the method below.
363 template< class Problem, class ElementVolumeVariables >
364 static typename Cache::HeatConductionTransmissibilityContainer
365 calculateTransmissibility(const Problem& problem,
366 const Element& element,
367 const FVElementGeometry& fvGeometry,
368 const ElementVolumeVariables& elemVolVars,
369 const SubControlVolumeFace& scvf)
370 {
371 typename Cache::HeatConductionTransmissibilityContainer tij;
372 if (!problem.couplingManager().isCoupled(element, scvf))
373 {
375 tij[Cache::insideTijIdx] = ParentType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
376 return tij;
377 }
378
380 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(), "FacetCoupling.Xi", 1.0);
381
382 // On surface grids only xi = 1.0 can be used, as the coupling condition
383 // for xi != 1.0 does not generalize for surface grids where the normal
384 // vectors of the inside/outside elements have different orientations.
385 if (Dune::FloatCmp::ne(xi, 1.0, 1e-6))
386 DUNE_THROW(Dune::InvalidStateException, "Xi != 1.0 cannot be used on surface grids");
387
388 const auto insideScvIdx = scvf.insideScvIdx();
389 const auto& insideScv = fvGeometry.scv(insideScvIdx);
390 const auto& insideVolVars = elemVolVars[insideScvIdx];
391 const auto wIn = scvf.area()*computeTpfaTransmissibility(scvf,
392 insideScv,
393 insideVolVars.effectiveThermalConductivity(),
394 insideVolVars.extrusionFactor());
395
396 // proceed depending on the interior BC types used
397 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
398
399 // neumann-coupling
400 if (iBcTypes.hasOnlyNeumann())
401 {
402 // Here we use the square root of the facet extrusion factor
403 // as an approximate average distance from scvf ip to facet center
404 using std::sqrt;
405 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
406 const auto wFacet = 2.0*scvf.area()*insideVolVars.extrusionFactor()
407 /sqrt(facetVolVars.extrusionFactor())
408 *vtmv(scvf.unitOuterNormal(),
409 facetVolVars.effectiveThermalConductivity(),
410 scvf.unitOuterNormal());
411
412 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
413 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
414 }
415 else if (iBcTypes.hasOnlyDirichlet())
416 {
417 tij[Cache::insideTijIdx] = wIn;
418 tij[Cache::facetTijIdx] = -wIn;
419 }
420 else
421 DUNE_THROW(Dune::NotImplemented, "Interior boundary types other than pure Dirichlet or Neumann");
422
423 return tij;
424 }
425};
426
427} // end namespace Dumux
428
429#endif
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
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:840
Definition: adapt.hh:29
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:47
TpfaFouriersLawCache Cache
export the type for the corresponding cache
Definition: flux/cctpfa/fourierslaw.hh:108
Forward declaration of the implementation.
Definition: multidomain/facet/cellcentered/tpfa/fourierslaw.hh:47
Specialization of CCTpfaFacetCouplingFouriersLawImpl for dim=dimWorld.
Definition: multidomain/facet/cellcentered/tpfa/fourierslaw.hh:67
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:175
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:142
Specialization of CCTpfaFacetCouplingFouriersLawImpl for dim<dimWorld.
Definition: multidomain/facet/cellcentered/tpfa/fourierslaw.hh:272
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:342
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:365
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.