3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
box/maxwellstefanslaw.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_MAXWELL_STEFAN_LAW_HH
26#define DUMUX_DISCRETIZATION_BOX_MAXWELL_STEFAN_LAW_HH
27
28#include <dune/common/float_cmp.hh>
29#include <dune/common/fmatrix.hh>
30
35
39
40namespace Dumux {
41
42// forward declaration
43template <class TypeTag, DiscretizationMethod discMethod, ReferenceSystemFormulation referenceSystem>
44class MaxwellStefansLawImplementation;
45
50template <class TypeTag, ReferenceSystemFormulation referenceSystem>
52{
57 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
58 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
59 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
60 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
62 using Element = typename GridView::template Codim<0>::Entity;
63
64 static constexpr auto numFluidPhases = GetPropType<TypeTag, Properties::ModelTraits>::numFluidPhases();
65 static constexpr auto numComponents = GetPropType<TypeTag, Properties::ModelTraits>::numFluidComponents();
66
67 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
68 using ReducedComponentVector = Dune::FieldVector<Scalar, numComponents-1>;
69 using ReducedComponentMatrix = Dune::FieldMatrix<Scalar, numComponents-1, numComponents-1>;
70
71 static_assert(referenceSystem == ReferenceSystemFormulation::massAveraged, "only the mass averaged reference system is supported for the Maxwell-Stefan formulation");
72
73public:
75
76 //return the reference system
78 { return referenceSystem; }
79
80 static ComponentFluxVector flux(const Problem& problem,
81 const Element& element,
82 const FVElementGeometry& fvGeometry,
83 const ElementVolumeVariables& elemVolVars,
84 const SubControlVolumeFace& scvf,
85 const int phaseIdx,
86 const ElementFluxVariablesCache& elemFluxVarsCache)
87 {
88 //this is to calculate the maxwellStefan diffusion in a multicomponent system.
89 //see: Multicomponent Mass Transfer. R. Taylor u. R. Krishna. J. Wiley & Sons, New York 1993
90 ComponentFluxVector componentFlux(0.0);
91 ReducedComponentMatrix reducedDiffusionMatrix(0.0);
92 ReducedComponentVector reducedFlux(0.0);
93 ComponentFluxVector moleFrac(0.0);
94 ReducedComponentVector normalX(0.0);
95
96 // evaluate gradX at integration point and interpolate density
97 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
98 const auto& shapeValues = fluxVarsCache.shapeValues();
99
100 Scalar rho(0.0);
101 Scalar avgMolarMass(0.0);
102 for (auto&& scv : scvs(fvGeometry))
103 {
104 const auto& volVars = elemVolVars[scv];
105
106 // density interpolation
107 rho += volVars.density(phaseIdx)*shapeValues[scv.indexInElement()][0];
108 //average molar mass interpolation
109 avgMolarMass += volVars.averageMolarMass(phaseIdx)*shapeValues[scv.indexInElement()][0];
110 //interpolate the mole fraction for the diffusion matrix
111 for (int compIdx = 0; compIdx < numComponents; compIdx++)
112 {
113 moleFrac[compIdx] += volVars.moleFraction(phaseIdx, compIdx)*shapeValues[scv.indexInElement()][0];
114 }
115 }
116
117 reducedDiffusionMatrix = setupMSMatrix_(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, moleFrac, avgMolarMass);
118
119 for (int compIdx = 0; compIdx < numComponents-1; compIdx++)
120 {
121 Dune::FieldVector<Scalar, GridView::dimensionworld> gradX(0.0);
122 for (auto&& scv : scvs(fvGeometry))
123 {
124 const auto& volVars = elemVolVars[scv];
125
126 // the mole/mass fraction gradient
127 gradX.axpy(volVars.moleFraction(phaseIdx, compIdx), fluxVarsCache.gradN(scv.indexInElement()));
128 }
129
130 normalX[compIdx] = gradX *scvf.unitOuterNormal();
131 }
132 reducedDiffusionMatrix.solve(reducedFlux,normalX);
133 reducedFlux *= -1.0*rho*scvf.area();
134
135 for (int compIdx = 0; compIdx < numComponents-1; compIdx++)
136 {
137 componentFlux[compIdx] = reducedFlux[compIdx];
138 componentFlux[numComponents-1] -= reducedFlux[compIdx];
139 }
140 return componentFlux ;
141 }
142
143private:
144
145 static ReducedComponentMatrix setupMSMatrix_(const Problem& problem,
146 const Element& element,
147 const FVElementGeometry& fvGeometry,
148 const ElementVolumeVariables& elemVolVars,
149 const SubControlVolumeFace& scvf,
150 const int phaseIdx,
151 const ComponentFluxVector moleFrac,
152 const Scalar avgMolarMass)
153 {
155 ReducedComponentMatrix reducedDiffusionMatrix(0.0);
156
157 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
158 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
159 const auto insideScvIdx = scvf.insideScvIdx();
160 const auto outsideScvIdx = scvf.outsideScvIdx();
161
162 //this is to not devide by 0 if the saturation in 0 and the effectiveDiffusionCoefficient becomes zero due to that
163 if(Dune::FloatCmp::eq<Scalar>(insideVolVars.saturation(phaseIdx), 0) || Dune::FloatCmp::eq<Scalar>(outsideVolVars.saturation(phaseIdx), 0))
164 return reducedDiffusionMatrix;
165
166 for (int compIIdx = 0; compIIdx < numComponents-1; compIIdx++)
167 {
168 //calculate diffusivity for i,numComponents
169 const auto xi = moleFrac[compIIdx];
170 const auto Mn = FluidSystem::molarMass(numComponents-1);
171
172 auto tinInside = Deprecated::template effectiveMSDiffusionCoefficient<EffDiffModel, FluidSystem>(insideVolVars, phaseIdx, compIIdx, numComponents-1, problem, element, fvGeometry.scv(insideScvIdx));
173 auto tinOutside = Deprecated::template effectiveMSDiffusionCoefficient<EffDiffModel, FluidSystem>(outsideVolVars, phaseIdx, compIIdx, numComponents-1, problem, element, fvGeometry.scv(outsideScvIdx));
174
175 // scale by extrusion factor
176 tinInside *= insideVolVars.extrusionFactor();
177 tinOutside *= outsideVolVars.extrusionFactor();
178
179 // the resulting averaged diffusion tensor
180 const auto tin = problem.spatialParams().harmonicMean(tinInside, tinOutside, scvf.unitOuterNormal());
181
182 //begin the entrys of the diffusion matrix of the diagonal
183 reducedDiffusionMatrix[compIIdx][compIIdx] += xi*avgMolarMass/(tin*Mn);
184
185 // now set the rest of the entries (off-diagonal and additional entries for diagonal)
186 for (int compJIdx = 0; compJIdx < numComponents; compJIdx++)
187 {
188 //we don't want to calculate e.g. water in water diffusion
189 if (compIIdx == compJIdx)
190 continue;
191
192 //calculate diffusivity for compIIdx, compJIdx
193 const auto xj = moleFrac[compJIdx];
194 const auto Mi = FluidSystem::molarMass(compIIdx);
195 const auto Mj = FluidSystem::molarMass(compJIdx);
196
197 // effective diffusion tensors
198 auto tijInside = Deprecated::template effectiveMSDiffusionCoefficient<EffDiffModel, FluidSystem>(insideVolVars, phaseIdx, compIIdx, compJIdx, problem, element, fvGeometry.scv(insideScvIdx));
199
200 auto tijOutside = Deprecated::template effectiveMSDiffusionCoefficient<EffDiffModel, FluidSystem>(outsideVolVars, phaseIdx, compIIdx, compJIdx, problem, element, fvGeometry.scv(outsideScvIdx));
201
202 // scale by extrusion factor
203 tijInside *= insideVolVars.extrusionFactor();
204 tijOutside *= outsideVolVars.extrusionFactor();
205
206 // the resulting averaged diffusion tensor
207 const auto tij = problem.spatialParams().harmonicMean(tijInside, tijOutside, scvf.unitOuterNormal());
208
209 reducedDiffusionMatrix[compIIdx][compIIdx] += xj*avgMolarMass/(tij*Mi);
210 if (compJIdx < numComponents-1)
211 reducedDiffusionMatrix[compIIdx][compJIdx] +=xi*(avgMolarMass/(tin*Mn) - avgMolarMass/(tij*Mj));
212 }
213 }
214 return reducedDiffusionMatrix;
215 }
216
217};
218} // end namespace Dumux
219
220#endif
Helpers for deprecation.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
Classes related to flux variables caching.
Container storing the diffusion coefficients required by the Maxwell- Stefan diffusion law....
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
ReferenceSystemFormulation
The formulations available for Fick's law related to the reference system.
Definition: referencesystemformulation.hh:45
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
Definition: maxwellstefanslaw.hh:37
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: box/maxwellstefanslaw.hh:77
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/maxwellstefanslaw.hh:80
Container storing the diffusion coefficients required by the Maxwell- Stefan diffusion law....
Definition: maxwellstefandiffusioncoefficients.hh:45
Declares all properties used in Dumux.