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