version 3.10-dev
flux/cctpfa/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_FOURIERS_LAW_HH
13#define DUMUX_DISCRETIZATION_CC_TPFA_FOURIERS_LAW_HH
14
17
21
22namespace Dumux {
23
24// forward declaration
25template<class TypeTag, class DiscretizationMethod>
26class FouriersLawImplementation;
27
32template <class TypeTag>
33class FouriersLawImplementation<TypeTag, DiscretizationMethods::CCTpfa>
34{
39 using FVElementGeometry = typename GridGeometry::LocalView;
40 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
41 using Extrusion = Extrusion_t<GridGeometry>;
42 using GridView = typename GridGeometry::GridView;
43 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
44 using VolumeVariables = typename ElementVolumeVariables::VolumeVariables;
45 using Element = typename GridView::template Codim<0>::Entity;
47 using ElementFluxVarsCache = typename GridFluxVariablesCache::LocalView;
48 using FluxVariablesCache = typename GridFluxVariablesCache::FluxVariablesCache;
49
50 static const int dim = GridView::dimension;
51 static const int dimWorld = GridView::dimensionworld;
52
54 class TpfaFouriersLawCacheFiller
55 {
56 public:
59 template<class FluxVariablesCacheFiller>
60 static void fill(FluxVariablesCache& scvfFluxVarsCache,
61 const Problem& problem,
62 const Element& element,
63 const FVElementGeometry& fvGeometry,
64 const ElementVolumeVariables& elemVolVars,
65 const SubControlVolumeFace& scvf,
66 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
67 {
68 scvfFluxVarsCache.updateHeatConduction(problem, element, fvGeometry, elemVolVars, scvf);
69 }
70 };
71
73 class TpfaFouriersLawCache
74 {
75 public:
76 using Filler = TpfaFouriersLawCacheFiller;
77
78 void updateHeatConduction(const Problem& problem,
79 const Element& element,
80 const FVElementGeometry& fvGeometry,
81 const ElementVolumeVariables& elemVolVars,
82 const SubControlVolumeFace &scvf)
83 {
84 tij_ = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
85 }
86
87 const Scalar& heatConductionTij() const
88 { return tij_; }
89
90 private:
91 Scalar tij_;
92 };
93
94public:
97 static constexpr DiscretizationMethod discMethod{};
98
100 using Cache = TpfaFouriersLawCache;
101
113 static Scalar flux(const Problem& problem,
114 const Element& element,
115 const FVElementGeometry& fvGeometry,
116 const VolumeVariables& insideVolVars,
117 const VolumeVariables& outsideVolVars,
118 const SubControlVolumeFace& scvf,
119 const ElementFluxVarsCache& elemFluxVarsCache)
120 {
121 if constexpr (isMixedDimensional_)
122 if (scvf.numOutsideScv() != 1)
123 DUNE_THROW(Dune::Exception, "This flux overload requires scvf.numOutsideScv() == 1");
124
125 // heat conductivities are always solution dependent (?)
126 Scalar tij = elemFluxVarsCache[scvf].heatConductionTij();
127
128 // get the inside/outside temperatures
129 const auto tInside = insideVolVars.temperature();
130 const auto tOutside = outsideVolVars.temperature();
131
132 return tij*(tInside - tOutside);
133 }
134
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 // heat conductivities are always solution dependent (?)
150 Scalar tij = elemFluxVarsCache[scvf].heatConductionTij();
151
152 // get the inside/outside temperatures
153 const auto tInside = elemVolVars[scvf.insideScvIdx()].temperature();
154 const auto tOutside = scvf.numOutsideScvs() == 1 ? elemVolVars[scvf.outsideScvIdx()].temperature()
155 : branchingFacetTemperature_(problem, element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf, tInside, tij);
156
157 return tij*(tInside - tOutside);
158 }
159
161 static Scalar calculateTransmissibility(const Problem& problem,
162 const Element& element,
163 const FVElementGeometry& fvGeometry,
164 const ElementVolumeVariables& elemVolVars,
165 const SubControlVolumeFace& scvf)
166 {
167 Scalar tij;
168
169 const auto insideScvIdx = scvf.insideScvIdx();
170 const auto& insideScv = fvGeometry.scv(insideScvIdx);
171 const auto& insideVolVars = elemVolVars[insideScvIdx];
172
173 const auto insideLambda = insideVolVars.effectiveThermalConductivity();
174 const Scalar ti = computeTpfaTransmissibility(fvGeometry, scvf, insideScv, insideLambda, insideVolVars.extrusionFactor());
175
176 // for the boundary (dirichlet) or at branching points we only need ti
177 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
178 {
179 tij = Extrusion::area(fvGeometry, scvf)*ti;
180 }
181 // otherwise we compute a tpfa harmonic mean
182 else
183 {
184 const auto outsideScvIdx = scvf.outsideScvIdx();
185 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
186 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
187
188 const auto outsideLambda = outsideVolVars.effectiveThermalConductivity();
189 Scalar tj;
190 if constexpr (dim == dimWorld)
191 // assume the normal vector from outside is anti parallel so we save flipping a vector
192 tj = -1.0*computeTpfaTransmissibility(fvGeometry, scvf, outsideScv, outsideLambda, outsideVolVars.extrusionFactor());
193 else
194 tj = computeTpfaTransmissibility(fvGeometry, fvGeometry.flipScvf(scvf.index()), outsideScv, outsideLambda, outsideVolVars.extrusionFactor());
195
196 // check for division by zero!
197 if (ti*tj <= 0.0)
198 tij = 0;
199 else
200 tij = Extrusion::area(fvGeometry, scvf)*(ti * tj)/(ti + tj);
201 }
202
203 return tij;
204 }
205
206private:
207
209 static Scalar branchingFacetTemperature_(const Problem& problem,
210 const Element& element,
211 const FVElementGeometry& fvGeometry,
212 const ElementVolumeVariables& elemVolVars,
213 const ElementFluxVarsCache& elemFluxVarsCache,
214 const SubControlVolumeFace& scvf,
215 Scalar insideTemperature,
216 Scalar insideTi)
217 {
218 Scalar sumTi(insideTi);
219 Scalar sumTempTi(insideTi*insideTemperature);
220
221 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
222 {
223 const auto outsideScvIdx = scvf.outsideScvIdx(i);
224 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
225 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
226 const auto& outsideFluxVarsCache = elemFluxVarsCache[flippedScvf];
227
228 auto outsideTi = outsideFluxVarsCache.heatConductionTij();
229 sumTi += outsideTi;
230 sumTempTi += outsideTi*outsideVolVars.temperature();
231 }
232 return sumTempTi/sumTi;
233 }
234
235 static constexpr bool isMixedDimensional_ = static_cast<int>(GridView::dimension) < static_cast<int>(GridView::dimensionworld);
236};
237
238} // end namespace Dumux
239
240#endif
Fourier's law for cell-centered finite volume schemes with two-point flux approximation.
Definition: flux/cctpfa/fourierslaw.hh:34
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVarsCache &elemFluxVarsCache)
Returns the heat flux within the porous medium (in J/s) across the given sub-control volume face.
Definition: flux/cctpfa/fourierslaw.hh:142
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const VolumeVariables &insideVolVars, const VolumeVariables &outsideVolVars, const SubControlVolumeFace &scvf, const ElementFluxVarsCache &elemFluxVarsCache)
Returns the heat flux within the porous medium (in J/s) across the given sub-control volume face.
Definition: flux/cctpfa/fourierslaw.hh:113
static Scalar calculateTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Compute transmissibilities.
Definition: flux/cctpfa/fourierslaw.hh:161
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.
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
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
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...