3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
flux/staggered/freeflow/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_DISCRETIZATION_STAGGERED_FICKS_LAW_HH
26#define DUMUX_DISCRETIZATION_STAGGERED_FICKS_LAW_HH
27
28#include <numeric>
29#include <dune/common/fvector.hh>
30
34#include <dumux/common/math.hh>
35
40
41
42namespace Dumux {
43
44// forward declaration
45template<class TypeTag, DiscretizationMethod discMethod, ReferenceSystemFormulation referenceSystem>
46class FicksLawImplementation;
47
52template <class TypeTag, ReferenceSystemFormulation referenceSystem>
53class FicksLawImplementation<TypeTag, DiscretizationMethod::staggered, referenceSystem>
54{
57 using FVElementGeometry = typename GridGeometry::LocalView;
58 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
60 using GridView = typename GridGeometry::GridView;
61 using Element = typename GridView::template Codim<0>::Entity;
63 using Indices = typename ModelTraits::Indices;
65
66 static constexpr int numComponents = ModelTraits::numFluidComponents();
67 static constexpr int numPhases = ModelTraits::numFluidPhases();
68
69 using NumEqVector = Dune::FieldVector<Scalar, numComponents>;
70
71 static_assert(ModelTraits::numFluidPhases() == 1, "Only one phase supported!");
72
73public:
74 // state the discretization method this implementation belongs to
76 //return the reference system
78 { return referenceSystem; }
79
83
85
86 template<class Problem, class ElementVolumeVariables>
87 static NumEqVector flux(const Problem& problem,
88 const Element& element,
89 const FVElementGeometry& fvGeometry,
90 const ElementVolumeVariables& elemVolVars,
91 const SubControlVolumeFace &scvf)
92 {
93 NumEqVector flux(0.0);
94
95 // There is no diffusion over outflow boundaries (grad x == 0).
96 // We assume that if an outflow BC is set for the first transported component, this
97 // also holds for all other components.
98 if (scvf.boundary() && problem.boundaryTypes(element, scvf).isOutflow(Indices::conti0EqIdx + 1))
99 return flux;
100
101 const int phaseIdx = 0;
102
103 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
104 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
105 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
106
107 const Scalar insideDistance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
108 const Scalar insideDensity = massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
109
110 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
111 {
112 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
113 continue;
114
115 const Scalar massOrMoleFractionInside = massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
116 const Scalar massOrMoleFractionOutside = massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
117 const Scalar insideD = getEffectiveDiffusionCoefficient_(insideVolVars, phaseIdx, compIdx) * insideVolVars.extrusionFactor();
118
119 if (scvf.boundary())
120 {
121 flux[compIdx] = insideDensity * insideD
122 * (massOrMoleFractionInside - massOrMoleFractionOutside) / insideDistance;
123 }
124 else
125 {
126 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
127 const Scalar outsideD = getEffectiveDiffusionCoefficient_(outsideVolVars, phaseIdx, compIdx)
128 * outsideVolVars.extrusionFactor();
129 const Scalar outsideDistance = (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm();
130 const Scalar outsideDensity = massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
131
132 const Scalar avgDensity = 0.5*(insideDensity + outsideDensity);
133 const Scalar avgD = harmonicMean(insideD, outsideD, insideDistance, outsideDistance);
134
135 flux[compIdx] = avgDensity * avgD
136 * (massOrMoleFractionInside - massOrMoleFractionOutside) / (insideDistance + outsideDistance);
137 }
138 }
139
140 // Fick's law (for binary systems) states that the net flux of mass within the bulk phase has to be zero:
141 const Scalar cumulativeFlux = std::accumulate(flux.begin(), flux.end(), 0.0);
142 flux[FluidSystem::getMainComponent(0)] = -cumulativeFlux;
143
144 flux *= scvf.area();
145
146 return flux;
147 }
148
149private:
150 static Scalar getEffectiveDiffusionCoefficient_(const VolumeVariables& volVars, const int phaseIdx, const int compIdx)
151 {
152 if constexpr (Dumux::Deprecated::hasEffDiffCoeff<VolumeVariables>)
153 return volVars.effectiveDiffusionCoefficient(phaseIdx,
154 VolumeVariables::FluidSystem::getMainComponent(phaseIdx), compIdx);
155 else
156 {
157 // TODO: remove this else clause after release 3.2!
158 return volVars.effectiveDiffusivity(phaseIdx, compIdx);
159 }
160 }
161};
162} // end namespace
163
164#endif
Define some often used mathematical functions.
Helpers for deprecation.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
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...
Classes related to flux variables caching.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
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:68
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
forward declaration of the method-specific implemetation
Definition: flux/box/fickslaw.hh:43
Container storing the diffusion coefficients required by Fick's law. Uses the minimal possible contai...
Definition: fickiandiffusioncoefficients.hh:44
Definition: fluxvariablescaching.hh:64
static NumEqVector flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Definition: flux/staggered/freeflow/fickslaw.hh:87
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: flux/staggered/freeflow/fickslaw.hh:77
Declares all properties used in Dumux.