version 3.8
box/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//
13#ifndef DUMUX_DISCRETIZATION_BOX_MAXWELL_STEFAN_LAW_HH
14#define DUMUX_DISCRETIZATION_BOX_MAXWELL_STEFAN_LAW_HH
15
16#include <dune/common/float_cmp.hh>
17#include <dune/common/fmatrix.hh>
18
23
28
29namespace Dumux {
30
31// forward declaration
32template <class TypeTag, class DiscretizationMethod, ReferenceSystemFormulation referenceSystem>
34
39template <class TypeTag, ReferenceSystemFormulation referenceSystem>
40class MaxwellStefansLawImplementation<TypeTag, DiscretizationMethods::Box, referenceSystem>
41{
47 using FVElementGeometry = typename GridGeometry::LocalView;
48 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
49 using Extrusion = Extrusion_t<GridGeometry>;
50 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
51 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
53 using Element = typename GridView::template Codim<0>::Entity;
54
55 static constexpr auto numFluidPhases = GetPropType<TypeTag, Properties::ModelTraits>::numFluidPhases();
56 static constexpr auto numComponents = GetPropType<TypeTag, Properties::ModelTraits>::numFluidComponents();
57
58 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
59 using ReducedComponentVector = Dune::FieldVector<Scalar, numComponents-1>;
60 using ReducedComponentMatrix = Dune::FieldMatrix<Scalar, numComponents-1, numComponents-1>;
61
62 static_assert(referenceSystem == ReferenceSystemFormulation::massAveraged, "only the mass averaged reference system is supported for the Maxwell-Stefan formulation");
63
64public:
66
67 //return the reference system
69 { return referenceSystem; }
70
76 static ComponentFluxVector flux(const Problem& problem,
77 const Element& element,
78 const FVElementGeometry& fvGeometry,
79 const ElementVolumeVariables& elemVolVars,
80 const SubControlVolumeFace& scvf,
81 const int phaseIdx,
82 const ElementFluxVariablesCache& elemFluxVarsCache)
83 {
84 //this is to calculate the maxwellStefan diffusion in a multicomponent system.
85 //see: Multicomponent Mass Transfer. R. Taylor u. R. Krishna. J. Wiley & Sons, New York 1993
86 ComponentFluxVector componentFlux(0.0);
87 ReducedComponentMatrix reducedDiffusionMatrix(0.0);
88 ReducedComponentVector reducedFlux(0.0);
89 ComponentFluxVector moleFrac(0.0);
90 ReducedComponentVector normalX(0.0);
91
92 // evaluate gradX at integration point and interpolate density
93 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
94 const auto& shapeValues = fluxVarsCache.shapeValues();
95
96 Scalar rho(0.0);
97 Scalar avgMolarMass(0.0);
98 for (auto&& scv : scvs(fvGeometry))
99 {
100 const auto& volVars = elemVolVars[scv];
101
102 // density interpolation
103 rho += volVars.density(phaseIdx)*shapeValues[scv.indexInElement()][0];
104 //average molar mass interpolation
105 avgMolarMass += volVars.averageMolarMass(phaseIdx)*shapeValues[scv.indexInElement()][0];
106 //interpolate the mole fraction for the diffusion matrix
107 for (int compIdx = 0; compIdx < numComponents; compIdx++)
108 {
109 moleFrac[compIdx] += volVars.moleFraction(phaseIdx, compIdx)*shapeValues[scv.indexInElement()][0];
110 }
111 }
112
113 reducedDiffusionMatrix = setupMSMatrix_(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, moleFrac, avgMolarMass);
114
115 for (int compIdx = 0; compIdx < numComponents-1; compIdx++)
116 {
117 Dune::FieldVector<Scalar, GridView::dimensionworld> gradX(0.0);
118 for (auto&& scv : scvs(fvGeometry))
119 {
120 const auto& volVars = elemVolVars[scv];
121
122 // the mole/mass fraction gradient
123 gradX.axpy(volVars.moleFraction(phaseIdx, compIdx), fluxVarsCache.gradN(scv.indexInElement()));
124 }
125
126 normalX[compIdx] = gradX *scvf.unitOuterNormal();
127 }
128 reducedDiffusionMatrix.solve(reducedFlux,normalX);
129 reducedFlux *= -1.0*rho*Extrusion::area(fvGeometry, scvf);
130
131 for (int compIdx = 0; compIdx < numComponents-1; compIdx++)
132 {
133 componentFlux[compIdx] = reducedFlux[compIdx];
134 componentFlux[numComponents-1] -= reducedFlux[compIdx];
135 }
136 return componentFlux ;
137 }
138
139private:
140
141 static ReducedComponentMatrix setupMSMatrix_(const Problem& problem,
142 const Element& element,
143 const FVElementGeometry& fvGeometry,
144 const ElementVolumeVariables& elemVolVars,
145 const SubControlVolumeFace& scvf,
146 const int phaseIdx,
147 const ComponentFluxVector moleFrac,
148 const Scalar avgMolarMass)
149 {
150 ReducedComponentMatrix reducedDiffusionMatrix(0.0);
151
152 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
153 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
154
155 //this is to not divide by 0 if the saturation in 0 and the effectiveDiffusionCoefficient becomes zero due to that
156 if(Dune::FloatCmp::eq<Scalar>(insideVolVars.saturation(phaseIdx), 0) || Dune::FloatCmp::eq<Scalar>(outsideVolVars.saturation(phaseIdx), 0))
157 return reducedDiffusionMatrix;
158
159 for (int compIIdx = 0; compIIdx < numComponents-1; compIIdx++)
160 {
161 //calculate diffusivity for i,numComponents
162 const auto xi = moleFrac[compIIdx];
163 const auto Mn = FluidSystem::molarMass(numComponents-1);
164
165 auto tinInside = insideVolVars.effectiveDiffusionCoefficient(phaseIdx, compIIdx, numComponents-1);
166 auto tinOutside = outsideVolVars.effectiveDiffusionCoefficient(phaseIdx, compIIdx, numComponents-1);
167
168 // scale by extrusion factor
169 tinInside *= insideVolVars.extrusionFactor();
170 tinOutside *= outsideVolVars.extrusionFactor();
171
172 // the resulting averaged diffusion tensor
173 const auto tin = faceTensorAverage(tinInside, tinOutside, scvf.unitOuterNormal());
174
175 // begin with the entries of the diffusion matrix of the diagonal
176 reducedDiffusionMatrix[compIIdx][compIIdx] += xi*avgMolarMass/(tin*Mn);
177
178 // now set the rest of the entries (off-diagonal and additional entries for diagonal)
179 for (int compJIdx = 0; compJIdx < numComponents; compJIdx++)
180 {
181 //we don't want to calculate e.g. water in water diffusion
182 if (compIIdx == compJIdx)
183 continue;
184
185 //calculate diffusivity for compIIdx, compJIdx
186 const auto xj = moleFrac[compJIdx];
187 const auto Mi = FluidSystem::molarMass(compIIdx);
188 const auto Mj = FluidSystem::molarMass(compJIdx);
189
190 // effective diffusion tensors
191 auto tijInside = insideVolVars.effectiveDiffusionCoefficient(phaseIdx, compIIdx, compJIdx);
192 auto tijOutside = outsideVolVars.effectiveDiffusionCoefficient(phaseIdx, compIIdx, compJIdx);
193
194 // scale by extrusion factor
195 tijInside *= insideVolVars.extrusionFactor();
196 tijOutside *= outsideVolVars.extrusionFactor();
197
198 // the resulting averaged diffusion tensor
199 const auto tij = faceTensorAverage(tijInside, tijOutside, scvf.unitOuterNormal());
200
201 reducedDiffusionMatrix[compIIdx][compIIdx] += xj*avgMolarMass/(tij*Mi);
202 if (compJIdx < numComponents-1)
203 reducedDiffusionMatrix[compIIdx][compJIdx] +=xi*(avgMolarMass/(tin*Mn) - avgMolarMass/(tij*Mj));
204 }
205 }
206 return reducedDiffusionMatrix;
207 }
208
209};
210} // end namespace Dumux
211
212#endif
Container storing the diffusion coefficients required by the Maxwell- Stefan diffusion law....
Definition: maxwellstefandiffusioncoefficients.hh:33
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: box/maxwellstefanslaw.hh:68
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: box/maxwellstefanslaw.hh:76
Definition: box/maxwellstefanslaw.hh:33
Defines all properties used in Dumux.
Helper classes to compute the integration elements.
A free function to average a Tensor at an interface.
Classes related to flux variables caching.
ReferenceSystemFormulation
The formulations available for Fick's law related to the reference system.
Definition: referencesystemformulation.hh:33
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.
CVFE< CVFEMethods::PQ1 > Box
Definition: method.hh:94
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
Scalar faceTensorAverage(const Scalar T1, const Scalar T2, const Dune::FieldVector< Scalar, dim > &normal)
Average of a discontinuous scalar field at discontinuity interface (for compatibility reasons with th...
Definition: facetensoraverage.hh:29
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...