3.2-git
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{
56 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
57 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
58 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
59 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
60 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
64 using Element = typename GridView::template Codim<0>::Entity;
66
67 enum { dim = GridView::dimension} ;
68 enum { dimWorld = GridView::dimensionworld} ;
69 enum
70 {
71 numPhases = ModelTraits::numFluidPhases(),
72 numComponents = ModelTraits::numFluidComponents()
73 };
74 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
75 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
76
77public:
79
80 //return the reference system
82 { return referenceSystem; }
83
84 static ComponentFluxVector flux(const Problem& problem,
85 const Element& element,
86 const FVElementGeometry& fvGeometry,
87 const ElementVolumeVariables& elemVolVars,
88 const SubControlVolumeFace& scvf,
89 const int phaseIdx,
90 const ElementFluxVariablesCache& elemFluxVarsCache)
91 {
92 // get inside and outside diffusion tensors and calculate the harmonic mean
93 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
94 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
95
96 // evaluate gradX at integration point and interpolate density
97 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
98 const auto& shapeValues = fluxVarsCache.shapeValues();
99
100 // density interpolation
101 Scalar rho(0.0);
102 for (auto&& scv : scvs(fvGeometry))
103 rho += massOrMolarDensity(elemVolVars[scv], referenceSystem, phaseIdx)*shapeValues[scv.indexInElement()][0];
104
105 ComponentFluxVector componentFlux(0.0);
106 for (int compIdx = 0; compIdx < numComponents; compIdx++)
107 {
108 if constexpr (!FluidSystem::isTracerFluidSystem())
109 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
110 continue;
111
112 const auto D = averageDiffusionCoefficient_(phaseIdx, compIdx, insideVolVars, outsideVolVars, problem, scvf);
113
114 // compute the diffusive flux
115 const auto massOrMoleFrac = [&](const SubControlVolume& scv){ return massOrMoleFraction(elemVolVars[scv], referenceSystem, phaseIdx, compIdx); };
116 componentFlux[compIdx] = discreteFlux_(fvGeometry, scvf, fluxVarsCache, massOrMoleFrac, D, rho);
117
118 // if the main component is balanced subtract the same flux from there (conservation)
119 if constexpr (!FluidSystem::isTracerFluidSystem())
120 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx))
121 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
122 }
123 return componentFlux;
124 }
125
126 // compute transmissibilities ti for analytical jacobians
127 static std::array<std::vector<Scalar>, numComponents>
128 calculateTransmissibilities(const Problem& problem,
129 const Element& element,
130 const FVElementGeometry& fvGeometry,
131 const ElementVolumeVariables& elemVolVars,
132 const SubControlVolumeFace& scvf,
133 const FluxVarCache& fluxVarCache,
134 const int phaseIdx)
135 {
136 Scalar rho(0.0);
137 const auto& shapeValues = fluxVarCache.shapeValues();
138 for (auto&& scv : scvs(fvGeometry))
139 rho += massOrMolarDensity(elemVolVars[scv], referenceSystem, phaseIdx)*shapeValues[scv.indexInElement()][0];
140
141 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
142 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
143
144 std::array<std::vector<Scalar>, numComponents> ti;
145 for (int compIdx = 0; compIdx < numComponents; compIdx++)
146 {
147 if constexpr (!FluidSystem::isTracerFluidSystem())
148 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
149 continue;
150
151 const auto D = averageDiffusionCoefficient_(phaseIdx, compIdx, insideVolVars, outsideVolVars, problem, scvf);
152
153 ti[compIdx].resize(fvGeometry.numScv());
154 for (auto&& scv : scvs(fvGeometry))
155 ti[compIdx][scv.indexInElement()] = -rho*vtmv(scvf.unitOuterNormal(), D, fluxVarCache.gradN(scv.indexInElement()))*scvf.area();
156 }
157
158 return ti;
159 }
160
161private:
162 static Scalar averageDiffusionCoefficient_(const int phaseIdx, const int compIdx,
163 const VolumeVariables& insideVV, const VolumeVariables& outsideVV,
164 const Problem& problem,
165 const SubControlVolumeFace& scvf)
166 {
167 // effective diffusion tensors
168 auto [insideD, outsideD] = diffusionCoefficientsAtInterface_(phaseIdx, compIdx, insideVV, outsideVV);
169
170 // scale by extrusion factor
171 insideD *= insideVV.extrusionFactor();
172 outsideD *= outsideVV.extrusionFactor();
173
174 // the resulting averaged diffusion tensor
175 return problem.spatialParams().harmonicMean(insideD, outsideD, scvf.unitOuterNormal());
176 }
177
178 static std::pair<Scalar, Scalar>
179 diffusionCoefficientsAtInterface_([[maybe_unused]] const int phaseIdx, const int compIdx,
180 const VolumeVariables& insideVV, const VolumeVariables& outsideVV)
181 {
182 if constexpr (!FluidSystem::isTracerFluidSystem())
183 {
184 using EffDiffModel = GetPropType<TypeTag, Properties::EffectiveDiffusivityModel>;
185 const auto mainCompIdx = FluidSystem::getMainComponent(phaseIdx);
186 const auto insideD = Deprecated::template effectiveDiffusionCoefficient<EffDiffModel>(insideVV, phaseIdx, mainCompIdx, compIdx);
187 const auto outsideD = Deprecated::template effectiveDiffusionCoefficient<EffDiffModel>(outsideVV, phaseIdx, mainCompIdx, compIdx);
188 return { std::move(insideD), std::move(outsideD) };
189 }
190 else
191 {
192 using EffDiffModel = GetPropType<TypeTag, Properties::EffectiveDiffusivityModel>;
193 const auto insideD = Deprecated::template effectiveDiffusionCoefficient<EffDiffModel>(insideVV, 0, 0, compIdx);
194 const auto outsideD = Deprecated::template effectiveDiffusionCoefficient<EffDiffModel>(outsideVV, 0, 0, compIdx);
195 return { std::move(insideD), std::move(outsideD) };
196 }
197 }
198
199 template<class EvaluateVariable, class Tensor>
200 static Scalar discreteFlux_(const FVElementGeometry& fvGeometry,
201 const SubControlVolumeFace& scvf,
202 const FluxVarCache& fluxVarsCache,
203 const EvaluateVariable& massOrMoleFraction,
204 const Tensor& D, const Scalar preFactor)
205 {
206 Dune::FieldVector<Scalar, dimWorld> gradX(0.0);
207 for (auto&& scv : scvs(fvGeometry))
208 gradX.axpy(massOrMoleFraction(scv), fluxVarsCache.gradN(scv.indexInElement()));
209 return -1.0*preFactor*vtmv(scvf.unitOuterNormal(), D, gradX)*scvf.area();
210 }
211};
212
213} // end namespace Dumux
214
215#endif
Define some often used mathematical functions.
Helpers for deprecation.
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...
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:840
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
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: flux/box/fickslaw.hh:81
static ComponentFluxVector flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVariablesCache &elemFluxVarsCache)
Definition: flux/box/fickslaw.hh:84
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:128
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.