3.3.0
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 Element = typename GridView::template Codim<0>::Entity;
59
60 static const int dim = GridView::dimension;
61 static const int dimWorld = GridView::dimensionworld;
62
64 class TpfaFouriersLawCacheFiller
65 {
66 public:
69 template<class FluxVariablesCacheFiller>
70 static void fill(FluxVariablesCache& scvfFluxVarsCache,
71 const Problem& problem,
72 const Element& element,
73 const FVElementGeometry& fvGeometry,
74 const ElementVolumeVariables& elemVolVars,
75 const SubControlVolumeFace& scvf,
76 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
77 {
78 scvfFluxVarsCache.updateHeatConduction(problem, element, fvGeometry, elemVolVars, scvf);
79 }
80 };
81
83 class TpfaFouriersLawCache
84 {
85 public:
86 using Filler = TpfaFouriersLawCacheFiller;
87
88 void updateHeatConduction(const Problem& problem,
89 const Element& element,
90 const FVElementGeometry& fvGeometry,
91 const ElementVolumeVariables& elemVolVars,
92 const SubControlVolumeFace &scvf)
93 {
94 tij_ = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
95 }
96
97 const Scalar& heatConductionTij() const
98 { return tij_; }
99
100 private:
101 Scalar tij_;
102 };
103
104public:
107
109 using Cache = TpfaFouriersLawCache;
110
118 static Scalar flux(const Problem& problem,
119 const Element& element,
120 const FVElementGeometry& fvGeometry,
121 const ElementVolumeVariables& elemVolVars,
122 const SubControlVolumeFace& scvf,
123 const ElementFluxVarsCache& elemFluxVarsCache)
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 = elemVolVars[scvf.insideScvIdx()].temperature();
130 const auto tOutside = scvf.numOutsideScvs() == 1 ? elemVolVars[scvf.outsideScvIdx()].temperature()
131 : branchingFacetTemperature_(problem, element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf, tInside, tij);
132
133 return tij*(tInside - tOutside);
134 }
135
137 static Scalar calculateTransmissibility(const Problem& problem,
138 const Element& element,
139 const FVElementGeometry& fvGeometry,
140 const ElementVolumeVariables& elemVolVars,
141 const SubControlVolumeFace& scvf)
142 {
143 Scalar tij;
144
145 const auto insideScvIdx = scvf.insideScvIdx();
146 const auto& insideScv = fvGeometry.scv(insideScvIdx);
147 const auto& insideVolVars = elemVolVars[insideScvIdx];
148
149 const auto insideLambda = insideVolVars.effectiveThermalConductivity();
150 const Scalar ti = computeTpfaTransmissibility(scvf, insideScv, insideLambda, insideVolVars.extrusionFactor());
151
152 // for the boundary (dirichlet) or at branching points we only need ti
153 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
154 {
155 tij = Extrusion::area(scvf)*ti;
156 }
157 // otherwise we compute a tpfa harmonic mean
158 else
159 {
160 const auto outsideScvIdx = scvf.outsideScvIdx();
161 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
162 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
163
164 const auto outsideLambda = outsideVolVars.effectiveThermalConductivity();
165 Scalar tj;
166 if (dim == dimWorld)
167 // assume the normal vector from outside is anti parallel so we save flipping a vector
168 tj = -1.0*computeTpfaTransmissibility(scvf, outsideScv, outsideLambda, outsideVolVars.extrusionFactor());
169 else
170 tj = computeTpfaTransmissibility(fvGeometry.flipScvf(scvf.index()), outsideScv, outsideLambda, outsideVolVars.extrusionFactor());
171
172 // check for division by zero!
173 if (ti*tj <= 0.0)
174 tij = 0;
175 else
176 tij = Extrusion::area(scvf)*(ti * tj)/(ti + tj);
177 }
178
179 return tij;
180 }
181
182private:
183
185 static Scalar branchingFacetTemperature_(const Problem& problem,
186 const Element& element,
187 const FVElementGeometry& fvGeometry,
188 const ElementVolumeVariables& elemVolVars,
189 const ElementFluxVarsCache& elemFluxVarsCache,
190 const SubControlVolumeFace& scvf,
191 Scalar insideTemperature,
192 Scalar insideTi)
193 {
194 Scalar sumTi(insideTi);
195 Scalar sumTempTi(insideTi*insideTemperature);
196
197 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
198 {
199 const auto outsideScvIdx = scvf.outsideScvIdx(i);
200 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
201 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
202 const auto& outsideFluxVarsCache = elemFluxVarsCache[flippedScvf];
203
204 auto outsideTi = outsideFluxVarsCache.heatConductionTij();
205 sumTi += outsideTi;
206 sumTempTi += outsideTi*outsideVolVars.temperature();
207 }
208 return sumTempTi/sumTi;
209 }
210};
211
212} // end namespace Dumux
213
214#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 (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: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:137
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:118
TpfaFouriersLawCache Cache
export the type for the corresponding cache
Definition: flux/cctpfa/fourierslaw.hh:109
Declares all properties used in Dumux.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...