25#ifndef DUMUX_FLUX_PNM_FICKS_LAW_HH
26#define DUMUX_FLUX_PNM_FICKS_LAW_HH
28#include <dune/common/fvector.hh>
39template<
class Scalar,
int numPhases,
int numComponents,
48 {
return referenceSystem; }
50 template<
class Problem,
class Element,
class FVElementGeometry,
51 class ElementVolumeVariables,
class ElementFluxVariablesCache>
52 static Dune::FieldVector<Scalar, numComponents>
53 flux(
const Problem& problem,
54 const Element& element,
55 const FVElementGeometry& fvGeometry,
56 const ElementVolumeVariables& elemVolVars,
57 const typename FVElementGeometry::SubControlVolumeFace& scvf,
59 const ElementFluxVariablesCache& elemFluxVarsCache)
61 Dune::FieldVector<Scalar, numComponents> componentFlux(0.0);
64 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
65 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
67 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
70 const Scalar throatLength = fluxVarsCache.throatLength();
71 const Scalar phaseCrossSectionalArea = fluxVarsCache.throatCrossSectionalArea(phaseIdx);
73 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
75 if(compIdx == phaseIdx)
78 auto insideD = getDiffusionCoefficient_(phaseIdx, compIdx, insideVolVars);
79 auto outsideD = getDiffusionCoefficient_(phaseIdx, compIdx, outsideVolVars);
82 insideD *= insideVolVars.extrusionFactor();
83 outsideD *= outsideVolVars.extrusionFactor();
88 const Scalar insideMoleFraction =
massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
89 const Scalar outsideMoleFraction =
massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
91 componentFlux[compIdx] =
density * (insideMoleFraction - outsideMoleFraction) / throatLength * D * phaseCrossSectionalArea;
92 componentFlux[phaseIdx] -= componentFlux[compIdx];
98 template<
class VolumeVariables>
99 static Scalar getDiffusionCoefficient_(
const int phaseIdx,
const int compIdx,
100 const VolumeVariables& volVars)
102 using FluidSystem =
typename VolumeVariables::FluidSystem;
104 if constexpr (!FluidSystem::isTracerFluidSystem())
106 const auto mainCompIdx = FluidSystem::getMainComponent(phaseIdx);
107 return volVars.diffusionCoefficient(phaseIdx, mainCompIdx, compIdx);
110 return volVars.diffusionCoefficient(0, 0, compIdx);
Define some often used mathematical functions.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
Container storing the diffusion coefficients required by Fick's law. Uses the minimal possible contai...
VolumeVariables::PrimaryVariables::value_type massOrMoleFraction(const VolumeVariables &volVars, ReferenceSystemFormulation referenceSys, const int phaseIdx, const int compIdx)
returns the mass or mole fraction to be used in Fick's law based on the reference system
Definition: referencesystemformulation.hh:66
VolumeVariables::PrimaryVariables::value_type massOrMolarDensity(const VolumeVariables &volVars, ReferenceSystemFormulation referenceSys, const int phaseIdx)
evaluates the density to be used in Fick's law based on the reference system
Definition: referencesystemformulation.hh:55
ReferenceSystemFormulation
The formulations available for Fick's law related to the reference system.
Definition: referencesystemformulation.hh:45
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:33
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
Container storing the diffusion coefficients required by Fick's law. Uses the minimal possible contai...
Definition: fickiandiffusioncoefficients.hh:44
Specialization of Fick's Law for the pore-network model.
Definition: flux/porenetwork/fickslaw.hh:42
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: flux/porenetwork/fickslaw.hh:47
static Dune::FieldVector< Scalar, numComponents > flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const typename FVElementGeometry::SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVariablesCache &elemFluxVarsCache)
Definition: flux/porenetwork/fickslaw.hh:53