3.6-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
39
40namespace Dumux {
41
42// forward declaration
43template<class TypeTag, class DiscretizationMethod, ReferenceSystemFormulation referenceSystem>
45
50template <class TypeTag, ReferenceSystemFormulation referenceSystem>
51class FicksLawImplementation<TypeTag, DiscretizationMethods::Box, referenceSystem>
52{
58 using FVElementGeometry = typename GridGeometry::LocalView;
59 using SubControlVolume = typename GridGeometry::SubControlVolume;
60 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
61 using Extrusion = Extrusion_t<GridGeometry>;
62 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
64 using ElementFluxVariablesCache = typename GridFluxVariablesCache::LocalView;
65 using FluxVarCache = typename GridFluxVariablesCache::FluxVariablesCache;
68 using Element = typename GridView::template Codim<0>::Entity;
70
71 enum { dim = GridView::dimension} ;
72 enum { dimWorld = GridView::dimensionworld} ;
73 enum
74 {
75 numPhases = ModelTraits::numFluidPhases(),
76 numComponents = ModelTraits::numFluidComponents()
77 };
78 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
79 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
80
81public:
83
84 //return the reference system
86 { return referenceSystem; }
87
94 static ComponentFluxVector flux(const Problem& problem,
95 const Element& element,
96 const FVElementGeometry& fvGeometry,
97 const ElementVolumeVariables& elemVolVars,
98 const SubControlVolumeFace& scvf,
99 const int phaseIdx,
100 const ElementFluxVariablesCache& elemFluxVarsCache)
101 {
102 // get inside and outside diffusion tensors and calculate the harmonic mean
103 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
104 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
105
106 // evaluate gradX at integration point and interpolate density
107 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
108 const auto& shapeValues = fluxVarsCache.shapeValues();
109
110 // density interpolation
111 Scalar rho(0.0);
112 for (auto&& scv : scvs(fvGeometry))
113 rho += massOrMolarDensity(elemVolVars[scv], referenceSystem, phaseIdx)*shapeValues[scv.indexInElement()][0];
114
115 ComponentFluxVector componentFlux(0.0);
116 for (int compIdx = 0; compIdx < numComponents; compIdx++)
117 {
118 if constexpr (!FluidSystem::isTracerFluidSystem())
119 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
120 continue;
121
122 const auto diffCoeff = averageDiffusionCoefficient_(phaseIdx, compIdx, insideVolVars, outsideVolVars, problem, scvf);
123
124 // compute the diffusive flux
125 const auto massOrMoleFrac = [&](const SubControlVolume& scv){ return massOrMoleFraction(elemVolVars[scv], referenceSystem, phaseIdx, compIdx); };
126 componentFlux[compIdx] = discreteFlux_(fvGeometry, scvf, fluxVarsCache, massOrMoleFrac, diffCoeff, rho);
127
128 // if the main component is balanced subtract the same flux from there (conservation)
129 if constexpr (!FluidSystem::isTracerFluidSystem())
130 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx))
131 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
132 }
133 return componentFlux;
134 }
135
136 // compute transmissibilities ti for analytical jacobians
137 static std::array<std::vector<Scalar>, numComponents>
138 calculateTransmissibilities(const Problem& problem,
139 const Element& element,
140 const FVElementGeometry& fvGeometry,
141 const ElementVolumeVariables& elemVolVars,
142 const SubControlVolumeFace& scvf,
143 const FluxVarCache& fluxVarCache,
144 const int phaseIdx)
145 {
146 Scalar rho(0.0);
147 const auto& shapeValues = fluxVarCache.shapeValues();
148 for (auto&& scv : scvs(fvGeometry))
149 rho += massOrMolarDensity(elemVolVars[scv], referenceSystem, phaseIdx)*shapeValues[scv.indexInElement()][0];
150
151 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
152 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
153
154 std::array<std::vector<Scalar>, numComponents> ti;
155 for (int compIdx = 0; compIdx < numComponents; compIdx++)
156 {
157 if constexpr (!FluidSystem::isTracerFluidSystem())
158 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
159 continue;
160
161 const auto diffCoeff = averageDiffusionCoefficient_(phaseIdx, compIdx, insideVolVars, outsideVolVars, problem, scvf);
162
163 ti[compIdx].resize(fvGeometry.numScv());
164 for (auto&& scv : scvs(fvGeometry))
165 ti[compIdx][scv.indexInElement()] = -rho*vtmv(scvf.unitOuterNormal(), diffCoeff, fluxVarCache.gradN(scv.indexInElement()))*Extrusion::area(fvGeometry, scvf);
166 }
167
168 return ti;
169 }
170
171private:
172 static Scalar averageDiffusionCoefficient_(const int phaseIdx, const int compIdx,
173 const VolumeVariables& insideVV, const VolumeVariables& outsideVV,
174 const Problem& problem,
175 const SubControlVolumeFace& scvf)
176 {
177 // effective diffusion tensors
178 auto [insideDiffCoeff, outsideDiffCoeff] = diffusionCoefficientsAtInterface_(phaseIdx, compIdx, insideVV, outsideVV);
179
180 // scale by extrusion factor
181 insideDiffCoeff *= insideVV.extrusionFactor();
182 outsideDiffCoeff *= outsideVV.extrusionFactor();
183
184 // the resulting averaged diffusion tensor
185 return faceTensorAverage(insideDiffCoeff, outsideDiffCoeff, scvf.unitOuterNormal());
186 }
187
188 static std::pair<Scalar, Scalar>
189 diffusionCoefficientsAtInterface_([[maybe_unused]] const int phaseIdx, const int compIdx,
190 const VolumeVariables& insideVV, const VolumeVariables& outsideVV)
191 {
192 if constexpr (!FluidSystem::isTracerFluidSystem())
193 {
194 const auto mainCompIdx = FluidSystem::getMainComponent(phaseIdx);
195 const auto insideDiffCoeff = insideVV.effectiveDiffusionCoefficient(phaseIdx, mainCompIdx, compIdx);
196 const auto outsideDiffCoeff = outsideVV.effectiveDiffusionCoefficient(phaseIdx, mainCompIdx, compIdx);
197 return { std::move(insideDiffCoeff), std::move(outsideDiffCoeff) };
198 }
199 else
200 {
201 const auto insideDiffCoeff = insideVV.effectiveDiffusionCoefficient(0, 0, compIdx);
202 const auto outsideDiffCoeff = outsideVV.effectiveDiffusionCoefficient(0, 0, compIdx);
203 return { std::move(insideDiffCoeff), std::move(outsideDiffCoeff) };
204 }
205 }
206
207 template<class EvaluateVariable, class Tensor>
208 static Scalar discreteFlux_(const FVElementGeometry& fvGeometry,
209 const SubControlVolumeFace& scvf,
210 const FluxVarCache& fluxVarsCache,
211 const EvaluateVariable& massOrMoleFraction,
212 const Tensor& D, const Scalar preFactor)
213 {
214 Dune::FieldVector<Scalar, dimWorld> gradX(0.0);
215 for (auto&& scv : scvs(fvGeometry))
216 gradX.axpy(massOrMoleFraction(scv), fluxVarsCache.gradN(scv.indexInElement()));
217 return -1.0*preFactor*vtmv(scvf.unitOuterNormal(), D, gradX)*Extrusion::area(fvGeometry, scvf);
218 }
219};
220
221} // end namespace Dumux
222
223#endif
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...
A free function to average a Tensor at an interface.
Define some often used mathematical functions.
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
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:863
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
Scalar faceTensorAverage(const Scalar T1, const Scalar T2, const Dune::FieldVector< Scalar, dim > &normal)
Average of a discontinuous scalar field at discontinuity interface (for compatibility reasons with th...
Definition: facetensoraverage.hh:41
CVFE< CVFEMethods::PQ1 > Box
Definition: method.hh:83
forward declaration of the method-specific implementation
Definition: flux/box/fickslaw.hh:44
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:138
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: flux/box/fickslaw.hh:85
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:94
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.