3.5-git
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, class DiscretizationMethod>
38class FouriersLawImplementation;
39
44template <class TypeTag>
45class FouriersLawImplementation<TypeTag, DiscretizationMethods::CCTpfa>
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 static constexpr DiscretizationMethod discMethod{};
109
111 using Cache = TpfaFouriersLawCache;
112
124 static Scalar flux(const Problem& problem,
125 const Element& element,
126 const FVElementGeometry& fvGeometry,
127 const VolumeVariables& insideVolVars,
128 const VolumeVariables& outsideVolVars,
129 const SubControlVolumeFace& scvf,
130 const ElementFluxVarsCache& elemFluxVarsCache)
131 {
132 if constexpr (isMixedDimensional_)
133 if (scvf.numOutsideScv() != 1)
134 DUNE_THROW(Dune::Exception, "This flux overload requires scvf.numOutsideScv() == 1");
135
136 // heat conductivities are always solution dependent (?)
137 Scalar tij = elemFluxVarsCache[scvf].heatConductionTij();
138
139 // get the inside/outside temperatures
140 const auto tInside = insideVolVars.temperature();
141 const auto tOutside = outsideVolVars.temperature();
142
143 return tij*(tInside - tOutside);
144 }
145
153 static Scalar flux(const Problem& problem,
154 const Element& element,
155 const FVElementGeometry& fvGeometry,
156 const ElementVolumeVariables& elemVolVars,
157 const SubControlVolumeFace& scvf,
158 const ElementFluxVarsCache& elemFluxVarsCache)
159 {
160 // heat conductivities are always solution dependent (?)
161 Scalar tij = elemFluxVarsCache[scvf].heatConductionTij();
162
163 // get the inside/outside temperatures
164 const auto tInside = elemVolVars[scvf.insideScvIdx()].temperature();
165 const auto tOutside = scvf.numOutsideScvs() == 1 ? elemVolVars[scvf.outsideScvIdx()].temperature()
166 : branchingFacetTemperature_(problem, element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf, tInside, tij);
167
168 return tij*(tInside - tOutside);
169 }
170
172 static Scalar calculateTransmissibility(const Problem& problem,
173 const Element& element,
174 const FVElementGeometry& fvGeometry,
175 const ElementVolumeVariables& elemVolVars,
176 const SubControlVolumeFace& scvf)
177 {
178 Scalar tij;
179
180 const auto insideScvIdx = scvf.insideScvIdx();
181 const auto& insideScv = fvGeometry.scv(insideScvIdx);
182 const auto& insideVolVars = elemVolVars[insideScvIdx];
183
184 const auto insideLambda = insideVolVars.effectiveThermalConductivity();
185 const Scalar ti = computeTpfaTransmissibility(scvf, insideScv, insideLambda, insideVolVars.extrusionFactor());
186
187 // for the boundary (dirichlet) or at branching points we only need ti
188 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
189 {
190 tij = Extrusion::area(scvf)*ti;
191 }
192 // otherwise we compute a tpfa harmonic mean
193 else
194 {
195 const auto outsideScvIdx = scvf.outsideScvIdx();
196 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
197 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
198
199 const auto outsideLambda = outsideVolVars.effectiveThermalConductivity();
200 Scalar tj;
201 if constexpr (dim == dimWorld)
202 // assume the normal vector from outside is anti parallel so we save flipping a vector
203 tj = -1.0*computeTpfaTransmissibility(scvf, outsideScv, outsideLambda, outsideVolVars.extrusionFactor());
204 else
205 tj = computeTpfaTransmissibility(fvGeometry.flipScvf(scvf.index()), outsideScv, outsideLambda, outsideVolVars.extrusionFactor());
206
207 // check for division by zero!
208 if (ti*tj <= 0.0)
209 tij = 0;
210 else
211 tij = Extrusion::area(scvf)*(ti * tj)/(ti + tj);
212 }
213
214 return tij;
215 }
216
217private:
218
220 static Scalar branchingFacetTemperature_(const Problem& problem,
221 const Element& element,
222 const FVElementGeometry& fvGeometry,
223 const ElementVolumeVariables& elemVolVars,
224 const ElementFluxVarsCache& elemFluxVarsCache,
225 const SubControlVolumeFace& scvf,
226 Scalar insideTemperature,
227 Scalar insideTi)
228 {
229 Scalar sumTi(insideTi);
230 Scalar sumTempTi(insideTi*insideTemperature);
231
232 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
233 {
234 const auto outsideScvIdx = scvf.outsideScvIdx(i);
235 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
236 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
237 const auto& outsideFluxVarsCache = elemFluxVarsCache[flippedScvf];
238
239 auto outsideTi = outsideFluxVarsCache.heatConductionTij();
240 sumTi += outsideTi;
241 sumTempTi += outsideTi*outsideVolVars.temperature();
242 }
243 return sumTempTi/sumTi;
244 }
245
246 static constexpr bool isMixedDimensional_ = static_cast<int>(GridView::dimension) < static_cast<int>(GridView::dimensionworld);
247};
248
249} // end namespace Dumux
250
251#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
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: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
Definition: method.hh:49
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
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:153
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:124
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:172
TpfaFouriersLawCache Cache
export the type for the corresponding cache
Definition: flux/cctpfa/fourierslaw.hh:111
Declares all properties used in Dumux.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...