3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
staggered/freeflow/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 *****************************************************************************/
24#ifndef DUMUX_DISCRETIZATION_STAGGERED_MAXWELL_STEFAN_LAW_HH
25#define DUMUX_DISCRETIZATION_STAGGERED_MAXWELL_STEFAN_LAW_HH
26
27#include <dune/common/fmatrix.hh>
28
33
37
38namespace Dumux {
39
40// forward declaration
41template <class TypeTag, class DiscretizationMethod, ReferenceSystemFormulation referenceSystem>
42class MaxwellStefansLawImplementation;
43
48template <class TypeTag, ReferenceSystemFormulation referenceSystem>
49class MaxwellStefansLawImplementation<TypeTag, DiscretizationMethods::Staggered, referenceSystem>
50{
54 using FVElementGeometry = typename GridGeometry::LocalView;
55 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
56 using Extrusion = Extrusion_t<GridGeometry>;
57 using GridView = typename GridGeometry::GridView;
58 using Element = typename GridView::template Codim<0>::Entity;
60 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
63
65
66 static const int numComponents = ModelTraits::numFluidComponents();
67 static const int numPhases = ModelTraits::numFluidPhases();
68 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
69
70 using ReducedComponentVector = Dune::FieldVector<Scalar, numComponents-1>;
71 using ReducedComponentMatrix = Dune::FieldMatrix<Scalar, numComponents-1, numComponents-1>;
72
73 static_assert(ModelTraits::numFluidPhases() == 1, "Only one phase allowed supported!");
74
75 static_assert(referenceSystem == ReferenceSystemFormulation::massAveraged, "only the mass averaged reference system is supported for the Maxwell-Stefan formulation");
76
77public:
79 // state the discretization method this implementation belongs to
80 static constexpr DiscretizationMethod discMethod{};
81
82 //return the reference system
84 { return referenceSystem; }
85
90
92
98 template<class ElementVolumeVariables>
99 static CellCenterPrimaryVariables flux(const Problem& problem,
100 const Element& element,
101 const FVElementGeometry& fvGeometry,
102 const ElementVolumeVariables& elemVolVars,
103 const SubControlVolumeFace& scvf)
104 {
105 //this is to calculate the maxwellStefan diffusion in a multicomponent system.
106 //see: Multicomponent Mass Transfer. R. Taylor u. R. Krishna. J. Wiley & Sons, New York 1993
107 CellCenterPrimaryVariables componentFlux(0.0);
108 ReducedComponentVector moleFracInside(0.0);
109 ReducedComponentVector moleFracOutside(0.0);
110 ReducedComponentVector reducedFlux(0.0);
111 ReducedComponentMatrix reducedDiffusionMatrixInside(0.0);
112 ReducedComponentMatrix reducedDiffusionMatrixOutside(0.0);
113
114 // get inside/outside volume variables
115 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
116 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
117 const auto rhoInside = insideVolVars.density();
118 const auto rhoOutside = outsideVolVars.density();
119
120 //to implement outflow boundaries correctly we need to loop over all components but the main component as only for the transported ones we implement the outflow boundary. diffusion then is 0.
121 for(int compIdx = 0; compIdx < numComponents; ++compIdx)
122 {
123 if(compIdx == FluidSystem::getMainComponent(0))
124 continue;
125
126 // get equation index
127 const auto eqIdx = Indices::conti0EqIdx + compIdx;
128 if(scvf.boundary())
129 {
130 const auto bcTypes = problem.boundaryTypes(element, scvf);
131 if((bcTypes.isOutflow(eqIdx)) || (bcTypes.isSymmetry()))
132 return componentFlux;
133 }
134 }
135
136 //calculate the mole fraction vectors
137 for (int compIdx = 0; compIdx < numComponents-1; compIdx++)
138 {
139 //calculate x_inside
140 const auto xInside = insideVolVars.moleFraction(compIdx);
141 //calculate outside molefraction with the respective transmissibility
142 const auto xOutside = outsideVolVars.moleFraction(compIdx);
143
144 moleFracInside[compIdx] = xInside;
145 moleFracOutside[compIdx] = xOutside;
146 }
147
148 //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
149 reducedDiffusionMatrixInside = setupMSMatrix_(problem, fvGeometry, insideVolVars, scvf);
150
151 reducedDiffusionMatrixOutside = setupMSMatrix_(problem, fvGeometry, outsideVolVars, scvf);
152
153 const auto insideScvIdx = scvf.insideScvIdx();
154 const auto& insideScv = fvGeometry.scv(insideScvIdx);
155 const auto outsideScvIdx = scvf.outsideScvIdx();
156 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
157
158 const Scalar insideDistance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
159
160 const Scalar omegai = calculateOmega_(insideDistance, insideVolVars.extrusionFactor());
161
162 //if on boundary
163 if (scvf.boundary())
164 {
165 moleFracOutside -= moleFracInside;
166 reducedDiffusionMatrixInside.solve(reducedFlux, moleFracOutside);
167 reducedFlux *= omegai*rhoInside;
168 }
169 else
170 {
171 const Scalar outsideDistance = (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm();
172 const Scalar omegaj = calculateOmega_(outsideDistance, outsideVolVars.extrusionFactor());
173
174 reducedDiffusionMatrixInside.invert();
175 reducedDiffusionMatrixOutside.invert();
176 reducedDiffusionMatrixInside *= omegai*rhoInside;
177 reducedDiffusionMatrixOutside *= omegaj*rhoOutside;
178
179 //in the helpervector we store the values for x*
180 ReducedComponentVector helperVector(0.0);
181 ReducedComponentVector gradientVectori(0.0);
182 ReducedComponentVector gradientVectorj(0.0);
183
184 reducedDiffusionMatrixInside.mv(moleFracInside, gradientVectori);
185 reducedDiffusionMatrixOutside.mv(moleFracOutside, gradientVectorj);
186
187 auto gradientVectorij = (gradientVectori + gradientVectorj);
188
189 //add the two matrixes to each other
190 reducedDiffusionMatrixOutside += reducedDiffusionMatrixInside;
191
192 reducedDiffusionMatrixOutside.solve(helperVector, gradientVectorij);
193
194 //Bi^-1 omegai rho(x*-xi)
195 helperVector -=moleFracInside;
196 reducedDiffusionMatrixInside.mv(helperVector, reducedFlux);
197 }
198
199 reducedFlux *= -Extrusion::area(fvGeometry, scvf);
200
201 for (int compIdx = 0; compIdx < numComponents-1; compIdx++)
202 {
203 componentFlux[compIdx] = reducedFlux[compIdx];
204 componentFlux[numComponents-1] -= reducedFlux[compIdx];
205 }
206
207 return componentFlux ;
208 }
209
210private:
211 static Scalar calculateOmega_(const Scalar distance,
212 const Scalar extrusionFactor)
213 {
214 Scalar omega = 1/distance;
215 omega *= extrusionFactor;
216
217 return omega;
218 }
219
220 static ReducedComponentMatrix setupMSMatrix_(const Problem& problem,
221 const FVElementGeometry& fvGeometry,
222 const VolumeVariables& volVars,
223 const SubControlVolumeFace& scvf)
224 {
225 ReducedComponentMatrix reducedDiffusionMatrix(0.0);
226
227 for (int compIIdx = 0; compIIdx < numComponents-1; compIIdx++)
228 {
229 const auto xi = volVars.moleFraction(compIIdx);
230 const auto avgMolarMass = volVars.averageMolarMass(0);
231 const auto Mn = FluidSystem::molarMass(numComponents-1);
232 const Scalar tin = getEffectiveDiffusionCoefficient_(volVars, compIIdx, numComponents-1);
233
234 // set the entries of the diffusion matrix of the diagonal
235 reducedDiffusionMatrix[compIIdx][compIIdx] += xi*avgMolarMass/(tin*Mn);
236
237 for (int compJIdx = 0; compJIdx < numComponents; compJIdx++)
238 {
239 // we don't want to calculate e.g. water in water diffusion
240 if (compJIdx == compIIdx)
241 continue;
242
243 const auto xj = volVars.moleFraction(compJIdx);
244 const auto Mi = FluidSystem::molarMass(compIIdx);
245 const auto Mj = FluidSystem::molarMass(compJIdx);
246 const Scalar tij = getEffectiveDiffusionCoefficient_(volVars, compIIdx, compJIdx);
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 static Scalar getEffectiveDiffusionCoefficient_(const VolumeVariables& volVars, const int phaseIdx, const int compIdx)
256 {
257 return volVars.effectiveDiffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
258 }
259};
260} // end namespace
261
262#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.
static ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:294
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: box/maxwellstefanslaw.hh:45
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
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: staggered/freeflow/maxwellstefanslaw.hh:83
static CellCenterPrimaryVariables flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Returns the diffusive fluxes of all components within a fluid phase across the given sub-control volu...
Definition: staggered/freeflow/maxwellstefanslaw.hh:99
Declares all properties used in Dumux.