3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 <dumux/common/math.hh>
31
33
34namespace Dumux {
35
36// forward declaration
37template<class TypeTag, DiscretizationMethod discMethod, ReferenceSystemFormulation referenceSystem>
39
44template <class TypeTag, ReferenceSystemFormulation referenceSystem>
45class FicksLawImplementation<TypeTag, DiscretizationMethod::box, referenceSystem>
46{
51 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
52 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
53 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
55 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
56 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
60 using Element = typename GridView::template Codim<0>::Entity;
62 using Indices = typename ModelTraits::Indices;
63
64 enum { dim = GridView::dimension} ;
65 enum { dimWorld = GridView::dimensionworld} ;
66 enum
67 {
68 numPhases = ModelTraits::numFluidPhases(),
69 numComponents = ModelTraits::numFluidComponents()
70 };
71 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
72 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
73
74public:
75 //return the reference system
77 { return referenceSystem; }
78
79 static ComponentFluxVector flux(const Problem& problem,
80 const Element& element,
81 const FVElementGeometry& fvGeometry,
82 const ElementVolumeVariables& elemVolVars,
83 const SubControlVolumeFace& scvf,
84 const int phaseIdx,
85 const ElementFluxVariablesCache& elemFluxVarsCache)
86 {
87 ComponentFluxVector componentFlux(0.0);
88
89 // get inside and outside diffusion tensors and calculate the harmonic mean
90 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
91 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
92
93 // evaluate gradX at integration point and interpolate density
94 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
95 const auto& shapeValues = fluxVarsCache.shapeValues();
96
97 // density interpolation
98 Scalar rho(0.0);
99 for (auto&& scv : scvs(fvGeometry))
100 rho += massOrMolarDensity(elemVolVars[scv], referenceSystem, phaseIdx)*shapeValues[scv.indexInElement()][0];
101
102 for (int compIdx = 0; compIdx < numComponents; compIdx++)
103 {
104 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
105 continue;
106
107 // effective diffusion tensors
109 auto insideD = EffDiffModel::effectiveDiffusivity(insideVolVars.porosity(),
110 insideVolVars.saturation(phaseIdx),
111 insideVolVars.diffusionCoefficient(phaseIdx, compIdx));
112 auto outsideD = EffDiffModel::effectiveDiffusivity(outsideVolVars.porosity(),
113 outsideVolVars.saturation(phaseIdx),
114 outsideVolVars.diffusionCoefficient(phaseIdx, compIdx));
115
116 // scale by extrusion factor
117 insideD *= insideVolVars.extrusionFactor();
118 outsideD *= outsideVolVars.extrusionFactor();
119
120 // the resulting averaged diffusion tensor
121 const auto D = problem.spatialParams().harmonicMean(insideD, outsideD, scvf.unitOuterNormal());
122
123 // the mole/mass fraction gradient
124 Dune::FieldVector<Scalar, dimWorld> gradX(0.0);
125 for (auto&& scv : scvs(fvGeometry))
126 gradX.axpy(massOrMoleFraction(elemVolVars[scv], referenceSystem, phaseIdx, compIdx), fluxVarsCache.gradN(scv.indexInElement()));
127
128 // compute the diffusive flux
129 componentFlux[compIdx] = -1.0*rho*vtmv(scvf.unitOuterNormal(), D, gradX)*scvf.area();
130 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx) && !FluidSystem::isTracerFluidSystem())
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(compIdx == FluidSystem::getMainComponent(phaseIdx))
158 continue;
159
160 // effective diffusion tensors
162 auto insideD = EffDiffModel::effectiveDiffusivity(insideVolVars.porosity(),
163 insideVolVars.saturation(phaseIdx),
164 insideVolVars.diffusionCoefficient(phaseIdx, compIdx));
165 auto outsideD = EffDiffModel::effectiveDiffusivity(outsideVolVars.porosity(),
166 outsideVolVars.saturation(phaseIdx),
167 outsideVolVars.diffusionCoefficient(phaseIdx, compIdx));
168
169 // scale by extrusion factor
170 insideD *= insideVolVars.extrusionFactor();
171 outsideD *= outsideVolVars.extrusionFactor();
172
173 // the resulting averaged diffusion tensor
174 const auto D = problem.spatialParams().harmonicMean(insideD, outsideD, scvf.unitOuterNormal());
175
176 ti[compIdx].resize(fvGeometry.numScv());
177 for (auto&& scv : scvs(fvGeometry))
178 ti[compIdx][scv.indexInElement()] = -rho*vtmv(scvf.unitOuterNormal(), D, fluxVarCache.gradN(scv.indexInElement()))*scvf.area();
179 }
180
181 return ti;
182 }
183};
184
185} // end namespace Dumux
186
187#endif
Define some often used mathematical functions.
The available discretization methods in Dumux.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
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
make the local view function available whenever we use the grid geometry
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
Property tag Indices
Definition: porousmediumflow/sequential/properties.hh:59
forward declaration of the method-specific implemetation
Definition: box/fickslaw.hh:38
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: box/fickslaw.hh:76
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: box/fickslaw.hh:79
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: box/fickslaw.hh:138
Declares all properties used in Dumux.