3.2-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/fmatrix.hh>
29#include <dune/common/float_cmp.hh>
30
34
39
40namespace Dumux {
41
42// forward declaration
43template <class TypeTag, DiscretizationMethod discMethod, ReferenceSystemFormulation referenceSystem>
44class MaxwellStefansLawImplementation;
45
50template <class TypeTag, ReferenceSystemFormulation referenceSystem>
52{
56 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
57 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
58 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
61 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
62 using Element = typename GridView::template Codim<0>::Entity;
63 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
66
67 static const int dim = GridView::dimension;
68 static const int dimWorld = GridView::dimensionworld;
71
72 static_assert(referenceSystem == ReferenceSystemFormulation::massAveraged, "only the mass averaged reference system is supported for the Maxwell-Stefan formulation");
73
74
75 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
76 using ReducedComponentVector = Dune::FieldVector<Scalar, numComponents-1>;
77 using ReducedComponentMatrix = Dune::FieldMatrix<Scalar, numComponents-1, numComponents-1>;
78
79public:
80 // state the discretization method this implementation belongs to
82 //return the reference system
84 { return referenceSystem; }
85
87
92
93 static ComponentFluxVector flux(const Problem& problem,
94 const Element& element,
95 const FVElementGeometry& fvGeometry,
96 const ElementVolumeVariables& elemVolVars,
97 const SubControlVolumeFace& scvf,
98 const int phaseIdx,
99 const ElementFluxVariablesCache& elemFluxVarsCache)
100 {
101 //this is to calculate the maxwellStefan diffusion in a multicomponent system.
102 //see: Multicomponent Mass Transfer. R. Taylor u. R. Krishna. J. Wiley & Sons, New York 1993
103 ComponentFluxVector componentFlux(0.0);
104 ReducedComponentVector moleFracInside(0.0);
105 ReducedComponentVector moleFracOutside(0.0);
106 ReducedComponentVector reducedFlux(0.0);
107 ReducedComponentMatrix reducedDiffusionMatrixInside(0.0);
108 ReducedComponentMatrix reducedDiffusionMatrixOutside(0.0);
109
110 // get inside/outside volume variables
111 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
112 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
113 const auto rhoInside = insideVolVars.density(phaseIdx);
114 const auto rhoOutside = outsideVolVars.density(phaseIdx);
115 //calculate the mole fraction vectors
116 for (int compIdx = 0; compIdx < numComponents-1; compIdx++)
117 {
118 //calculate x_inside
119 const auto xInside = insideVolVars.moleFraction(phaseIdx, compIdx);
120 //calculate outside molefraction with the respective transmissibility
121 const auto xOutside = outsideVolVars.moleFraction(phaseIdx, compIdx);
122
123 moleFracInside[compIdx] = xInside;
124 moleFracOutside[compIdx] = xOutside;
125 }
126
127 //we cannot solve that if the matrix is 0 everywhere
128 if(!(Dune::FloatCmp::eq<Scalar>(insideVolVars.saturation(phaseIdx), 0) || Dune::FloatCmp::eq<Scalar>(outsideVolVars.saturation(phaseIdx), 0)))
129 {
130 const auto insideScvIdx = scvf.insideScvIdx();
131 const auto& insideScv = fvGeometry.scv(insideScvIdx);
132 const Scalar omegai = calculateOmega_(scvf,
133 insideScv,
134 insideVolVars.extrusionFactor());
135
136 //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
137 reducedDiffusionMatrixInside = setupMSMatrix_(problem, element, fvGeometry, insideVolVars, insideScv, phaseIdx);
138
139 //if on boundary
140 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
141 {
142 moleFracOutside -= moleFracInside;
143 reducedDiffusionMatrixInside.solve(reducedFlux, moleFracOutside);
144 reducedFlux *= omegai*rhoInside;
145 }
146 else //we need outside cells as well if we are not on the boundary
147 {
148 Scalar omegaj;
149 const auto outsideScvIdx = scvf.outsideScvIdx();
150 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
151
152 reducedDiffusionMatrixOutside = setupMSMatrix_(problem, element, fvGeometry, outsideVolVars, outsideScv, phaseIdx);
153
154 if (dim == dimWorld)
155 // assume the normal vector from outside is anti parallel so we save flipping a vector
156 omegaj = -1.0*calculateOmega_(scvf,
157 outsideScv,
158 outsideVolVars.extrusionFactor());
159 else
160 omegaj = calculateOmega_(fvGeometry.flipScvf(scvf.index()),
161 outsideScv,
162 outsideVolVars.extrusionFactor());
163
164 reducedDiffusionMatrixInside.invert();
165 reducedDiffusionMatrixOutside.invert();
166 reducedDiffusionMatrixInside *= omegai*rhoInside;
167 reducedDiffusionMatrixOutside *= omegaj*rhoOutside;
168
169 //in the helpervector we store the values for x*
170 ReducedComponentVector helperVector(0.0);
171 ReducedComponentVector gradientVectori(0.0);
172 ReducedComponentVector gradientVectorj(0.0);
173
174 reducedDiffusionMatrixInside.mv(moleFracInside, gradientVectori);
175 reducedDiffusionMatrixOutside.mv(moleFracOutside, gradientVectorj);
176
177 auto gradientVectorij = (gradientVectori + gradientVectorj);
178
179 //add the two matrixes to each other
180 reducedDiffusionMatrixOutside += reducedDiffusionMatrixInside;
181
182 reducedDiffusionMatrixOutside.solve(helperVector, gradientVectorij);
183
184 //Bi^-1 omegai rho(x*-xi)
185 helperVector -=moleFracInside;
186 reducedDiffusionMatrixInside.mv(helperVector, reducedFlux);
187 }
188
189 reducedFlux *= -scvf.area();
190 for (int compIdx = 0; compIdx < numComponents-1; compIdx++)
191 {
192 componentFlux[compIdx] = reducedFlux[compIdx];
193 componentFlux[numComponents-1] -=reducedFlux[compIdx];
194 }
195 }
196 return componentFlux ;
197 }
198
199private:
200 static Scalar calculateOmega_(const SubControlVolumeFace& scvf,
201 const SubControlVolume &scv,
202 const Scalar extrusionFactor)
203 {
204 auto distanceVector = scvf.ipGlobal();
205 distanceVector -= scv.center();
206 distanceVector /= distanceVector.two_norm2();
207
208 Scalar omega = (distanceVector * scvf.unitOuterNormal());
209 omega *= extrusionFactor;
210
211 return omega;
212 }
213
214 static ReducedComponentMatrix setupMSMatrix_(const Problem& problem,
215 const Element& element,
216 const FVElementGeometry& fvGeometry,
217 const VolumeVariables& volVars,
218 const SubControlVolume& scv,
219 const int phaseIdx)
220 {
221 using EffDiffModel = GetPropType<TypeTag, Properties::EffectiveDiffusivityModel>;
222 ReducedComponentMatrix reducedDiffusionMatrix(0.0);
223
224 //this is to not devide by 0 if the saturation in 0 and the effectiveDiffusionCoefficient becomes zero due to that
225 if(Dune::FloatCmp::eq<Scalar>(volVars.saturation(phaseIdx), 0))
226 return reducedDiffusionMatrix;
227
228 const auto avgMolarMass = volVars.averageMolarMass(phaseIdx);
229
230 for (int compIIdx = 0; compIIdx < numComponents-1; compIIdx++)
231 {
232 const auto xi = volVars.moleFraction(phaseIdx, compIIdx);
233
234 const auto Mn = FluidSystem::molarMass(numComponents-1);
235 Scalar tin = Deprecated::template effectiveMSDiffusionCoefficient<EffDiffModel, FluidSystem>(volVars, phaseIdx, compIIdx, numComponents-1, problem, element, scv);
236
237 // set the entries of the diffusion matrix of the diagonal
238 reducedDiffusionMatrix[compIIdx][compIIdx] += xi*avgMolarMass/(tin*Mn);
239
240 // now set the rest of the entries (off-diagonal and additional entries for diagonal)
241 for (int compJIdx = 0; compJIdx < numComponents; compJIdx++)
242 {
243 // we don't want to calculate e.g. water in water diffusion
244 if (compJIdx == compIIdx)
245 continue;
246
247 const auto xj = volVars.moleFraction(phaseIdx, compJIdx);
248 const auto Mi = FluidSystem::molarMass(compIIdx);
249 const auto Mj = FluidSystem::molarMass(compJIdx);
250 Scalar tij = Deprecated::template effectiveMSDiffusionCoefficient<EffDiffModel, FluidSystem>(volVars, phaseIdx, compIIdx, compJIdx, problem, element, scv);
251 reducedDiffusionMatrix[compIIdx][compIIdx] += xj*avgMolarMass/(tij*Mi);
252
253 if (compJIdx < numComponents-1)
254 reducedDiffusionMatrix[compIIdx][compJIdx] += xi*(avgMolarMass/(tin*Mn) - avgMolarMass/(tij*Mj));
255 }
256 }
257 return reducedDiffusionMatrix;
258 }
259
260};
261} // end namespace
262
263#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: cctpfa/maxwellstefanslaw.hh:83
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:93
The empty filler class corresponding to EmptyCache.
Definition: fluxvariablescaching.hh:32
Definition: fluxvariablescaching.hh:64
Container storing the diffusion coefficients required by the Maxwell- Stefan diffusion law....
Definition: maxwellstefandiffusioncoefficients.hh:45
Declares all properties used in Dumux.