3.1-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/float_cmp.hh>
28
29#include <dumux/common/math.hh>
35
36namespace Dumux {
37
38// forward declaration
39template <class TypeTag, DiscretizationMethod discMethod, ReferenceSystemFormulation referenceSystem>
40class MaxwellStefansLawImplementation;
41
46template <class TypeTag, ReferenceSystemFormulation referenceSystem>
48{
52 using FVElementGeometry = typename GridGeometry::LocalView;
53 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
54 using GridView = typename GridGeometry::GridView;
55 using Element = typename GridView::template Codim<0>::Entity;
56 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
58 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
61
63
64 static const int numComponents = ModelTraits::numFluidComponents();
65 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
66
67 using ReducedComponentVector = Dune::FieldVector<Scalar, numComponents-1>;
68 using ReducedComponentMatrix = Dune::FieldMatrix<Scalar, numComponents-1, numComponents-1>;
69
70 static_assert(ModelTraits::numFluidPhases() == 1, "Only one phase allowed supported!");
71
72 static_assert(referenceSystem == ReferenceSystemFormulation::massAveraged, "only the mass averaged reference system is supported for the Maxwell-Stefan formulation");
73
74public:
75 // state the discretization method this implementation belongs to
77 //return the reference system
79 { return referenceSystem; }
80
85
86 static CellCenterPrimaryVariables flux(const Problem& problem,
87 const Element& element,
88 const FVElementGeometry& fvGeometry,
89 const ElementVolumeVariables& elemVolVars,
90 const SubControlVolumeFace& scvf)
91 {
92 //this is to calculate the maxwellStefan diffusion in a multicomponent system.
93 //see: Multicomponent Mass Transfer. R. Taylor u. R. Krishna. J. Wiley & Sons, New York 1993
94 CellCenterPrimaryVariables componentFlux(0.0);
95 ReducedComponentVector moleFracInside(0.0);
96 ReducedComponentVector moleFracOutside(0.0);
97 ReducedComponentVector reducedFlux(0.0);
98 ReducedComponentMatrix reducedDiffusionMatrixInside(0.0);
99 ReducedComponentMatrix reducedDiffusionMatrixOutside(0.0);
100
101 // get inside/outside volume variables
102 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
103 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
104 const auto rhoInside = insideVolVars.density();
105 const auto rhoOutside = outsideVolVars.density();
106
107 //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.
108 for(int compIdx = 0; compIdx < numComponents; ++compIdx)
109 {
110 if(compIdx == FluidSystem::getMainComponent(0))
111 continue;
112
113 // get equation index
114 const auto eqIdx = Indices::conti0EqIdx + compIdx;
115 if(scvf.boundary())
116 {
117 const auto bcTypes = problem.boundaryTypes(element, scvf);
118 if((bcTypes.isOutflow(eqIdx)) || (bcTypes.isSymmetry()))
119 return componentFlux;
120 }
121 }
122
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(compIdx);
128 //calculate outside molefraction with the respective transmissibility
129 const auto xOutside = outsideVolVars.moleFraction(compIdx);
130
131 moleFracInside[compIdx] = xInside;
132 moleFracOutside[compIdx] = xOutside;
133 }
134
135 //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
136 reducedDiffusionMatrixInside = setupMSMatrix_(problem, fvGeometry, insideVolVars, scvf);
137
138 reducedDiffusionMatrixOutside = setupMSMatrix_(problem, fvGeometry, outsideVolVars, scvf);
139
140 const auto insideScvIdx = scvf.insideScvIdx();
141 const auto& insideScv = fvGeometry.scv(insideScvIdx);
142 const auto outsideScvIdx = scvf.outsideScvIdx();
143 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
144
145 const Scalar insideDistance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
146
147 const Scalar omegai = calculateOmega_(insideDistance, insideVolVars.extrusionFactor());
148
149 //if on boundary
150 if (scvf.boundary())
151 {
152 moleFracOutside -= moleFracInside;
153 reducedDiffusionMatrixInside.solve(reducedFlux, moleFracOutside);
154 reducedFlux *= omegai*rhoInside;
155 }
156 else
157 {
158 const Scalar outsideDistance = (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm();
159 const Scalar omegaj = calculateOmega_(outsideDistance, outsideVolVars.extrusionFactor());
160
161 reducedDiffusionMatrixInside.invert();
162 reducedDiffusionMatrixOutside.invert();
163 reducedDiffusionMatrixInside *= omegai*rhoInside;
164 reducedDiffusionMatrixOutside *= omegaj*rhoOutside;
165
166 //in the helpervector we store the values for x*
167 ReducedComponentVector helperVector(0.0);
168 ReducedComponentVector gradientVectori(0.0);
169 ReducedComponentVector gradientVectorj(0.0);
170
171 reducedDiffusionMatrixInside.mv(moleFracInside, gradientVectori);
172 reducedDiffusionMatrixOutside.mv(moleFracOutside, gradientVectorj);
173
174 auto gradientVectorij = (gradientVectori + gradientVectorj);
175
176 //add the two matrixes to each other
177 reducedDiffusionMatrixOutside += reducedDiffusionMatrixInside;
178
179 reducedDiffusionMatrixOutside.solve(helperVector, gradientVectorij);
180
181 //Bi^-1 omegai rho(x*-xi)
182 helperVector -=moleFracInside;
183 reducedDiffusionMatrixInside.mv(helperVector, reducedFlux);
184 }
185
186 reducedFlux *= -scvf.area();
187
188 for (int compIdx = 0; compIdx < numComponents-1; compIdx++)
189 {
190 componentFlux[compIdx] = reducedFlux[compIdx];
191 componentFlux[numComponents-1] -= reducedFlux[compIdx];
192 }
193
194 return componentFlux ;
195 }
196
197private:
198 static Scalar calculateOmega_(const Scalar distance,
199 const Scalar extrusionFactor)
200 {
201 Scalar omega = 1/distance;
202 omega *= extrusionFactor;
203
204 return omega;
205 }
206
207 static ReducedComponentMatrix setupMSMatrix_(const Problem& problem,
208 const FVElementGeometry& fvGeometry,
209 const VolumeVariables& volVars,
210 const SubControlVolumeFace& scvf)
211 {
212 ReducedComponentMatrix reducedDiffusionMatrix(0.0);
213
214 for (int compIIdx = 0; compIIdx < numComponents-1; compIIdx++)
215 {
216 const auto xi = volVars.moleFraction(compIIdx);
217 const auto avgMolarMass = volVars.averageMolarMass(0);
218 const auto Mn = FluidSystem::molarMass(numComponents-1);
219 const Scalar tin = volVars.effectiveDiffusivity(compIIdx, numComponents-1);
220
221 // set the entries of the diffusion matrix of the diagonal
222 reducedDiffusionMatrix[compIIdx][compIIdx] += xi*avgMolarMass/(tin*Mn);
223
224 for (int compJIdx = 0; compJIdx < numComponents; compJIdx++)
225 {
226 // we don't want to calculate e.g. water in water diffusion
227 if (compJIdx == compIIdx)
228 continue;
229
230 const auto xj = volVars.moleFraction(compJIdx);
231 const auto Mi = FluidSystem::molarMass(compIIdx);
232 const auto Mj = FluidSystem::molarMass(compJIdx);
233 const Scalar tij = volVars.effectiveDiffusivity(compIIdx, compJIdx);
234 reducedDiffusionMatrix[compIIdx][compIIdx] += xj*avgMolarMass/(tij*Mi);
235 if (compJIdx < numComponents-1)
236 reducedDiffusionMatrix[compIIdx][compJIdx] += xi*(avgMolarMass/(tin*Mn) - avgMolarMass/(tij*Mj));
237 }
238 }
239 return reducedDiffusionMatrix;
240 }
241};
242} // end namespace
243
244#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
The empty filler class corresponding to EmptyCache.
Definition: fluxvariablescaching.hh:32
Definition: fluxvariablescaching.hh:64
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: staggered/freeflow/maxwellstefanslaw.hh:78
static CellCenterPrimaryVariables flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Definition: staggered/freeflow/maxwellstefanslaw.hh:86
Declares all properties used in Dumux.