24#ifndef DUMUX_DISCRETIZATION_STAGGERED_FOURIERS_LAW_HH
25#define DUMUX_DISCRETIZATION_STAGGERED_FOURIERS_LAW_HH
36template<
class TypeTag, DiscretizationMethod discMethod>
37class FouriersLawImplementation;
43template <
class TypeTag>
48 using FVElementGeometry =
typename GridGeometry::LocalView;
49 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
51 using Element =
typename GridGeometry::GridView::template Codim<0>::Entity;
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)
73 if (scvf.boundary() && problem.boundaryTypes(element, scvf).isOutflow(Indices::energyEqIdx))
76 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
77 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
78 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
80 const Scalar insideTemperature = insideVolVars.temperature();
81 const Scalar outsideTemperature = outsideVolVars.temperature();
83 const Scalar insideLambda = insideVolVars.effectiveThermalConductivity() * insideVolVars.extrusionFactor();
84 const Scalar insideDistance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
88 flux = insideLambda * (insideTemperature - outsideTemperature) / insideDistance;
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);
97 flux = avgLambda * (insideTemperature - outsideTemperature) / (insideDistance + outsideDistance);
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
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.