version 3.8
flux/porenetwork/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//
13#ifndef DUMUX_FLUX_PNM_FOURIERS_LAW_HH
14#define DUMUX_FLUX_PNM_FOURIERS_LAW_HH
15
16#include <dumux/common/math.hh>
18#include <type_traits>
19
20namespace Dumux::PoreNetwork {
21
22namespace Detail {
23
25
26} // end namespace Detail
27
32template<class MolecularDiffusionType = Detail::NoDiffusionType>
34{
35
36 template<class Problem, class Element, class FVElementGeometry,
37 class ElementVolumeVariables, class ElementFluxVariablesCache>
38 static auto flux(const Problem& problem,
39 const Element& element,
40 const FVElementGeometry& fvGeometry,
41 const ElementVolumeVariables& elemVolVars,
42 const typename FVElementGeometry::SubControlVolumeFace& scvf,
43 const ElementFluxVariablesCache& elemFluxVarsCache)
44 {
45 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
46 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
47 const auto& insideVolVars = elemVolVars[insideScv];
48 const auto& outsideVolVars = elemVolVars[outsideScv];
49 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
50
51 static constexpr auto numPhases = ElementVolumeVariables::VolumeVariables::numFluidPhases();
52 using Scalar = typename ElementVolumeVariables::VolumeVariables::PrimaryVariables::value_type;
53 Scalar heatflux = 0;
54
55 for (int phaseIdx = 0; phaseIdx < numPhases; phaseIdx++)
56 {
57 auto insideThermalConducitivity = insideVolVars.fluidThermalConductivity(phaseIdx);
58 auto outsideThermalConducitivity = outsideVolVars.fluidThermalConductivity(phaseIdx);
59
60 auto thermalConductivity = Dumux::harmonicMean(insideThermalConducitivity, outsideThermalConducitivity);
61 auto area = fluxVarsCache.throatCrossSectionalArea(phaseIdx);
62
63 // calculate the temperature gradient
64 const Scalar deltaT = insideVolVars.temperature() - outsideVolVars.temperature();
65 const Scalar gradT = deltaT/fluxVarsCache.throatLength();
66
67 heatflux += thermalConductivity*gradT*area;
68
69 if constexpr (!std::is_same_v<MolecularDiffusionType, Detail::NoDiffusionType>)
70 heatflux += componentEnthalpyFlux_(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache, phaseIdx);
71 }
72
73 return heatflux;
74 }
75
76private:
77 template<class Problem, class Element, class FVElementGeometry,
78 class ElementVolumeVariables, class ElementFluxVariablesCache>
79 static auto componentEnthalpyFlux_(const Problem& problem,
80 const Element& element,
81 const FVElementGeometry& fvGeometry,
82 const ElementVolumeVariables& elemVolVars,
83 const typename FVElementGeometry::SubControlVolumeFace& scvf,
84 const ElementFluxVariablesCache& elemFluxVarsCache,
85 const int phaseIdx)
86 {
87 using Scalar = typename ElementVolumeVariables::VolumeVariables::PrimaryVariables::value_type;
88 Scalar heatflux = 0.0;
89 using FluidSystem = typename ElementVolumeVariables::VolumeVariables::FluidSystem;
90 const auto diffusiveFlux = MolecularDiffusionType::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
91 for (int compIdx = 0; compIdx < ElementVolumeVariables::VolumeVariables::numFluidComponents(); ++compIdx)
92 {
93 const bool insideIsUpstream = diffusiveFlux[compIdx] > 0.0;
94 const auto& upstreamVolVars = insideIsUpstream ? elemVolVars[scvf.insideScvIdx()] : elemVolVars[scvf.outsideScvIdx()];
95 const Scalar componentEnthalpy = FluidSystem::componentEnthalpy(upstreamVolVars.fluidState(), phaseIdx, compIdx);
96
97 if (MolecularDiffusionType::referenceSystemFormulation() == ReferenceSystemFormulation::massAveraged)
98 heatflux += diffusiveFlux[compIdx] * componentEnthalpy;
99 else
100 heatflux += diffusiveFlux[compIdx] * FluidSystem::molarMass(compIdx) * componentEnthalpy;
101 }
102
103 return heatflux;
104 }
105};
106
107} // end namespace Dumux::PoreNetwork
108
109#endif
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
Define some often used mathematical functions.
Definition: discretization/porenetwork/fvelementgeometry.hh:24
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
Definition: flux/porenetwork/fourierslaw.hh:24
Specialization of Fourier's Law for the pore-network model.
Definition: flux/porenetwork/fourierslaw.hh:34
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:38