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