25#ifndef DUMUX_FLUX_PNM_FOURIERS_LAW_HH
26#define DUMUX_FLUX_PNM_FOURIERS_LAW_HH
44template<
class MolecularDiffusionType = Detail::NoDiffusionType>
48 template<
class Problem,
class Element,
class FVElementGeometry,
49 class ElementVolumeVariables,
class ElementFluxVariablesCache>
50 static auto flux(
const Problem& problem,
51 const Element& element,
52 const FVElementGeometry& fvGeometry,
53 const ElementVolumeVariables& elemVolVars,
54 const typename FVElementGeometry::SubControlVolumeFace& scvf,
55 const ElementFluxVariablesCache& elemFluxVarsCache)
57 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
58 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
59 const auto& insideVolVars = elemVolVars[insideScv];
60 const auto& outsideVolVars = elemVolVars[outsideScv];
61 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
63 static constexpr auto numPhases = ElementVolumeVariables::VolumeVariables::numFluidPhases();
64 using Scalar =
typename ElementVolumeVariables::VolumeVariables::PrimaryVariables::value_type;
67 for (
int phaseIdx = 0; phaseIdx < numPhases; phaseIdx++)
69 auto insideThermalConducitivity = insideVolVars.fluidThermalConductivity(phaseIdx);
70 auto outsideThermalConducitivity = outsideVolVars.fluidThermalConductivity(phaseIdx);
72 auto thermalConductivity =
Dumux::harmonicMean(insideThermalConducitivity, outsideThermalConducitivity);
73 auto area = fluxVarsCache.throatCrossSectionalArea(phaseIdx);
76 const Scalar deltaT = insideVolVars.temperature() - outsideVolVars.temperature();
77 const Scalar gradT = deltaT/fluxVarsCache.throatLength();
79 heatflux += thermalConductivity*gradT*area;
81 if constexpr (!std::is_same_v<MolecularDiffusionType, Detail::NoDiffusionType>)
82 heatflux += componentEnthalpyFlux_(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache, phaseIdx);
89 template<
class Problem,
class Element,
class FVElementGeometry,
90 class ElementVolumeVariables,
class ElementFluxVariablesCache>
91 static auto componentEnthalpyFlux_(
const Problem& problem,
92 const Element& element,
93 const FVElementGeometry& fvGeometry,
94 const ElementVolumeVariables& elemVolVars,
95 const typename FVElementGeometry::SubControlVolumeFace& scvf,
96 const ElementFluxVariablesCache& elemFluxVarsCache,
99 using Scalar =
typename ElementVolumeVariables::VolumeVariables::PrimaryVariables::value_type;
100 Scalar heatflux = 0.0;
101 using FluidSystem =
typename ElementVolumeVariables::VolumeVariables::FluidSystem;
102 const auto diffusiveFlux = MolecularDiffusionType::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
103 for (
int compIdx = 0; compIdx < ElementVolumeVariables::VolumeVariables::numFluidComponents(); ++compIdx)
105 const bool insideIsUpstream = diffusiveFlux[compIdx] > 0.0;
106 const auto& upstreamVolVars = insideIsUpstream ? elemVolVars[scvf.insideScvIdx()] : elemVolVars[scvf.outsideScvIdx()];
107 const Scalar componentEnthalpy = FluidSystem::componentEnthalpy(upstreamVolVars.fluidState(), phaseIdx, compIdx);
110 heatflux += diffusiveFlux[compIdx] * componentEnthalpy;
112 heatflux += diffusiveFlux[compIdx] * FluidSystem::molarMass(compIdx) * componentEnthalpy;
Define some often used mathematical functions.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
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: discretization/porenetwork/fvelementgeometry.hh:34
Definition: flux/porenetwork/fourierslaw.hh:36
Specialization of Fourier's Law for the pore-network model.
Definition: flux/porenetwork/fourierslaw.hh:46
static auto flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElementFluxVariablesCache &elemFluxVarsCache)
Definition: flux/porenetwork/fourierslaw.hh:50