3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * 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_FOURIERS_LAW_HH
25#define DUMUX_DISCRETIZATION_CC_TPFA_FOURIERS_LAW_HH
26
29
33
34namespace Dumux {
35
36// forward declaration
37template<class TypeTag, DiscretizationMethod discMethod>
38class FouriersLawImplementation;
39
44template <class TypeTag>
46{
51 using FVElementGeometry = typename GridGeometry::LocalView;
52 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
53 using Extrusion = Extrusion_t<GridGeometry>;
54 using GridView = typename GridGeometry::GridView;
55 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
56 using VolumeVariables = typename ElementVolumeVariables::VolumeVariables;
57 using Element = typename GridView::template Codim<0>::Entity;
60
61 static const int dim = GridView::dimension;
62 static const int dimWorld = GridView::dimensionworld;
63
65 class TpfaFouriersLawCacheFiller
66 {
67 public:
70 template<class FluxVariablesCacheFiller>
71 static void fill(FluxVariablesCache& scvfFluxVarsCache,
72 const Problem& problem,
73 const Element& element,
74 const FVElementGeometry& fvGeometry,
75 const ElementVolumeVariables& elemVolVars,
76 const SubControlVolumeFace& scvf,
77 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
78 {
79 scvfFluxVarsCache.updateHeatConduction(problem, element, fvGeometry, elemVolVars, scvf);
80 }
81 };
82
84 class TpfaFouriersLawCache
85 {
86 public:
87 using Filler = TpfaFouriersLawCacheFiller;
88
89 void updateHeatConduction(const Problem& problem,
90 const Element& element,
91 const FVElementGeometry& fvGeometry,
92 const ElementVolumeVariables& elemVolVars,
93 const SubControlVolumeFace &scvf)
94 {
95 tij_ = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
96 }
97
98 const Scalar& heatConductionTij() const
99 { return tij_; }
100
101 private:
102 Scalar tij_;
103 };
104
105public:
108
110 using Cache = TpfaFouriersLawCache;
111
123 static Scalar flux(const Problem& problem,
124 const Element& element,
125 const FVElementGeometry& fvGeometry,
126 const VolumeVariables& insideVolVars,
127 const VolumeVariables& outsideVolVars,
128 const SubControlVolumeFace& scvf,
129 const ElementFluxVarsCache& elemFluxVarsCache)
130 {
131 if constexpr (isMixedDimensional_)
132 if (scvf.numOutsideScv() != 1)
133 DUNE_THROW(Dune::Exception, "This flux overload requires scvf.numOutsideScv() == 1");
134
135 // heat conductivities are always solution dependent (?)
136 Scalar tij = elemFluxVarsCache[scvf].heatConductionTij();
137
138 // get the inside/outside temperatures
139 const auto tInside = insideVolVars.temperature();
140 const auto tOutside = outsideVolVars.temperature();
141
142 return tij*(tInside - tOutside);
143 }
144
152 static Scalar flux(const Problem& problem,
153 const Element& element,
154 const FVElementGeometry& fvGeometry,
155 const ElementVolumeVariables& elemVolVars,
156 const SubControlVolumeFace& scvf,
157 const ElementFluxVarsCache& elemFluxVarsCache)
158 {
159 // heat conductivities are always solution dependent (?)
160 Scalar tij = elemFluxVarsCache[scvf].heatConductionTij();
161
162 // get the inside/outside temperatures
163 const auto tInside = elemVolVars[scvf.insideScvIdx()].temperature();
164 const auto tOutside = scvf.numOutsideScvs() == 1 ? elemVolVars[scvf.outsideScvIdx()].temperature()
165 : branchingFacetTemperature_(problem, element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf, tInside, tij);
166
167 return tij*(tInside - tOutside);
168 }
169
171 static Scalar calculateTransmissibility(const Problem& problem,
172 const Element& element,
173 const FVElementGeometry& fvGeometry,
174 const ElementVolumeVariables& elemVolVars,
175 const SubControlVolumeFace& scvf)
176 {
177 Scalar tij;
178
179 const auto insideScvIdx = scvf.insideScvIdx();
180 const auto& insideScv = fvGeometry.scv(insideScvIdx);
181 const auto& insideVolVars = elemVolVars[insideScvIdx];
182
183 const auto insideLambda = insideVolVars.effectiveThermalConductivity();
184 const Scalar ti = computeTpfaTransmissibility(scvf, insideScv, insideLambda, insideVolVars.extrusionFactor());
185
186 // for the boundary (dirichlet) or at branching points we only need ti
187 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
188 {
189 tij = Extrusion::area(scvf)*ti;
190 }
191 // otherwise we compute a tpfa harmonic mean
192 else
193 {
194 const auto outsideScvIdx = scvf.outsideScvIdx();
195 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
196 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
197
198 const auto outsideLambda = outsideVolVars.effectiveThermalConductivity();
199 Scalar tj;
200 if constexpr (dim == dimWorld)
201 // assume the normal vector from outside is anti parallel so we save flipping a vector
202 tj = -1.0*computeTpfaTransmissibility(scvf, outsideScv, outsideLambda, outsideVolVars.extrusionFactor());
203 else
204 tj = computeTpfaTransmissibility(fvGeometry.flipScvf(scvf.index()), outsideScv, outsideLambda, outsideVolVars.extrusionFactor());
205
206 // check for division by zero!
207 if (ti*tj <= 0.0)
208 tij = 0;
209 else
210 tij = Extrusion::area(scvf)*(ti * tj)/(ti + tj);
211 }
212
213 return tij;
214 }
215
216private:
217
219 static Scalar branchingFacetTemperature_(const Problem& problem,
220 const Element& element,
221 const FVElementGeometry& fvGeometry,
222 const ElementVolumeVariables& elemVolVars,
223 const ElementFluxVarsCache& elemFluxVarsCache,
224 const SubControlVolumeFace& scvf,
225 Scalar insideTemperature,
226 Scalar insideTi)
227 {
228 Scalar sumTi(insideTi);
229 Scalar sumTempTi(insideTi*insideTemperature);
230
231 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
232 {
233 const auto outsideScvIdx = scvf.outsideScvIdx(i);
234 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
235 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
236 const auto& outsideFluxVarsCache = elemFluxVarsCache[flippedScvf];
237
238 auto outsideTi = outsideFluxVarsCache.heatConductionTij();
239 sumTi += outsideTi;
240 sumTempTi += outsideTi*outsideVolVars.temperature();
241 }
242 return sumTempTi/sumTi;
243 }
244
245 static constexpr bool isMixedDimensional_ = static_cast<int>(GridView::dimension) < static_cast<int>(GridView::dimensionworld);
246};
247
248} // end namespace Dumux
249
250#endif
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
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
Definition: propertysystem.hh:150
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
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:171
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:152
TpfaFouriersLawCache Cache
export the type for the corresponding cache
Definition: flux/cctpfa/fourierslaw.hh:110
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:123
Declares all properties used in Dumux.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...