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