Line | Branch | Exec | Source |
---|---|---|---|
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 | // | ||
7 | /*! | ||
8 | * \file | ||
9 | * \ingroup BoxFlux | ||
10 | * \brief This file contains the data which is required to calculate | ||
11 | * diffusive molar fluxes due to molecular diffusion with Maxwell Stefan | ||
12 | */ | ||
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 | |||
19 | #include <dumux/common/properties.hh> | ||
20 | #include <dumux/common/parameters.hh> | ||
21 | #include <dumux/discretization/method.hh> | ||
22 | #include <dumux/discretization/extrusion.hh> | ||
23 | |||
24 | #include <dumux/flux/maxwellstefandiffusioncoefficients.hh> | ||
25 | #include <dumux/flux/fluxvariablescaching.hh> | ||
26 | #include <dumux/flux/referencesystemformulation.hh> | ||
27 | #include <dumux/flux/facetensoraverage.hh> | ||
28 | |||
29 | namespace Dumux { | ||
30 | |||
31 | // forward declaration | ||
32 | template <class TypeTag, class DiscretizationMethod, ReferenceSystemFormulation referenceSystem> | ||
33 | class MaxwellStefansLawImplementation; | ||
34 | |||
35 | /*! | ||
36 | * \ingroup BoxFlux | ||
37 | * \brief Specialization of Maxwell Stefan's Law for the Box method. | ||
38 | */ | ||
39 | template <class TypeTag, ReferenceSystemFormulation referenceSystem> | ||
40 | class MaxwellStefansLawImplementation<TypeTag, DiscretizationMethods::Box, referenceSystem> | ||
41 | { | ||
42 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
43 | using Problem = GetPropType<TypeTag, Properties::Problem>; | ||
44 | using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>; | ||
45 | using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>; | ||
46 | using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>; | ||
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; | ||
52 | using GridView = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView; | ||
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 | |||
64 | public: | ||
65 | using DiffusionCoefficientsContainer = MaxwellStefanDiffusionCoefficients<Scalar, numFluidPhases, numComponents>; | ||
66 | |||
67 | //return the reference system | ||
68 | static constexpr ReferenceSystemFormulation referenceSystemFormulation() | ||
69 | { return referenceSystem; } | ||
70 | |||
71 | /*! | ||
72 | * \brief Returns the diffusive fluxes of all components within | ||
73 | * a fluid phase across the given sub-control volume face. | ||
74 | * The computed fluxes are given in kg/s. | ||
75 | */ | ||
76 | 764400 | 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 | 764400 | ComponentFluxVector componentFlux(0.0); | |
87 | 764400 | ReducedComponentMatrix reducedDiffusionMatrix(0.0); | |
88 | 764400 | ReducedComponentVector reducedFlux(0.0); | |
89 | 764400 | ComponentFluxVector moleFrac(0.0); | |
90 | 764400 | ReducedComponentVector normalX(0.0); | |
91 | |||
92 | // evaluate gradX at integration point and interpolate density | ||
93 | 764400 | const auto& fluxVarsCache = elemFluxVarsCache[scvf]; | |
94 | 764400 | const auto& shapeValues = fluxVarsCache.shapeValues(); | |
95 | |||
96 | 764400 | Scalar rho(0.0); | |
97 | 764400 | Scalar avgMolarMass(0.0); | |
98 |
4/4✓ Branch 0 taken 3057600 times.
✓ Branch 1 taken 764400 times.
✓ Branch 2 taken 3057600 times.
✓ Branch 3 taken 764400 times.
|
4586400 | for (auto&& scv : scvs(fvGeometry)) |
99 | { | ||
100 | 3057600 | const auto& volVars = elemVolVars[scv]; | |
101 | |||
102 | // density interpolation | ||
103 | 6115200 | rho += volVars.density(phaseIdx)*shapeValues[scv.indexInElement()][0]; | |
104 | //average molar mass interpolation | ||
105 | 3057600 | avgMolarMass += volVars.averageMolarMass(phaseIdx)*shapeValues[scv.indexInElement()][0]; | |
106 | //interpolate the mole fraction for the diffusion matrix | ||
107 |
2/2✓ Branch 0 taken 9172800 times.
✓ Branch 1 taken 3057600 times.
|
12230400 | for (int compIdx = 0; compIdx < numComponents; compIdx++) |
108 | { | ||
109 | 27518400 | moleFrac[compIdx] += volVars.moleFraction(phaseIdx, compIdx)*shapeValues[scv.indexInElement()][0]; | |
110 | } | ||
111 | } | ||
112 | |||
113 | 764400 | reducedDiffusionMatrix = setupMSMatrix_(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, moleFrac, avgMolarMass); | |
114 | |||
115 |
2/2✓ Branch 0 taken 1528800 times.
✓ Branch 1 taken 764400 times.
|
2293200 | for (int compIdx = 0; compIdx < numComponents-1; compIdx++) |
116 | { | ||
117 | 1528800 | Dune::FieldVector<Scalar, GridView::dimensionworld> gradX(0.0); | |
118 |
4/4✓ Branch 0 taken 6115200 times.
✓ Branch 1 taken 1528800 times.
✓ Branch 2 taken 6115200 times.
✓ Branch 3 taken 1528800 times.
|
15288000 | for (auto&& scv : scvs(fvGeometry)) |
119 | { | ||
120 | 6115200 | const auto& volVars = elemVolVars[scv]; | |
121 | |||
122 | // the mole/mass fraction gradient | ||
123 | 24460800 | gradX.axpy(volVars.moleFraction(phaseIdx, compIdx), fluxVarsCache.gradN(scv.indexInElement())); | |
124 | } | ||
125 | |||
126 | 6115200 | normalX[compIdx] = gradX *scvf.unitOuterNormal(); | |
127 | } | ||
128 | 764400 | reducedDiffusionMatrix.solve(reducedFlux,normalX); | |
129 | 1528800 | reducedFlux *= -1.0*rho*Extrusion::area(fvGeometry, scvf); | |
130 | |||
131 |
2/2✓ Branch 0 taken 1528800 times.
✓ Branch 1 taken 764400 times.
|
2293200 | for (int compIdx = 0; compIdx < numComponents-1; compIdx++) |
132 | { | ||
133 | 3057600 | componentFlux[compIdx] = reducedFlux[compIdx]; | |
134 | 4586400 | componentFlux[numComponents-1] -= reducedFlux[compIdx]; | |
135 | } | ||
136 | 764400 | return componentFlux ; | |
137 | } | ||
138 | |||
139 | private: | ||
140 | |||
141 | 764400 | 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 | 764400 | ReducedComponentMatrix reducedDiffusionMatrix(0.0); | |
151 | |||
152 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 764400 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 764400 times.
|
1528800 | const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()]; |
153 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 764400 times.
|
764400 | 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 | 2293200 | if(Dune::FloatCmp::eq<Scalar>(insideVolVars.saturation(phaseIdx), 0) || Dune::FloatCmp::eq<Scalar>(outsideVolVars.saturation(phaseIdx), 0)) | |
157 | return reducedDiffusionMatrix; | ||
158 | |||
159 |
2/2✓ Branch 0 taken 1528800 times.
✓ Branch 1 taken 764400 times.
|
2293200 | for (int compIIdx = 0; compIIdx < numComponents-1; compIIdx++) |
160 | { | ||
161 | //calculate diffusivity for i,numComponents | ||
162 | 1528800 | const auto xi = moleFrac[compIIdx]; | |
163 | 1528800 | const auto Mn = FluidSystem::molarMass(numComponents-1); | |
164 | |||
165 | 1528800 | auto tinInside = insideVolVars.effectiveDiffusionCoefficient(phaseIdx, compIIdx, numComponents-1); | |
166 | 1528800 | auto tinOutside = outsideVolVars.effectiveDiffusionCoefficient(phaseIdx, compIIdx, numComponents-1); | |
167 | |||
168 | // scale by extrusion factor | ||
169 | 1528800 | tinInside *= insideVolVars.extrusionFactor(); | |
170 | 1528800 | tinOutside *= outsideVolVars.extrusionFactor(); | |
171 | |||
172 | // the resulting averaged diffusion tensor | ||
173 |
2/4✓ Branch 0 taken 1528800 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 1528800 times.
✗ Branch 3 not taken.
|
3057600 | const auto tin = faceTensorAverage(tinInside, tinOutside, scvf.unitOuterNormal()); |
174 | |||
175 | // begin with the entries of the diffusion matrix of the diagonal | ||
176 | 3057600 | reducedDiffusionMatrix[compIIdx][compIIdx] += xi*avgMolarMass/(tin*Mn); | |
177 | |||
178 | // now set the rest of the entries (off-diagonal and additional entries for diagonal) | ||
179 |
2/2✓ Branch 0 taken 4586400 times.
✓ Branch 1 taken 1528800 times.
|
6115200 | for (int compJIdx = 0; compJIdx < numComponents; compJIdx++) |
180 | { | ||
181 | //we don't want to calculate e.g. water in water diffusion | ||
182 |
2/2✓ Branch 0 taken 3057600 times.
✓ Branch 1 taken 1528800 times.
|
4586400 | if (compIIdx == compJIdx) |
183 | continue; | ||
184 | |||
185 | //calculate diffusivity for compIIdx, compJIdx | ||
186 | 3057600 | const auto xj = moleFrac[compJIdx]; | |
187 | 3057600 | const auto Mi = FluidSystem::molarMass(compIIdx); | |
188 | 3057600 | const auto Mj = FluidSystem::molarMass(compJIdx); | |
189 | |||
190 | // effective diffusion tensors | ||
191 | 3057600 | auto tijInside = insideVolVars.effectiveDiffusionCoefficient(phaseIdx, compIIdx, compJIdx); | |
192 | 3057600 | auto tijOutside = outsideVolVars.effectiveDiffusionCoefficient(phaseIdx, compIIdx, compJIdx); | |
193 | |||
194 | // scale by extrusion factor | ||
195 | 3057600 | tijInside *= insideVolVars.extrusionFactor(); | |
196 | 3057600 | tijOutside *= outsideVolVars.extrusionFactor(); | |
197 | |||
198 | // the resulting averaged diffusion tensor | ||
199 |
1/2✓ Branch 0 taken 3057600 times.
✗ Branch 1 not taken.
|
3057600 | const auto tij = faceTensorAverage(tijInside, tijOutside, scvf.unitOuterNormal()); |
200 | |||
201 |
4/4✓ Branch 0 taken 1528800 times.
✓ Branch 1 taken 1528800 times.
✓ Branch 2 taken 1528800 times.
✓ Branch 3 taken 1528800 times.
|
6115200 | reducedDiffusionMatrix[compIIdx][compIIdx] += xj*avgMolarMass/(tij*Mi); |
202 |
2/2✓ Branch 0 taken 1528800 times.
✓ Branch 1 taken 1528800 times.
|
3057600 | if (compJIdx < numComponents-1) |
203 | 4586400 | 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 | ||
213 |