3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * 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_FICKS_LAW_HH
26#define DUMUX_FLUX_PNM_FICKS_LAW_HH
27
28#include <dune/common/fvector.hh>
29#include <dumux/common/math.hh>
32
33namespace Dumux::PoreNetwork {
34
39template<class Scalar, int numPhases, int numComponents,
42{
43public:
45
46 //return the reference system
48 { return referenceSystem; }
49
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,
58 const int phaseIdx,
59 const ElementFluxVariablesCache& elemFluxVarsCache)
60 {
61 Dune::FieldVector<Scalar, numComponents> componentFlux(0.0);
62
63 // get inside and outside diffusion tensors and calculate the harmonic mean
64 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
65 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
66
67 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
68
69 const Scalar density = 0.5 * (massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx) + massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx));
70 const Scalar throatLength = fluxVarsCache.throatLength();
71 const Scalar phaseCrossSectionalArea = fluxVarsCache.throatCrossSectionalArea(phaseIdx);
72
73 for (int compIdx = 0; compIdx < numComponents; compIdx++)
74 {
75 if(compIdx == phaseIdx)
76 continue;
77
78 auto insideD = getDiffusionCoefficient_(phaseIdx, compIdx, insideVolVars);
79 auto outsideD = getDiffusionCoefficient_(phaseIdx, compIdx, outsideVolVars);
80
81 // scale by extrusion factor
82 insideD *= insideVolVars.extrusionFactor();
83 outsideD *= outsideVolVars.extrusionFactor();
84
85 // the resulting averaged diffusion coefficient
86 const auto D = harmonicMean(insideD, outsideD);
87
88 const Scalar insideMoleFraction = massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
89 const Scalar outsideMoleFraction = massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
90
91 componentFlux[compIdx] = density * (insideMoleFraction - outsideMoleFraction) / throatLength * D * phaseCrossSectionalArea;
92 componentFlux[phaseIdx] -= componentFlux[compIdx];
93 }
94 return componentFlux;
95 }
96private:
97
98 template<class VolumeVariables>
99 static Scalar getDiffusionCoefficient_(const int phaseIdx, const int compIdx,
100 const VolumeVariables& volVars)
101 {
102 using FluidSystem = typename VolumeVariables::FluidSystem;
103
104 if constexpr (!FluidSystem::isTracerFluidSystem())
105 {
106 const auto mainCompIdx = FluidSystem::getMainComponent(phaseIdx);
107 return volVars.diffusionCoefficient(phaseIdx, mainCompIdx, compIdx);
108 }
109 else
110 return volVars.diffusionCoefficient(0, 0, compIdx);
111 }
112};
113} // end namespace Dumux::PoreNetwork
114
115#endif
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:34
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