3.5-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
33
34namespace Dumux {
35
36// forward declaration
37template<class TypeTag, class DiscretizationMethod>
38class FouriersLawImplementation;
39
44template <class TypeTag>
45class FouriersLawImplementation<TypeTag, DiscretizationMethods::Staggered>
46{
49 using FVElementGeometry = typename GridGeometry::LocalView;
50 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
51 using Extrusion = Extrusion_t<GridGeometry>;
52 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
53 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
55
56public:
58 // state the discretization method this implementation belongs to
59 static constexpr DiscretizationMethod discMethod{};
60
64
72 template<class Problem>
73 static Scalar flux(const Problem& problem,
74 const Element& element,
75 const FVElementGeometry& fvGeometry,
76 const ElementVolumeVariables& elemVolVars,
77 const SubControlVolumeFace &scvf)
78 {
79 Scalar flux(0.0);
80
81 // conductive energy flux is zero for outflow boundary conditions
82 if (scvf.boundary() && problem.boundaryTypes(element, scvf).isOutflow(Indices::energyEqIdx))
83 return flux;
84
85 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
86 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
87 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
88
89 const Scalar insideTemperature = insideVolVars.temperature();
90 const Scalar outsideTemperature = outsideVolVars.temperature();
91
92 const Scalar insideLambda = insideVolVars.effectiveThermalConductivity() * insideVolVars.extrusionFactor();
93 const Scalar insideDistance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
94
95 if (scvf.boundary())
96 {
97 flux = insideLambda * (insideTemperature - outsideTemperature) / insideDistance;
98 }
99 else
100 {
101 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
102 const Scalar outsideLambda = outsideVolVars.effectiveThermalConductivity() * outsideVolVars.extrusionFactor();
103 const Scalar outsideDistance = (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm();
104 const Scalar avgLambda = harmonicMean(insideLambda, outsideLambda, insideDistance, outsideDistance);
105
106 flux = avgLambda * (insideTemperature - outsideTemperature) / (insideDistance + outsideDistance);
107 }
108
109 flux *= Extrusion::area(scvf);
110 return flux;
111 }
112};
113
114} // end namespace Dumux
115
116#endif
Define some often used mathematical functions.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
Classes related to flux variables caching.
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:69
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
forward declaration of the method-specific implementation
Definition: flux/box/fourierslaw.hh:38
Definition: fluxvariablescaching.hh:67
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Returns the heat flux within the porous medium (in J/s) across the given sub-control volume face.
Definition: flux/staggered/freeflow/fourierslaw.hh:73
Declares all properties used in Dumux.