25#ifndef DUMUX_DISCRETIZATION_BOX_FOURIERS_LAW_NONEQUILIBRIUM_HH
26#define DUMUX_DISCRETIZATION_BOX_FOURIERS_LAW_NONEQUILIBRIUM_HH
28#include <dune/common/fvector.hh>
40template <
class TypeTag, DiscretizationMethod discMethod>
41class FouriersLawNonEquilibriumImplementation;
47template <
class TypeTag>
53 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
59 using Element =
typename GridView::template Codim<0>::Entity;
61 static constexpr auto numEnergyEqSolid = getPropValue<TypeTag, Properties::NumEnergyEqSolid>();
62 static constexpr auto numEnergyEqFluid = getPropValue<TypeTag, Properties::NumEnergyEqFluid>();
63 static constexpr auto numEnergyEq = numEnergyEqSolid + numEnergyEqFluid;
64 static constexpr auto sPhaseIdx = ModelTraits::numFluidPhases();
67 static Scalar
flux(
const Problem& problem,
68 const Element& element,
69 const FVElementGeometry& fvGeometry,
70 const ElementVolumeVariables& elemVolVars,
71 const SubControlVolumeFace& scvf,
73 const ElementFluxVariablesCache& elemFluxVarsCache)
76 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
77 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
78 const auto& insideVolVars = elemVolVars[insideScv];
79 const auto& outsideVolVars = elemVolVars[outsideScv];
80 const auto computeLambda = [&](
const auto& v){
81 if constexpr (numEnergyEq == 1)
82 return v.effectiveThermalConductivity();
83 else if constexpr (numEnergyEqFluid == 1)
84 return (phaseIdx != sPhaseIdx)
85 ? v.effectiveFluidThermalConductivity()
86 : v.effectiveSolidThermalConductivity();
88 return v.effectivePhaseThermalConductivity(phaseIdx);
91 auto insideLambda = computeLambda(insideVolVars);
92 auto outsideLambda = computeLambda(outsideVolVars);
95 insideLambda *= insideVolVars.extrusionFactor();
96 outsideLambda *= outsideVolVars.extrusionFactor();
99 const auto lambda = problem.spatialParams().harmonicMean(insideLambda, outsideLambda, scvf.unitOuterNormal());
102 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
104 Dune::FieldVector<Scalar, GridView::dimensionworld> gradTemp(0.0);
105 for (
auto&& scv : scvs(fvGeometry))
108 if (phaseIdx < numEnergyEqFluid)
109 gradTemp.axpy(elemVolVars[scv].temperatureFluid(phaseIdx), fluxVarsCache.gradN(scv.indexInElement()));
111 gradTemp.axpy(elemVolVars[scv].temperatureSolid(), fluxVarsCache.gradN(scv.indexInElement()));
115 return -1.0*
vtmv(scvf.unitOuterNormal(), lambda, gradTemp)*scvf.area();
Define some often used mathematical functions.
The available discretization methods in Dumux.
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:840
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
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)
Definition: box/fourierslawnonequilibrium.hh:67
Declares all properties used in Dumux.