3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
flux/staggered/freeflow/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_STAGGERED_FOURIERS_LAW_HH
25#define DUMUX_DISCRETIZATION_STAGGERED_FOURIERS_LAW_HH
26
28#include <dumux/common/math.hh>
29
32
33namespace Dumux {
34
35// forward declaration
36template<class TypeTag, DiscretizationMethod discMethod>
37class FouriersLawImplementation;
38
43template <class TypeTag>
45{
48 using FVElementGeometry = typename GridGeometry::LocalView;
49 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
50 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
51 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
53
54public:
55 // state the discretization method this implementation belongs to
57
61
63 template<class Problem>
64 static Scalar flux(const Problem& problem,
65 const Element& element,
66 const FVElementGeometry& fvGeometry,
67 const ElementVolumeVariables& elemVolVars,
68 const SubControlVolumeFace &scvf)
69 {
70 Scalar flux(0.0);
71
72 // conductive energy flux is zero for outflow boundary conditions
73 if (scvf.boundary() && problem.boundaryTypes(element, scvf).isOutflow(Indices::energyEqIdx))
74 return flux;
75
76 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
77 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
78 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
79
80 const Scalar insideTemperature = insideVolVars.temperature();
81 const Scalar outsideTemperature = outsideVolVars.temperature();
82
83 const Scalar insideLambda = insideVolVars.effectiveThermalConductivity() * insideVolVars.extrusionFactor();
84 const Scalar insideDistance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
85
86 if (scvf.boundary())
87 {
88 flux = insideLambda * (insideTemperature - outsideTemperature) / insideDistance;
89 }
90 else
91 {
92 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
93 const Scalar outsideLambda = outsideVolVars.effectiveThermalConductivity() * outsideVolVars.extrusionFactor();
94 const Scalar outsideDistance = (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm();
95 const Scalar avgLambda = harmonicMean(insideLambda, outsideLambda, insideDistance, outsideDistance);
96
97 flux = avgLambda * (insideTemperature - outsideTemperature) / (insideDistance + outsideDistance);
98 }
99
100 flux *= scvf.area();
101 return flux;
102 }
103
104
105};
106} // end namespace
107
108#endif
Define some often used mathematical functions.
The available discretization methods in Dumux.
Classes related to flux variables caching.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
constexpr Scalar harmonicMean(Scalar x, Scalar y, Scalar wx=1.0, Scalar wy=1.0) noexcept
Calculate the (weighted) harmonic mean of two scalar values.
Definition: math.hh:68
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
Definition: fluxvariablescaching.hh:64
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
calculate the diffusive energy fluxes
Definition: flux/staggered/freeflow/fourierslaw.hh:64
Declares all properties used in Dumux.