25#ifndef DUMUX_DISCRETIZATION_BOX_FOURIERS_LAW_NONEQUILIBRIUM_HH
26#define DUMUX_DISCRETIZATION_BOX_FOURIERS_LAW_NONEQUILIBRIUM_HH
28#include <dune/common/fvector.hh>
39template <
class TypeTag, DiscretizationMethod discMethod>
40class FouriersLawNonEquilibriumImplementation;
46template <
class TypeTag>
52 using FVElementGeometry =
typename GridGeometry::LocalView;
53 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
60 using Element =
typename GridView::template Codim<0>::Entity;
62 static constexpr auto numEnergyEqSolid = getPropValue<TypeTag, Properties::NumEnergyEqSolid>();
63 static constexpr auto numEnergyEqFluid = getPropValue<TypeTag, Properties::NumEnergyEqFluid>();
64 static constexpr auto numEnergyEq = numEnergyEqSolid + numEnergyEqFluid;
65 static constexpr auto sPhaseIdx = ModelTraits::numFluidPhases();
72 static Scalar
flux(
const Problem& problem,
73 const Element& element,
74 const FVElementGeometry& fvGeometry,
75 const ElementVolumeVariables& elemVolVars,
76 const SubControlVolumeFace& scvf,
78 const ElementFluxVariablesCache& elemFluxVarsCache)
81 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
82 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
83 const auto& insideVolVars = elemVolVars[insideScv];
84 const auto& outsideVolVars = elemVolVars[outsideScv];
85 const auto computeLambda = [&](
const auto& v){
86 if constexpr (numEnergyEq == 1)
87 return v.effectiveThermalConductivity();
88 else if constexpr (numEnergyEqFluid == 1)
89 return (phaseIdx != sPhaseIdx)
90 ? v.effectiveFluidThermalConductivity()
91 : v.effectiveSolidThermalConductivity();
93 return v.effectivePhaseThermalConductivity(phaseIdx);
96 auto insideLambda = computeLambda(insideVolVars);
97 auto outsideLambda = computeLambda(outsideVolVars);
100 insideLambda *= insideVolVars.extrusionFactor();
101 outsideLambda *= outsideVolVars.extrusionFactor();
104 const auto lambda = problem.spatialParams().harmonicMean(insideLambda, outsideLambda, scvf.unitOuterNormal());
107 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
109 Dune::FieldVector<Scalar, GridView::dimensionworld> gradTemp(0.0);
110 for (
auto&& scv : scvs(fvGeometry))
113 if (phaseIdx < numEnergyEqFluid)
114 gradTemp.axpy(elemVolVars[scv].temperatureFluid(phaseIdx), fluxVarsCache.gradN(scv.indexInElement()));
116 gradTemp.axpy(elemVolVars[scv].temperatureSolid(), fluxVarsCache.gradN(scv.indexInElement()));
120 return -1.0*
vtmv(scvf.unitOuterNormal(), lambda, gradTemp)*Extrusion::area(scvf);
Define some often used mathematical functions.
The available discretization methods in Dumux.
Helper classes to compute the integration elements.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
Dune::DenseMatrix< MAT >::value_type vtmv(const Dune::DenseVector< V1 > &v1, const Dune::DenseMatrix< MAT > &M, const Dune::DenseVector< V2 > &v2)
Evaluates the scalar product of a vector v2, projected by a matrix M, with a vector v1.
Definition: math.hh:863
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
Definition: fourierslawnonequilibrium.hh:36
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVariablesCache &elemFluxVarsCache)
Returns the heat flux within a fluid or solid phase (in J/s) across the given sub-control volume face...
Definition: box/fourierslawnonequilibrium.hh:72
Declares all properties used in Dumux.