3.6-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
33#include <dumux/common/math.hh>
34
40
41
42namespace Dumux {
43
44// forward declaration
45template<class TypeTag, class DiscretizationMethod, ReferenceSystemFormulation referenceSystem>
46class FicksLawImplementation;
47
52template <class TypeTag, ReferenceSystemFormulation referenceSystem>
53class FicksLawImplementation<TypeTag, DiscretizationMethods::Staggered, referenceSystem>
54{
57 using FVElementGeometry = typename GridGeometry::LocalView;
58 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
59 using Extrusion = Extrusion_t<GridGeometry>;
61 using GridView = typename GridGeometry::GridView;
62 using Element = typename GridView::template Codim<0>::Entity;
64 using Indices = typename ModelTraits::Indices;
66
67 static constexpr int numComponents = ModelTraits::numFluidComponents();
68 static constexpr int numPhases = ModelTraits::numFluidPhases();
69
70 using NumEqVector = Dune::FieldVector<Scalar, numComponents>;
71
72 static_assert(ModelTraits::numFluidPhases() == 1, "Only one phase supported!");
73
74public:
76 // state the discretization method this implementation belongs to
77 static constexpr DiscretizationMethod discMethod{};
78
79 //return the reference system
81 { return referenceSystem; }
82
86
88
95 template<class Problem, class ElementVolumeVariables>
96 static NumEqVector flux(const Problem& problem,
97 const Element& element,
98 const FVElementGeometry& fvGeometry,
99 const ElementVolumeVariables& elemVolVars,
100 const SubControlVolumeFace &scvf)
101 {
102 NumEqVector flux(0.0);
103
104 // There is no diffusion over outflow boundaries (grad x == 0).
105 // We assume that if an outflow BC is set for the first transported component, this
106 // also holds for all other components.
107 if (scvf.boundary() && problem.boundaryTypes(element, scvf).isOutflow(Indices::conti0EqIdx + 1))
108 return flux;
109
110 const int phaseIdx = 0;
111
112 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
113 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
114 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
115
116 const Scalar insideDistance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
117 const Scalar insideDensity = massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
118
119 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
120 {
121 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
122 continue;
123
124 const Scalar massOrMoleFractionInside = massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
125 const Scalar massOrMoleFractionOutside = massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
126 const Scalar insideDiffCoeff = getEffectiveDiffusionCoefficient_(insideVolVars, phaseIdx, compIdx) * insideVolVars.extrusionFactor();
127
128 if (scvf.boundary())
129 {
130 flux[compIdx] = insideDensity * insideDiffCoeff
131 * (massOrMoleFractionInside - massOrMoleFractionOutside) / insideDistance;
132 }
133 else
134 {
135 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
136 const Scalar outsideDiffCoeff = getEffectiveDiffusionCoefficient_(outsideVolVars, phaseIdx, compIdx)
137 * outsideVolVars.extrusionFactor();
138 const Scalar outsideDistance = (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm();
139 const Scalar outsideDensity = massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
140
141 const Scalar avgDensity = 0.5*(insideDensity + outsideDensity);
142 const Scalar avgDiffCoeff = harmonicMean(insideDiffCoeff, outsideDiffCoeff, insideDistance, outsideDistance);
143
144 flux[compIdx] = avgDensity * avgDiffCoeff
145 * (massOrMoleFractionInside - massOrMoleFractionOutside) / (insideDistance + outsideDistance);
146 }
147 }
148
149 // Fick's law (for binary systems) states that the net flux of mass within the bulk phase has to be zero:
150 const Scalar cumulativeFlux = std::accumulate(flux.begin(), flux.end(), 0.0);
151 flux[FluidSystem::getMainComponent(0)] = -cumulativeFlux;
152
153 flux *= Extrusion::area(fvGeometry, scvf);
154
155 return flux;
156 }
157
158private:
159 static Scalar getEffectiveDiffusionCoefficient_(const VolumeVariables& volVars, const int phaseIdx, const int compIdx)
160 {
161 return volVars.effectiveDiffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
162 }
163};
164} // end namespace Dumux
165
166#endif
Classes related to flux variables caching.
Container storing the diffusion coefficients required by Fick's law. Uses the minimal possible contai...
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
Helper classes to compute the integration elements.
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
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
forward declaration of the method-specific implementation
Definition: flux/box/fickslaw.hh:44
Container storing the diffusion coefficients required by Fick's law. Uses the minimal possible contai...
Definition: fickiandiffusioncoefficients.hh:44
Definition: fluxvariablescaching.hh:67
static NumEqVector flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Returns the diffusive fluxes of all components within a fluid phase across the given sub-control volu...
Definition: flux/staggered/freeflow/fickslaw.hh:96
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: flux/staggered/freeflow/fickslaw.hh:80
Declares all properties used in Dumux.