3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
cctpfa/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_CC_TPFA_MAXWELL_STEFAN_LAW_HH
26#define DUMUX_DISCRETIZATION_CC_TPFA_MAXWELL_STEFAN_LAW_HH
27
28#include <dune/common/float_cmp.hh>
29
30#include <dumux/common/math.hh>
36
37namespace Dumux {
38
39// forward declaration
40template <class TypeTag, DiscretizationMethod discMethod, ReferenceSystemFormulation referenceSystem>
41class MaxwellStefansLawImplementation;
42
47template <class TypeTag, ReferenceSystemFormulation referenceSystem>
49{
53 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
54 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
55 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
58 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
59 using Element = typename GridView::template Codim<0>::Entity;
60 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
63
64 static const int dim = GridView::dimension;
65 static const int dimWorld = GridView::dimensionworld;
68
69 static_assert(referenceSystem == ReferenceSystemFormulation::massAveraged, "only the mass averaged reference system is supported for the Maxwell-Stefan formulation");
70
71
72 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
73 using ReducedComponentVector = Dune::FieldVector<Scalar, numComponents-1>;
74 using ReducedComponentMatrix = Dune::FieldMatrix<Scalar, numComponents-1, numComponents-1>;
75
76public:
77 // state the discretization method this implementation belongs to
79 //return the reference system
81 { return referenceSystem; }
82
87
88 static ComponentFluxVector flux(const Problem& problem,
89 const Element& element,
90 const FVElementGeometry& fvGeometry,
91 const ElementVolumeVariables& elemVolVars,
92 const SubControlVolumeFace& scvf,
93 const int phaseIdx,
94 const ElementFluxVariablesCache& elemFluxVarsCache)
95 {
96 //this is to calculate the maxwellStefan diffusion in a multicomponent system.
97 //see: Multicomponent Mass Transfer. R. Taylor u. R. Krishna. J. Wiley & Sons, New York 1993
98 ComponentFluxVector componentFlux(0.0);
99 ReducedComponentVector moleFracInside(0.0);
100 ReducedComponentVector moleFracOutside(0.0);
101 ReducedComponentVector reducedFlux(0.0);
102 ReducedComponentMatrix reducedDiffusionMatrixInside(0.0);
103 ReducedComponentMatrix reducedDiffusionMatrixOutside(0.0);
104
105 // get inside/outside volume variables
106 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
107 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
108 const auto rhoInside = insideVolVars.density(phaseIdx);
109 const auto rhoOutside = outsideVolVars.density(phaseIdx);
110 //calculate the mole fraction vectors
111 for (int compIdx = 0; compIdx < numComponents-1; compIdx++)
112 {
113 //calculate x_inside
114 const auto xInside = insideVolVars.moleFraction(phaseIdx, compIdx);
115 //calculate outside molefraction with the respective transmissibility
116 const auto xOutside = outsideVolVars.moleFraction(phaseIdx, compIdx);
117
118 moleFracInside[compIdx] = xInside;
119 moleFracOutside[compIdx] = xOutside;
120 }
121
122 //we cannot solve that if the matrix is 0 everywhere
123 if(!(Dune::FloatCmp::eq<Scalar>(insideVolVars.saturation(phaseIdx), 0) || Dune::FloatCmp::eq<Scalar>(outsideVolVars.saturation(phaseIdx), 0)))
124 {
125 const auto insideScvIdx = scvf.insideScvIdx();
126 const auto& insideScv = fvGeometry.scv(insideScvIdx);
127 const Scalar omegai = calculateOmega_(scvf,
128 insideScv,
129 insideVolVars.extrusionFactor());
130
131 //now we have to do the tpfa: J_i = J_j which leads to: tij(xi -xj) = -rho Bi^-1 omegai(x*-xi) with x* = (omegai Bi^-1 + omegaj Bj^-1)^-1 (xi omegai Bi^-1 + xj omegaj Bj^-1) with i inside and j outside
132 reducedDiffusionMatrixInside = setupMSMatrix_(problem, element, fvGeometry, insideVolVars, insideScv, phaseIdx);
133
134 //if on boundary
135 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
136 {
137 moleFracOutside -= moleFracInside;
138 reducedDiffusionMatrixInside.solve(reducedFlux, moleFracOutside);
139 reducedFlux *= omegai*rhoInside;
140 }
141 else //we need outside cells as well if we are not on the boundary
142 {
143 Scalar omegaj;
144 const auto outsideScvIdx = scvf.outsideScvIdx();
145 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
146
147 reducedDiffusionMatrixOutside = setupMSMatrix_(problem, element, fvGeometry, outsideVolVars, outsideScv, phaseIdx);
148
149 if (dim == dimWorld)
150 // assume the normal vector from outside is anti parallel so we save flipping a vector
151 omegaj = -1.0*calculateOmega_(scvf,
152 outsideScv,
153 outsideVolVars.extrusionFactor());
154 else
155 omegaj = calculateOmega_(fvGeometry.flipScvf(scvf.index()),
156 outsideScv,
157 outsideVolVars.extrusionFactor());
158
159 reducedDiffusionMatrixInside.invert();
160 reducedDiffusionMatrixOutside.invert();
161 reducedDiffusionMatrixInside *= omegai*rhoInside;
162 reducedDiffusionMatrixOutside *= omegaj*rhoOutside;
163
164 //in the helpervector we store the values for x*
165 ReducedComponentVector helperVector(0.0);
166 ReducedComponentVector gradientVectori(0.0);
167 ReducedComponentVector gradientVectorj(0.0);
168
169 reducedDiffusionMatrixInside.mv(moleFracInside, gradientVectori);
170 reducedDiffusionMatrixOutside.mv(moleFracOutside, gradientVectorj);
171
172 auto gradientVectorij = (gradientVectori + gradientVectorj);
173
174 //add the two matrixes to each other
175 reducedDiffusionMatrixOutside += reducedDiffusionMatrixInside;
176
177 reducedDiffusionMatrixOutside.solve(helperVector, gradientVectorij);
178
179 //Bi^-1 omegai rho(x*-xi)
180 helperVector -=moleFracInside;
181 reducedDiffusionMatrixInside.mv(helperVector, reducedFlux);
182 }
183
184 reducedFlux *= -scvf.area();
185 for (int compIdx = 0; compIdx < numComponents-1; compIdx++)
186 {
187 componentFlux[compIdx] = reducedFlux[compIdx];
188 componentFlux[numComponents-1] -=reducedFlux[compIdx];
189 }
190 }
191 return componentFlux ;
192 }
193
194private:
195 static Scalar calculateOmega_(const SubControlVolumeFace& scvf,
196 const SubControlVolume &scv,
197 const Scalar extrusionFactor)
198 {
199 auto distanceVector = scvf.ipGlobal();
200 distanceVector -= scv.center();
201 distanceVector /= distanceVector.two_norm2();
202
203 Scalar omega = (distanceVector * scvf.unitOuterNormal());
204 omega *= extrusionFactor;
205
206 return omega;
207 }
208
209 static ReducedComponentMatrix setupMSMatrix_(const Problem& problem,
210 const Element& element,
211 const FVElementGeometry& fvGeometry,
212 const VolumeVariables& volVars,
213 const SubControlVolume& scv,
214 const int phaseIdx)
215 {
216 ReducedComponentMatrix reducedDiffusionMatrix(0.0);
217
218 //this is to not devide by 0 if the saturation in 0 and the effectiveDiffusivity becomes zero due to that
219 if(Dune::FloatCmp::eq<Scalar>(volVars.saturation(phaseIdx), 0))
220 return reducedDiffusionMatrix;
221
222 const auto avgMolarMass = volVars.averageMolarMass(phaseIdx);
223
224 for (int compIIdx = 0; compIIdx < numComponents-1; compIIdx++)
225 {
226 const auto xi = volVars.moleFraction(phaseIdx, compIIdx);
227
228 const auto Mn = FluidSystem::molarMass(numComponents-1);
229 Scalar tin = getDiffusionCoefficient(phaseIdx, compIIdx, numComponents-1, problem, element, volVars, scv);
230 tin = EffDiffModel::effectiveDiffusivity(volVars.porosity(), volVars.saturation(phaseIdx), tin);
231
232 // set the entries of the diffusion matrix of the diagonal
233 reducedDiffusionMatrix[compIIdx][compIIdx] += xi*avgMolarMass/(tin*Mn);
234
235 // now set the rest of the entries (off-diagonal and additional entries for diagonal)
236 for (int compJIdx = 0; compJIdx < numComponents; compJIdx++)
237 {
238 // we don't want to calculate e.g. water in water diffusion
239 if (compJIdx == compIIdx)
240 continue;
241
242 const auto xj = volVars.moleFraction(phaseIdx, compJIdx);
243 const auto Mi = FluidSystem::molarMass(compIIdx);
244 const auto Mj = FluidSystem::molarMass(compJIdx);
245 Scalar tij = getDiffusionCoefficient(phaseIdx, compIIdx, compJIdx, problem, element, volVars, scv);
246 tij = EffDiffModel::effectiveDiffusivity(volVars.porosity(), volVars.saturation(phaseIdx), tij);
247 reducedDiffusionMatrix[compIIdx][compIIdx] += xj*avgMolarMass/(tij*Mi);
248 if (compJIdx < numComponents-1)
249 reducedDiffusionMatrix[compIIdx][compJIdx] += xi*(avgMolarMass/(tin*Mn) - avgMolarMass/(tij*Mj));
250 }
251 }
252 return reducedDiffusionMatrix;
253 }
254
255 template <class T = TypeTag, typename std::enable_if_t<GetPropType<T, Properties::FluidSystem>::isTracerFluidSystem(), int> =0 >
256 static Scalar getDiffusionCoefficient(const int phaseIdx,
257 const int compIIdx,
258 const int compJIdx,
259 const Problem& problem,
260 const Element& element,
261 const VolumeVariables& volVars,
262 const SubControlVolume& scv)
263 {
264 return FluidSystem::binaryDiffusionCoefficient(compIIdx,
265 compJIdx,
266 problem,
267 element,
268 scv);
269 }
270
271 template <class T = TypeTag, typename std::enable_if_t<!GetPropType<T, Properties::FluidSystem>::isTracerFluidSystem(), int> =0 >
272 static Scalar getDiffusionCoefficient(const int phaseIdx,
273 const int compIIdx,
274 const int compJIdx,
275 const Problem& problem,
276 const Element& element,
277 const VolumeVariables& volVars,
278 const SubControlVolume& scv)
279 {
280 auto fluidState = volVars.fluidState();
281 typename FluidSystem::ParameterCache paramCache;
282 paramCache.updateAll(fluidState);
283 return FluidSystem::binaryDiffusionCoefficient(fluidState,
284 paramCache,
285 phaseIdx,
286 compIIdx,
287 compJIdx);
288 }
289
290
291
292};
293} // end namespace
294
295#endif
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
Classes related to flux variables caching.
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
ReferenceSystemFormulation
The formulations available for Fick's law related to the reference system.
Definition: referencesystemformulation.hh:45
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
Definition: maxwellstefanslaw.hh:37
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: cctpfa/maxwellstefanslaw.hh:80
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: cctpfa/maxwellstefanslaw.hh:88
The empty filler class corresponding to EmptyCache.
Definition: fluxvariablescaching.hh:32
Definition: fluxvariablescaching.hh:64
Declares all properties used in Dumux.