3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_FLUX_PNM_FOURIERS_LAW_HH
26#define DUMUX_FLUX_PNM_FOURIERS_LAW_HH
27
28#include <dumux/common/math.hh>
30#include <type_traits>
31
32namespace Dumux::PoreNetwork {
33
34namespace Detail {
35
37
38} // end namespace Detail
39
44template<class MolecularDiffusionType = Detail::NoDiffusionType>
46{
47
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)
56 {
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];
62
63 static constexpr auto numPhases = ElementVolumeVariables::VolumeVariables::numFluidPhases();
64 using Scalar = typename ElementVolumeVariables::VolumeVariables::PrimaryVariables::value_type;
65 Scalar heatflux = 0;
66
67 for (int phaseIdx = 0; phaseIdx < numPhases; phaseIdx++)
68 {
69 auto insideThermalConducitivity = insideVolVars.fluidThermalConductivity(phaseIdx);
70 auto outsideThermalConducitivity = outsideVolVars.fluidThermalConductivity(phaseIdx);
71
72 auto thermalConductivity = Dumux::harmonicMean(insideThermalConducitivity, outsideThermalConducitivity);
73 auto area = fluxVarsCache.throatCrossSectionalArea(phaseIdx);
74
75 // calculate the temperature gradient
76 const Scalar deltaT = insideVolVars.temperature() - outsideVolVars.temperature();
77 const Scalar gradT = deltaT/fluxVarsCache.throatLength();
78
79 heatflux += thermalConductivity*gradT*area;
80
81 if constexpr (!std::is_same_v<MolecularDiffusionType, Detail::NoDiffusionType>)
82 heatflux += componentEnthalpyFlux_(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache, phaseIdx);
83 }
84
85 return heatflux;
86 }
87
88private:
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,
97 const int phaseIdx)
98 {
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)
104 {
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);
108
109 if (MolecularDiffusionType::referenceSystemFormulation() == ReferenceSystemFormulation::massAveraged)
110 heatflux += diffusiveFlux[compIdx] * componentEnthalpy;
111 else
112 heatflux += diffusiveFlux[compIdx] * FluidSystem::molarMass(compIdx) * componentEnthalpy;
113 }
114
115 return heatflux;
116 }
117};
118
119} // end namespace Dumux::PoreNetwork
120
121#endif
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