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