version 3.8
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_DISCRETIZATION_STAGGERED_MAXWELL_STEFAN_LAW_HH
13#define DUMUX_DISCRETIZATION_STAGGERED_MAXWELL_STEFAN_LAW_HH
14
15#include <dune/common/fmatrix.hh>
16
21
25
26namespace Dumux {
27
28// forward declaration
29template <class TypeTag, class DiscretizationMethod, ReferenceSystemFormulation referenceSystem>
30class MaxwellStefansLawImplementation;
31
36template <class TypeTag, ReferenceSystemFormulation referenceSystem>
37class MaxwellStefansLawImplementation<TypeTag, DiscretizationMethods::Staggered, referenceSystem>
38{
42 using FVElementGeometry = typename GridGeometry::LocalView;
43 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
44 using Extrusion = Extrusion_t<GridGeometry>;
45 using GridView = typename GridGeometry::GridView;
46 using Element = typename GridView::template Codim<0>::Entity;
48 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
51
53
54 static const int numComponents = ModelTraits::numFluidComponents();
55 static const int numPhases = ModelTraits::numFluidPhases();
56 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
57
58 using ReducedComponentVector = Dune::FieldVector<Scalar, numComponents-1>;
59 using ReducedComponentMatrix = Dune::FieldMatrix<Scalar, numComponents-1, numComponents-1>;
60
61 static_assert(ModelTraits::numFluidPhases() == 1, "Only one phase allowed supported!");
62
63 static_assert(referenceSystem == ReferenceSystemFormulation::massAveraged, "only the mass averaged reference system is supported for the Maxwell-Stefan formulation");
64
65public:
67 // state the discretization method this implementation belongs to
68 static constexpr DiscretizationMethod discMethod{};
69
70 //return the reference system
72 { return referenceSystem; }
73
78
80
86 template<class ElementVolumeVariables>
87 static CellCenterPrimaryVariables flux(const Problem& problem,
88 const Element& element,
89 const FVElementGeometry& fvGeometry,
90 const ElementVolumeVariables& elemVolVars,
91 const SubControlVolumeFace& scvf)
92 {
93 //this is to calculate the maxwellStefan diffusion in a multicomponent system.
94 //see: Multicomponent Mass Transfer. R. Taylor u. R. Krishna. J. Wiley & Sons, New York 1993
95 CellCenterPrimaryVariables componentFlux(0.0);
96 ReducedComponentVector moleFracInside(0.0);
97 ReducedComponentVector moleFracOutside(0.0);
98 ReducedComponentVector reducedFlux(0.0);
99 ReducedComponentMatrix reducedDiffusionMatrixInside(0.0);
100 ReducedComponentMatrix reducedDiffusionMatrixOutside(0.0);
101
102 // get inside/outside volume variables
103 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
104 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
105 const auto rhoInside = insideVolVars.density();
106 const auto rhoOutside = outsideVolVars.density();
107
108 //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.
109 for(int compIdx = 0; compIdx < numComponents; ++compIdx)
110 {
111 if(compIdx == FluidSystem::getMainComponent(0))
112 continue;
113
114 // get equation index
115 const auto eqIdx = Indices::conti0EqIdx + compIdx;
116 if(scvf.boundary())
117 {
118 const auto bcTypes = problem.boundaryTypes(element, scvf);
119 if((bcTypes.isOutflow(eqIdx)) || (bcTypes.isSymmetry()))
120 return componentFlux;
121 }
122 }
123
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(compIdx);
129 //calculate outside molefraction with the respective transmissibility
130 const auto xOutside = outsideVolVars.moleFraction(compIdx);
131
132 moleFracInside[compIdx] = xInside;
133 moleFracOutside[compIdx] = xOutside;
134 }
135
136 //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
137 reducedDiffusionMatrixInside = setupMSMatrix_(problem, fvGeometry, insideVolVars, scvf);
138
139 reducedDiffusionMatrixOutside = setupMSMatrix_(problem, fvGeometry, outsideVolVars, scvf);
140
141 const auto insideScvIdx = scvf.insideScvIdx();
142 const auto& insideScv = fvGeometry.scv(insideScvIdx);
143 const auto outsideScvIdx = scvf.outsideScvIdx();
144 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
145
146 const Scalar insideDistance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
147
148 const Scalar omegai = calculateOmega_(insideDistance, insideVolVars.extrusionFactor());
149
150 //if on boundary
151 if (scvf.boundary())
152 {
153 moleFracOutside -= moleFracInside;
154 reducedDiffusionMatrixInside.solve(reducedFlux, moleFracOutside);
155 reducedFlux *= omegai*rhoInside;
156 }
157 else
158 {
159 const Scalar outsideDistance = (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm();
160 const Scalar omegaj = calculateOmega_(outsideDistance, outsideVolVars.extrusionFactor());
161
162 reducedDiffusionMatrixInside.invert();
163 reducedDiffusionMatrixOutside.invert();
164 reducedDiffusionMatrixInside *= omegai*rhoInside;
165 reducedDiffusionMatrixOutside *= omegaj*rhoOutside;
166
167 //in the helpervector we store the values for x*
168 ReducedComponentVector helperVector(0.0);
169 ReducedComponentVector gradientVectori(0.0);
170 ReducedComponentVector gradientVectorj(0.0);
171
172 reducedDiffusionMatrixInside.mv(moleFracInside, gradientVectori);
173 reducedDiffusionMatrixOutside.mv(moleFracOutside, gradientVectorj);
174
175 auto gradientVectorij = (gradientVectori + gradientVectorj);
176
177 //add the two matrixes to each other
178 reducedDiffusionMatrixOutside += reducedDiffusionMatrixInside;
179
180 reducedDiffusionMatrixOutside.solve(helperVector, gradientVectorij);
181
182 //Bi^-1 omegai rho(x*-xi)
183 helperVector -=moleFracInside;
184 reducedDiffusionMatrixInside.mv(helperVector, reducedFlux);
185 }
186
187 reducedFlux *= -Extrusion::area(fvGeometry, scvf);
188
189 for (int compIdx = 0; compIdx < numComponents-1; compIdx++)
190 {
191 componentFlux[compIdx] = reducedFlux[compIdx];
192 componentFlux[numComponents-1] -= reducedFlux[compIdx];
193 }
194
195 return componentFlux ;
196 }
197
198private:
199 static Scalar calculateOmega_(const Scalar distance,
200 const Scalar extrusionFactor)
201 {
202 Scalar omega = 1/distance;
203 omega *= extrusionFactor;
204
205 return omega;
206 }
207
208 static ReducedComponentMatrix setupMSMatrix_(const Problem& problem,
209 const FVElementGeometry& fvGeometry,
210 const VolumeVariables& volVars,
211 const SubControlVolumeFace& scvf)
212 {
213 ReducedComponentMatrix reducedDiffusionMatrix(0.0);
214
215 for (int compIIdx = 0; compIIdx < numComponents-1; compIIdx++)
216 {
217 const auto xi = volVars.moleFraction(compIIdx);
218 const auto avgMolarMass = volVars.averageMolarMass(0);
219 const auto Mn = FluidSystem::molarMass(numComponents-1);
220 const Scalar tin = getEffectiveDiffusionCoefficient_(volVars, compIIdx, numComponents-1);
221
222 // set the entries of the diffusion matrix of the diagonal
223 reducedDiffusionMatrix[compIIdx][compIIdx] += xi*avgMolarMass/(tin*Mn);
224
225 for (int compJIdx = 0; compJIdx < numComponents; compJIdx++)
226 {
227 // we don't want to calculate e.g. water in water diffusion
228 if (compJIdx == compIIdx)
229 continue;
230
231 const auto xj = volVars.moleFraction(compJIdx);
232 const auto Mi = FluidSystem::molarMass(compIIdx);
233 const auto Mj = FluidSystem::molarMass(compJIdx);
234 const Scalar tij = getEffectiveDiffusionCoefficient_(volVars, compIIdx, compJIdx);
235 reducedDiffusionMatrix[compIIdx][compIIdx] += xj*avgMolarMass/(tij*Mi);
236 if (compJIdx < numComponents-1)
237 reducedDiffusionMatrix[compIIdx][compJIdx] += xi*(avgMolarMass/(tin*Mn) - avgMolarMass/(tij*Mj));
238 }
239 }
240 return reducedDiffusionMatrix;
241 }
242
243 static Scalar getEffectiveDiffusionCoefficient_(const VolumeVariables& volVars, const int phaseIdx, const int compIdx)
244 {
245 return volVars.effectiveDiffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
246 }
247};
248} // end namespace
249
250#endif
Container storing the diffusion coefficients required by the Maxwell- Stefan diffusion law....
Definition: maxwellstefandiffusioncoefficients.hh:33
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: staggered/freeflow/maxwellstefanslaw.hh:71
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:87
Definition: box/maxwellstefanslaw.hh:33
Defines all properties used in Dumux.
Helper classes to compute the integration elements.
Classes related to flux variables caching.
ReferenceSystemFormulation
The formulations available for Fick's law related to the reference system.
Definition: referencesystemformulation.hh:33
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:282
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Container storing the diffusion coefficients required by the Maxwell- Stefan diffusion law....
The available discretization methods in Dumux.
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
The empty filler class corresponding to EmptyCache.
Definition: fluxvariablescaching.hh:20
Definition: fluxvariablescaching.hh:55