version 3.8
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_DISCRETIZATION_STAGGERED_FOURIERS_LAW_HH
13#define DUMUX_DISCRETIZATION_STAGGERED_FOURIERS_LAW_HH
14
16#include <dumux/common/math.hh>
17
21
22namespace Dumux {
23
24// forward declaration
25template<class TypeTag, class DiscretizationMethod>
26class FouriersLawImplementation;
27
32template <class TypeTag>
33class FouriersLawImplementation<TypeTag, DiscretizationMethods::Staggered>
34{
37 using FVElementGeometry = typename GridGeometry::LocalView;
38 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
39 using Extrusion = Extrusion_t<GridGeometry>;
40 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
41 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
43
44public:
46 // state the discretization method this implementation belongs to
47 static constexpr DiscretizationMethod discMethod{};
48
52
60 template<class Problem>
61 static Scalar flux(const Problem& problem,
62 const Element& element,
63 const FVElementGeometry& fvGeometry,
64 const ElementVolumeVariables& elemVolVars,
65 const SubControlVolumeFace &scvf)
66 {
67 Scalar flux(0.0);
68
69 // conductive energy flux is zero for outflow boundary conditions
70 if (scvf.boundary() && problem.boundaryTypes(element, scvf).isOutflow(Indices::energyEqIdx))
71 return flux;
72
73 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
74 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
75 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
76
77 const Scalar insideTemperature = insideVolVars.temperature();
78 const Scalar outsideTemperature = outsideVolVars.temperature();
79
80 const Scalar insideLambda = insideVolVars.effectiveThermalConductivity() * insideVolVars.extrusionFactor();
81 const Scalar insideDistance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
82
83 if (scvf.boundary())
84 {
85 flux = insideLambda * (insideTemperature - outsideTemperature) / insideDistance;
86 }
87 else
88 {
89 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
90 const Scalar outsideLambda = outsideVolVars.effectiveThermalConductivity() * outsideVolVars.extrusionFactor();
91 const Scalar outsideDistance = (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm();
92 const Scalar avgLambda = harmonicMean(insideLambda, outsideLambda, insideDistance, outsideDistance);
93
94 flux = avgLambda * (insideTemperature - outsideTemperature) / (insideDistance + outsideDistance);
95 }
96
97 flux *= Extrusion::area(fvGeometry, scvf);
98 return flux;
99 }
100};
101
102} // end namespace Dumux
103
104#endif
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:61
forward declaration of the method-specific implementation
Definition: flux/box/fourierslaw.hh:26
Defines all properties used in Dumux.
Helper classes to compute the integration elements.
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:57
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Define some often used mathematical functions.
The available discretization methods in Dumux.
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
Definition: fluxvariablescaching.hh:55