3.2-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
32
33#include <dumux/common/deprecated.hh> // effective thermal conductivity interface
34
35namespace Dumux {
36
37// forward declaration
38template<class TypeTag, DiscretizationMethod discMethod>
39class FouriersLawImplementation;
40
45template <class TypeTag>
47{
51 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
52 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
54 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
55 using Element = typename GridView::template Codim<0>::Entity;
58
59 static const int dim = GridView::dimension;
60 static const int dimWorld = GridView::dimensionworld;
61
63 class TpfaFouriersLawCacheFiller
64 {
65 public:
68 template<class FluxVariablesCacheFiller>
69 static void fill(FluxVariablesCache& scvfFluxVarsCache,
70 const Problem& problem,
71 const Element& element,
72 const FVElementGeometry& fvGeometry,
73 const ElementVolumeVariables& elemVolVars,
74 const SubControlVolumeFace& scvf,
75 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
76 {
77 scvfFluxVarsCache.updateHeatConduction(problem, element, fvGeometry, elemVolVars, scvf);
78 }
79 };
80
82 class TpfaFouriersLawCache
83 {
84 public:
85 using Filler = TpfaFouriersLawCacheFiller;
86
87 void updateHeatConduction(const Problem& problem,
88 const Element& element,
89 const FVElementGeometry& fvGeometry,
90 const ElementVolumeVariables& elemVolVars,
91 const SubControlVolumeFace &scvf)
92 {
93 tij_ = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
94 }
95
96 const Scalar& heatConductionTij() const
97 { return tij_; }
98
99 private:
100 Scalar tij_;
101 };
102
103public:
106
108 using Cache = TpfaFouriersLawCache;
109
111 static Scalar flux(const Problem& problem,
112 const Element& element,
113 const FVElementGeometry& fvGeometry,
114 const ElementVolumeVariables& elemVolVars,
115 const SubControlVolumeFace& scvf,
116 const ElementFluxVarsCache& elemFluxVarsCache)
117 {
118 // heat conductivities are always solution dependent (?)
119 Scalar tij = elemFluxVarsCache[scvf].heatConductionTij();
120
121 // get the inside/outside temperatures
122 const auto tInside = elemVolVars[scvf.insideScvIdx()].temperature();
123 const auto tOutside = scvf.numOutsideScvs() == 1 ? elemVolVars[scvf.outsideScvIdx()].temperature()
124 : branchingFacetTemperature_(problem, element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf, tInside, tij);
125
126 return tij*(tInside - tOutside);
127 }
128
130 static Scalar calculateTransmissibility(const Problem& problem,
131 const Element& element,
132 const FVElementGeometry& fvGeometry,
133 const ElementVolumeVariables& elemVolVars,
134 const SubControlVolumeFace& scvf)
135 {
136 Scalar tij;
137
138 const auto insideScvIdx = scvf.insideScvIdx();
139 const auto& insideScv = fvGeometry.scv(insideScvIdx);
140 const auto& insideVolVars = elemVolVars[insideScvIdx];
141
142 const auto insideLambda = insideVolVars.effectiveThermalConductivity();
143 const Scalar ti = computeTpfaTransmissibility(scvf, insideScv, insideLambda, insideVolVars.extrusionFactor());
144
145 // for the boundary (dirichlet) or at branching points we only need ti
146 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
147 {
148 tij = scvf.area()*ti;
149 }
150 // otherwise we compute a tpfa harmonic mean
151 else
152 {
153 const auto outsideScvIdx = scvf.outsideScvIdx();
154 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
155 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
156
157 const auto outsideLambda = outsideVolVars.effectiveThermalConductivity();
158 Scalar tj;
159 if (dim == dimWorld)
160 // assume the normal vector from outside is anti parallel so we save flipping a vector
161 tj = -1.0*computeTpfaTransmissibility(scvf, outsideScv, outsideLambda, outsideVolVars.extrusionFactor());
162 else
163 tj = computeTpfaTransmissibility(fvGeometry.flipScvf(scvf.index()), outsideScv, outsideLambda, outsideVolVars.extrusionFactor());
164
165 // check for division by zero!
166 if (ti*tj <= 0.0)
167 tij = 0;
168 else
169 tij = scvf.area()*(ti * tj)/(ti + tj);
170 }
171
172 return tij;
173 }
174
175private:
176
178 static Scalar branchingFacetTemperature_(const Problem& problem,
179 const Element& element,
180 const FVElementGeometry& fvGeometry,
181 const ElementVolumeVariables& elemVolVars,
182 const ElementFluxVarsCache& elemFluxVarsCache,
183 const SubControlVolumeFace& scvf,
184 Scalar insideTemperature,
185 Scalar insideTi)
186 {
187 Scalar sumTi(insideTi);
188 Scalar sumTempTi(insideTi*insideTemperature);
189
190 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
191 {
192 const auto outsideScvIdx = scvf.outsideScvIdx(i);
193 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
194 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
195 const auto& outsideFluxVarsCache = elemFluxVarsCache[flippedScvf];
196
197 auto outsideTi = outsideFluxVarsCache.heatConductionTij();
198 sumTi += outsideTi;
199 sumTempTi += outsideTi*outsideVolVars.temperature();
200 }
201 return sumTempTi/sumTi;
202 }
203};
204
205} // end namespace Dumux
206
207#endif
Helpers for deprecation.
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
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
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:130
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVarsCache &elemFluxVarsCache)
Compute the heat condution flux assuming thermal equilibrium.
Definition: flux/cctpfa/fourierslaw.hh:111
TpfaFouriersLawCache Cache
export the type for the corresponding cache
Definition: flux/cctpfa/fourierslaw.hh:108
Declares all properties used in Dumux.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...