version 3.8
flux/porenetwork/fickslaw.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_FICKS_LAW_HH
14#define DUMUX_FLUX_PNM_FICKS_LAW_HH
15
16#include <dune/common/fvector.hh>
17#include <dumux/common/math.hh>
20
21namespace Dumux::PoreNetwork {
22
27template<class Scalar, int numPhases, int numComponents,
30{
31public:
33
34 //return the reference system
36 { return referenceSystem; }
37
38 template<class Problem, class Element, class FVElementGeometry,
39 class ElementVolumeVariables, class ElementFluxVariablesCache>
40 static Dune::FieldVector<Scalar, numComponents>
41 flux(const Problem& problem,
42 const Element& element,
43 const FVElementGeometry& fvGeometry,
44 const ElementVolumeVariables& elemVolVars,
45 const typename FVElementGeometry::SubControlVolumeFace& scvf,
46 const int phaseIdx,
47 const ElementFluxVariablesCache& elemFluxVarsCache)
48 {
49 Dune::FieldVector<Scalar, numComponents> componentFlux(0.0);
50
51 // get inside and outside diffusion tensors and calculate the harmonic mean
52 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
53 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
54
55 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
56
57 const Scalar density = 0.5 * (massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx) + massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx));
58 const Scalar throatLength = fluxVarsCache.throatLength();
59 const Scalar phaseCrossSectionalArea = fluxVarsCache.throatCrossSectionalArea(phaseIdx);
60
61 for (int compIdx = 0; compIdx < numComponents; compIdx++)
62 {
63 if(compIdx == phaseIdx)
64 continue;
65
66 auto insideDiffCoeff = getDiffusionCoefficient_(phaseIdx, compIdx, insideVolVars);
67 auto outsideDiffCoeff = getDiffusionCoefficient_(phaseIdx, compIdx, outsideVolVars);
68
69 // scale by extrusion factor
70 insideDiffCoeff *= insideVolVars.extrusionFactor();
71 outsideDiffCoeff *= outsideVolVars.extrusionFactor();
72
73 // the resulting averaged diffusion coefficient
74 const auto diffCoeff = harmonicMean(insideDiffCoeff, outsideDiffCoeff);
75
76 const Scalar insideMoleFraction = massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
77 const Scalar outsideMoleFraction = massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
78
79 componentFlux[compIdx] = density * (insideMoleFraction - outsideMoleFraction) / throatLength * diffCoeff * phaseCrossSectionalArea;
80 componentFlux[phaseIdx] -= componentFlux[compIdx];
81 }
82 return componentFlux;
83 }
84private:
85
86 template<class VolumeVariables>
87 static Scalar getDiffusionCoefficient_(const int phaseIdx, const int compIdx,
88 const VolumeVariables& volVars)
89 {
90 using FluidSystem = typename VolumeVariables::FluidSystem;
91
92 if constexpr (!FluidSystem::isTracerFluidSystem())
93 {
94 const auto mainCompIdx = FluidSystem::getMainComponent(phaseIdx);
95 return volVars.diffusionCoefficient(phaseIdx, mainCompIdx, compIdx);
96 }
97 else
98 return volVars.diffusionCoefficient(0, 0, compIdx);
99 }
100};
101} // end namespace Dumux::PoreNetwork
102
103#endif
Container storing the diffusion coefficients required by Fick's law. Uses the minimal possible contai...
Definition: fickiandiffusioncoefficients.hh:32
Specialization of Fick's Law for the pore-network model.
Definition: flux/porenetwork/fickslaw.hh:30
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: flux/porenetwork/fickslaw.hh:35
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:41
Container storing the diffusion coefficients required by Fick's law. Uses the minimal possible contai...
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
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:54
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:43
ReferenceSystemFormulation
The formulations available for Fick's law related to the reference system.
Definition: referencesystemformulation.hh:33
Define some often used mathematical functions.
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:53
Definition: discretization/porenetwork/fvelementgeometry.hh:24
The reference frameworks and formulations available for splitting total fluxes into a advective and d...