3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 SubControlVolume = typename FVElementGeometry::SubControlVolume;
53 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
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
63 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
64
66
68 class TpfaFouriersLawCacheFiller
69 {
70 public:
73 template<class FluxVariablesCacheFiller>
74 static void fill(FluxVariablesCache& scvfFluxVarsCache,
75 const Problem& problem,
76 const Element& element,
77 const FVElementGeometry& fvGeometry,
78 const ElementVolumeVariables& elemVolVars,
79 const SubControlVolumeFace& scvf,
80 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
81 {
82 scvfFluxVarsCache.updateHeatConduction(problem, element, fvGeometry, elemVolVars, scvf);
83 }
84 };
85
87 class TpfaFouriersLawCache
88 {
89 public:
90 using Filler = TpfaFouriersLawCacheFiller;
91
92 void updateHeatConduction(const Problem& problem,
93 const Element& element,
94 const FVElementGeometry& fvGeometry,
95 const ElementVolumeVariables& elemVolVars,
96 const SubControlVolumeFace &scvf)
97 {
98 tij_ = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
99 }
100
101 const Scalar& heatConductionTij() const
102 { return tij_; }
103
104 private:
105 Scalar tij_;
106 };
107
108public:
111
113 using Cache = TpfaFouriersLawCache;
114
116 static Scalar flux(const Problem& problem,
117 const Element& element,
118 const FVElementGeometry& fvGeometry,
119 const ElementVolumeVariables& elemVolVars,
120 const SubControlVolumeFace& scvf,
121 const ElementFluxVarsCache& elemFluxVarsCache)
122 {
123 // heat conductivities are always solution dependent (?)
124 Scalar tij = elemFluxVarsCache[scvf].heatConductionTij();
125
126 // get the inside/outside temperatures
127 const auto tInside = elemVolVars[scvf.insideScvIdx()].temperature();
128 const auto tOutside = scvf.numOutsideScvs() == 1 ? elemVolVars[scvf.outsideScvIdx()].temperature()
129 : branchingFacetTemperature_(problem, element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf, tInside, tij);
130
131 return tij*(tInside - tOutside);
132 }
133
135 static Scalar calculateTransmissibility(const Problem& problem,
136 const Element& element,
137 const FVElementGeometry& fvGeometry,
138 const ElementVolumeVariables& elemVolVars,
139 const SubControlVolumeFace& scvf)
140 {
141 Scalar tij;
142
143 const auto insideScvIdx = scvf.insideScvIdx();
144 const auto& insideScv = fvGeometry.scv(insideScvIdx);
145 const auto& insideVolVars = elemVolVars[insideScvIdx];
146
147 const auto insideLambda = Deprecated::template effectiveThermalConductivity<ThermalConductivityModel>(
148 insideVolVars, problem.spatialParams(), element, fvGeometry, insideScv);
149 const Scalar ti = computeTpfaTransmissibility(scvf, insideScv, insideLambda, insideVolVars.extrusionFactor());
150
151 // for the boundary (dirichlet) or at branching points we only need ti
152 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
153 {
154 tij = scvf.area()*ti;
155 }
156 // otherwise we compute a tpfa harmonic mean
157 else
158 {
159 const auto outsideScvIdx = scvf.outsideScvIdx();
160 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
161 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
162 const auto outsideElement = fvGeometry.gridGeometry().element(outsideScvIdx);
163
164 const auto outsideLambda = Deprecated::template effectiveThermalConductivity<ThermalConductivityModel>(
165 outsideVolVars, problem.spatialParams(), outsideElement, fvGeometry, outsideScv);
166 Scalar tj;
167 if (dim == dimWorld)
168 // assume the normal vector from outside is anti parallel so we save flipping a vector
169 tj = -1.0*computeTpfaTransmissibility(scvf, outsideScv, outsideLambda, outsideVolVars.extrusionFactor());
170 else
171 tj = computeTpfaTransmissibility(fvGeometry.flipScvf(scvf.index()), outsideScv, outsideLambda, outsideVolVars.extrusionFactor());
172
173 // check for division by zero!
174 if (ti*tj <= 0.0)
175 tij = 0;
176 else
177 tij = scvf.area()*(ti * tj)/(ti + tj);
178 }
179
180 return tij;
181 }
182
183private:
184
186 static Scalar branchingFacetTemperature_(const Problem& problem,
187 const Element& element,
188 const FVElementGeometry& fvGeometry,
189 const ElementVolumeVariables& elemVolVars,
190 const ElementFluxVarsCache& elemFluxVarsCache,
191 const SubControlVolumeFace& scvf,
192 Scalar insideTemperature,
193 Scalar insideTi)
194 {
195 Scalar sumTi(insideTi);
196 Scalar sumTempTi(insideTi*insideTemperature);
197
198 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
199 {
200 const auto outsideScvIdx = scvf.outsideScvIdx(i);
201 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
202 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
203 const auto& outsideFluxVarsCache = elemFluxVarsCache[flippedScvf];
204
205 auto outsideTi = outsideFluxVarsCache.heatConductionTij();
206 sumTi += outsideTi;
207 sumTempTi += outsideTi*outsideVolVars.temperature();
208 }
209 return sumTempTi/sumTi;
210 }
211};
212
213} // end namespace Dumux
214
215#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
make the local view function available whenever we use the grid geometry
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: fourierslaw.hh:37
Fourier's law for cell-centered finite volume schemes with two-point flux approximation.
Definition: 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: cctpfa/fourierslaw.hh:135
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: cctpfa/fourierslaw.hh:116
TpfaFouriersLawCache Cache
export the type for the corresponding cache
Definition: cctpfa/fourierslaw.hh:113
Declares all properties used in Dumux.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...