3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
flux/box/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_BOX_FICKS_LAW_HH
26#define DUMUX_DISCRETIZATION_BOX_FICKS_LAW_HH
27
28#include <dune/common/fvector.hh>
29#include <dune/common/fmatrix.hh>
30
31#include <dumux/common/math.hh>
35
38
39namespace Dumux {
40
41// forward declaration
42template<class TypeTag, DiscretizationMethod discMethod, ReferenceSystemFormulation referenceSystem>
44
49template <class TypeTag, ReferenceSystemFormulation referenceSystem>
50class FicksLawImplementation<TypeTag, DiscretizationMethod::box, referenceSystem>
51{
57 using FVElementGeometry = typename GridGeometry::LocalView;
58 using SubControlVolume = typename GridGeometry::SubControlVolume;
59 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
60 using Extrusion = Extrusion_t<GridGeometry>;
61 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
62 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
66 using Element = typename GridView::template Codim<0>::Entity;
68
69 enum { dim = GridView::dimension} ;
70 enum { dimWorld = GridView::dimensionworld} ;
71 enum
72 {
73 numPhases = ModelTraits::numFluidPhases(),
74 numComponents = ModelTraits::numFluidComponents()
75 };
76 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
77 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
78
79public:
81
82 //return the reference system
84 { return referenceSystem; }
85
92 static ComponentFluxVector flux(const Problem& problem,
93 const Element& element,
94 const FVElementGeometry& fvGeometry,
95 const ElementVolumeVariables& elemVolVars,
96 const SubControlVolumeFace& scvf,
97 const int phaseIdx,
98 const ElementFluxVariablesCache& elemFluxVarsCache)
99 {
100 // get inside and outside diffusion tensors and calculate the harmonic mean
101 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
102 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
103
104 // evaluate gradX at integration point and interpolate density
105 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
106 const auto& shapeValues = fluxVarsCache.shapeValues();
107
108 // density interpolation
109 Scalar rho(0.0);
110 for (auto&& scv : scvs(fvGeometry))
111 rho += massOrMolarDensity(elemVolVars[scv], referenceSystem, phaseIdx)*shapeValues[scv.indexInElement()][0];
112
113 ComponentFluxVector componentFlux(0.0);
114 for (int compIdx = 0; compIdx < numComponents; compIdx++)
115 {
116 if constexpr (!FluidSystem::isTracerFluidSystem())
117 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
118 continue;
119
120 const auto D = averageDiffusionCoefficient_(phaseIdx, compIdx, insideVolVars, outsideVolVars, problem, scvf);
121
122 // compute the diffusive flux
123 const auto massOrMoleFrac = [&](const SubControlVolume& scv){ return massOrMoleFraction(elemVolVars[scv], referenceSystem, phaseIdx, compIdx); };
124 componentFlux[compIdx] = discreteFlux_(fvGeometry, scvf, fluxVarsCache, massOrMoleFrac, D, rho);
125
126 // if the main component is balanced subtract the same flux from there (conservation)
127 if constexpr (!FluidSystem::isTracerFluidSystem())
128 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx))
129 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
130 }
131 return componentFlux;
132 }
133
134 // compute transmissibilities ti for analytical jacobians
135 static std::array<std::vector<Scalar>, numComponents>
136 calculateTransmissibilities(const Problem& problem,
137 const Element& element,
138 const FVElementGeometry& fvGeometry,
139 const ElementVolumeVariables& elemVolVars,
140 const SubControlVolumeFace& scvf,
141 const FluxVarCache& fluxVarCache,
142 const int phaseIdx)
143 {
144 Scalar rho(0.0);
145 const auto& shapeValues = fluxVarCache.shapeValues();
146 for (auto&& scv : scvs(fvGeometry))
147 rho += massOrMolarDensity(elemVolVars[scv], referenceSystem, phaseIdx)*shapeValues[scv.indexInElement()][0];
148
149 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
150 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
151
152 std::array<std::vector<Scalar>, numComponents> ti;
153 for (int compIdx = 0; compIdx < numComponents; compIdx++)
154 {
155 if constexpr (!FluidSystem::isTracerFluidSystem())
156 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
157 continue;
158
159 const auto D = averageDiffusionCoefficient_(phaseIdx, compIdx, insideVolVars, outsideVolVars, problem, scvf);
160
161 ti[compIdx].resize(fvGeometry.numScv());
162 for (auto&& scv : scvs(fvGeometry))
163 ti[compIdx][scv.indexInElement()] = -rho*vtmv(scvf.unitOuterNormal(), D, fluxVarCache.gradN(scv.indexInElement()))*Extrusion::area(scvf);
164 }
165
166 return ti;
167 }
168
169private:
170 static Scalar averageDiffusionCoefficient_(const int phaseIdx, const int compIdx,
171 const VolumeVariables& insideVV, const VolumeVariables& outsideVV,
172 const Problem& problem,
173 const SubControlVolumeFace& scvf)
174 {
175 // effective diffusion tensors
176 auto [insideD, outsideD] = diffusionCoefficientsAtInterface_(phaseIdx, compIdx, insideVV, outsideVV);
177
178 // scale by extrusion factor
179 insideD *= insideVV.extrusionFactor();
180 outsideD *= outsideVV.extrusionFactor();
181
182 // the resulting averaged diffusion tensor
183 return problem.spatialParams().harmonicMean(insideD, outsideD, scvf.unitOuterNormal());
184 }
185
186 static std::pair<Scalar, Scalar>
187 diffusionCoefficientsAtInterface_([[maybe_unused]] const int phaseIdx, const int compIdx,
188 const VolumeVariables& insideVV, const VolumeVariables& outsideVV)
189 {
190 if constexpr (!FluidSystem::isTracerFluidSystem())
191 {
192 const auto mainCompIdx = FluidSystem::getMainComponent(phaseIdx);
193 const auto insideD = insideVV.effectiveDiffusionCoefficient(phaseIdx, mainCompIdx, compIdx);
194 const auto outsideD = outsideVV.effectiveDiffusionCoefficient(phaseIdx, mainCompIdx, compIdx);
195 return { std::move(insideD), std::move(outsideD) };
196 }
197 else
198 {
199 const auto insideD = insideVV.effectiveDiffusionCoefficient(0, 0, compIdx);
200 const auto outsideD = outsideVV.effectiveDiffusionCoefficient(0, 0, compIdx);
201 return { std::move(insideD), std::move(outsideD) };
202 }
203 }
204
205 template<class EvaluateVariable, class Tensor>
206 static Scalar discreteFlux_(const FVElementGeometry& fvGeometry,
207 const SubControlVolumeFace& scvf,
208 const FluxVarCache& fluxVarsCache,
209 const EvaluateVariable& massOrMoleFraction,
210 const Tensor& D, const Scalar preFactor)
211 {
212 Dune::FieldVector<Scalar, dimWorld> gradX(0.0);
213 for (auto&& scv : scvs(fvGeometry))
214 gradX.axpy(massOrMoleFraction(scv), fluxVarsCache.gradN(scv.indexInElement()));
215 return -1.0*preFactor*vtmv(scvf.unitOuterNormal(), D, gradX)*Extrusion::area(scvf);
216 }
217};
218
219} // end namespace Dumux
220
221#endif
Define some often used mathematical functions.
The available discretization methods in Dumux.
Helper classes to compute the integration elements.
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...
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
Dune::DenseMatrix< MAT >::value_type vtmv(const Dune::DenseVector< V1 > &v1, const Dune::DenseMatrix< MAT > &M, const Dune::DenseVector< V2 > &v2)
Evaluates the scalar product of a vector v2, projected by a matrix M, with a vector v1.
Definition: math.hh:849
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
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
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: flux/box/fickslaw.hh:83
static ComponentFluxVector flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVariablesCache &elemFluxVarsCache)
Returns the diffusive fluxes of all components within a fluid phase across the given sub-control volu...
Definition: flux/box/fickslaw.hh:92
static std::array< std::vector< Scalar >, numComponents > calculateTransmissibilities(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const FluxVarCache &fluxVarCache, const int phaseIdx)
Definition: flux/box/fickslaw.hh:136
Container storing the diffusion coefficients required by Fick's law. Uses the minimal possible contai...
Definition: fickiandiffusioncoefficients.hh:44
Declares all properties used in Dumux.