GCC Code Coverage Report


Directory: ../../../builds/dumux-repositories/
File: /builds/dumux-repositories/dumux/dumux/flux/cctpfa/maxwellstefanslaw.hh
Date: 2024-09-21 20:52:54
Exec Total Coverage
Lines: 72 73 98.6%
Functions: 15 15 100.0%
Branches: 46 50 92.0%

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 CCTpfaFlux
10 * \brief This file contains the data which is required to calculate
11 * diffusive mass fluxes due to molecular diffusion with Maxwell Stefan
12 */
13 #ifndef DUMUX_DISCRETIZATION_CC_TPFA_MAXWELL_STEFAN_LAW_HH
14 #define DUMUX_DISCRETIZATION_CC_TPFA_MAXWELL_STEFAN_LAW_HH
15
16 #include <dune/common/fmatrix.hh>
17 #include <dune/common/float_cmp.hh>
18
19 #include <dumux/common/properties.hh>
20 #include <dumux/common/parameters.hh>
21
22 #include <dumux/discretization/method.hh>
23 #include <dumux/discretization/extrusion.hh>
24 #include <dumux/flux/fluxvariablescaching.hh>
25 #include <dumux/flux/referencesystemformulation.hh>
26 #include <dumux/flux/maxwellstefandiffusioncoefficients.hh>
27
28 namespace Dumux {
29
30 // forward declaration
31 template <class TypeTag, class DiscretizationMethod, ReferenceSystemFormulation referenceSystem>
32 class MaxwellStefansLawImplementation;
33
34 /*!
35 * \ingroup CCTpfaFlux
36 * \brief Specialization of Maxwell Stefan's Law for the CCTpfa method.
37 */
38 template <class TypeTag, ReferenceSystemFormulation referenceSystem>
39 class MaxwellStefansLawImplementation<TypeTag, DiscretizationMethods::CCTpfa, referenceSystem >
40 {
41 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
42 using Problem = GetPropType<TypeTag, Properties::Problem>;
43 using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>;
44 using FVElementGeometry = typename GridGeometry::LocalView;
45 using SubControlVolume = typename GridGeometry::SubControlVolume;
46 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
47 using Extrusion = Extrusion_t<GridGeometry>;
48 using GridView = typename GridGeometry::GridView;
49 using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>;
50 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
51 using Element = typename GridView::template Codim<0>::Entity;
52 using GridFluxVariablesCache = GetPropType<TypeTag, Properties::GridFluxVariablesCache>;
53 using ElementFluxVariablesCache = typename GridFluxVariablesCache::LocalView;
54 using FluxVariablesCache = typename GridFluxVariablesCache::FluxVariablesCache;
55 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
56
57 static const int dim = GridView::dimension;
58 static const int dimWorld = GridView::dimensionworld;
59 static const int numPhases = GetPropType<TypeTag, Properties::ModelTraits>::numFluidPhases();
60 static const int numComponents = GetPropType<TypeTag, Properties::ModelTraits>::numFluidComponents();
61
62 static_assert(referenceSystem == ReferenceSystemFormulation::massAveraged, "only the mass averaged reference system is supported for the Maxwell-Stefan formulation");
63
64
65 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
66 using ReducedComponentVector = Dune::FieldVector<Scalar, numComponents-1>;
67 using ReducedComponentMatrix = Dune::FieldMatrix<Scalar, numComponents-1, numComponents-1>;
68
69 public:
70 using DiscretizationMethod = DiscretizationMethods::CCTpfa;
71 // state the discretization method this implementation belongs to
72 static constexpr DiscretizationMethod discMethod{};
73
74 //return the reference system
75 static constexpr ReferenceSystemFormulation referenceSystemFormulation()
76 { return referenceSystem; }
77
78 using DiffusionCoefficientsContainer = MaxwellStefanDiffusionCoefficients<Scalar, numPhases, numComponents>;
79
80 //! state the type for the corresponding cache and its filler
81 //! We don't cache anything for this law
82 using Cache = FluxVariablesCaching::EmptyDiffusionCache;
83 using CacheFiller = FluxVariablesCaching::EmptyCacheFiller;
84
85 /*!
86 * \brief Returns the diffusive fluxes of all components within
87 * a fluid phase across the given sub-control volume face.
88 * The computed fluxes are given in kg/s.
89 */
90 3058604 static ComponentFluxVector flux(const Problem& problem,
91 const Element& element,
92 const FVElementGeometry& fvGeometry,
93 const ElementVolumeVariables& elemVolVars,
94 const SubControlVolumeFace& scvf,
95 const int phaseIdx,
96 const ElementFluxVariablesCache& elemFluxVarsCache)
97 {
98 //this is to calculate the maxwellStefan diffusion in a multicomponent system.
99 //see: Multicomponent Mass Transfer. R. Taylor u. R. Krishna. J. Wiley & Sons, New York 1993
100 3058604 ComponentFluxVector componentFlux(0.0);
101 3058604 ReducedComponentVector moleFracInside(0.0);
102 3058604 ReducedComponentVector moleFracOutside(0.0);
103 3058604 ReducedComponentVector reducedFlux(0.0);
104 3058604 ReducedComponentMatrix reducedDiffusionMatrixInside(0.0);
105 3058604 ReducedComponentMatrix reducedDiffusionMatrixOutside(0.0);
106
107 // get inside/outside volume variables
108 6117208 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
109 6117208 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
110 3058604 const auto rhoInside = insideVolVars.density(phaseIdx);
111 3058604 const auto rhoOutside = outsideVolVars.density(phaseIdx);
112 //calculate the mole fraction vectors
113
2/2
✓ Branch 0 taken 4867628 times.
✓ Branch 1 taken 3058604 times.
7926232 for (int compIdx = 0; compIdx < numComponents-1; compIdx++)
114 {
115 //calculate x_inside
116 4867628 const auto xInside = insideVolVars.moleFraction(phaseIdx, compIdx);
117 //calculate outside molefraction with the respective transmissibility
118 4867628 const auto xOutside = outsideVolVars.moleFraction(phaseIdx, compIdx);
119
120 4867628 moleFracInside[compIdx] = xInside;
121 8485676 moleFracOutside[compIdx] = xOutside;
122 }
123
124 //we cannot solve that if the matrix is 0 everywhere
125
12/12
✓ Branch 0 taken 1078988 times.
✓ Branch 1 taken 1033056 times.
✓ Branch 2 taken 1078988 times.
✓ Branch 3 taken 1033056 times.
✓ Branch 4 taken 444708 times.
✓ Branch 5 taken 428472 times.
✓ Branch 6 taken 442266 times.
✓ Branch 7 taken 2442 times.
✓ Branch 8 taken 442266 times.
✓ Branch 9 taken 2442 times.
✓ Branch 10 taken 442266 times.
✓ Branch 11 taken 2442 times.
6990388 if(!(Dune::FloatCmp::eq<Scalar>(insideVolVars.saturation(phaseIdx), 0) || Dune::FloatCmp::eq<Scalar>(outsideVolVars.saturation(phaseIdx), 0)))
126 {
127
2/2
✓ Branch 0 taken 868118 times.
✓ Branch 1 taken 813012 times.
2627690 const auto insideScvIdx = scvf.insideScvIdx();
128
2/2
✓ Branch 0 taken 868118 times.
✓ Branch 1 taken 813012 times.
2627690 const auto& insideScv = fvGeometry.scv(insideScvIdx);
129 2627690 const Scalar omegai = calculateOmega_(scvf,
130 insideScv,
131 insideVolVars.extrusionFactor());
132
133 //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
134 2627690 reducedDiffusionMatrixInside = setupMSMatrix_(problem, element, fvGeometry, insideVolVars, insideScv, phaseIdx);
135
136 //if on boundary
137
6/6
✓ Branch 0 taken 2602280 times.
✓ Branch 1 taken 25410 times.
✓ Branch 2 taken 416856 times.
✓ Branch 3 taken 2185424 times.
✓ Branch 4 taken 416856 times.
✓ Branch 5 taken 2185424 times.
2627690 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
138 {
139 25410 moleFracOutside -= moleFracInside;
140 25410 reducedDiffusionMatrixInside.solve(reducedFlux, moleFracOutside);
141 25410 reducedFlux *= omegai*rhoInside;
142 }
143 else //we need outside cells as well if we are not on the boundary
144 {
145 Scalar omegaj;
146
2/2
✓ Branch 0 taken 813012 times.
✓ Branch 1 taken 842708 times.
2602280 const auto outsideScvIdx = scvf.outsideScvIdx();
147
2/2
✓ Branch 0 taken 813012 times.
✓ Branch 1 taken 842708 times.
2602280 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
148
149 2602280 reducedDiffusionMatrixOutside = setupMSMatrix_(problem, element, fvGeometry, outsideVolVars, outsideScv, phaseIdx);
150
151 if (dim == dimWorld)
152 // assume the normal vector from outside is anti parallel so we save flipping a vector
153 2602280 omegaj = -1.0*calculateOmega_(scvf,
154 outsideScv,
155 outsideVolVars.extrusionFactor());
156 else
157 omegaj = calculateOmega_(fvGeometry.flipScvf(scvf.index()),
158 outsideScv,
159 outsideVolVars.extrusionFactor());
160
161 2602280 reducedDiffusionMatrixInside.invert();
162 2602280 reducedDiffusionMatrixOutside.invert();
163 2602280 reducedDiffusionMatrixInside *= omegai*rhoInside;
164 2602280 reducedDiffusionMatrixOutside *= omegaj*rhoOutside;
165
166 //in the helpervector we store the values for x*
167 2602280 ReducedComponentVector helperVector(0.0);
168 2602280 ReducedComponentVector gradientVectori(0.0);
169 2602280 ReducedComponentVector gradientVectorj(0.0);
170
171 793256 reducedDiffusionMatrixInside.mv(moleFracInside, gradientVectori);
172 2602280 reducedDiffusionMatrixOutside.mv(moleFracOutside, gradientVectorj);
173
174 2602280 auto gradientVectorij = (gradientVectori + gradientVectorj);
175
176 //add the two matrixes to each other
177 3395536 reducedDiffusionMatrixOutside += reducedDiffusionMatrixInside;
178
179 2602280 reducedDiffusionMatrixOutside.solve(helperVector, gradientVectorij);
180
181 //Bi^-1 omegai rho(x*-xi)
182 helperVector -=moleFracInside;
183 1809024 reducedDiffusionMatrixInside.mv(helperVector, reducedFlux);
184 }
185
186 5255380 reducedFlux *= -Extrusion::area(fvGeometry, scvf);
187
2/2
✓ Branch 0 taken 4436714 times.
✓ Branch 1 taken 2627690 times.
7064404 for (int compIdx = 0; compIdx < numComponents-1; compIdx++)
188 {
189 8054762 componentFlux[compIdx] = reducedFlux[compIdx];
190 8873428 componentFlux[numComponents-1] -=reducedFlux[compIdx];
191 }
192 }
193 3058604 return componentFlux ;
194 }
195
196 private:
197 5229970 static Scalar calculateOmega_(const SubControlVolumeFace& scvf,
198 const SubControlVolume &scv,
199 const Scalar extrusionFactor)
200 {
201 5229970 auto distanceVector = scvf.ipGlobal();
202 10459940 distanceVector -= scv.center();
203 distanceVector /= distanceVector.two_norm2();
204
205 5229970 Scalar omega = (distanceVector * scvf.unitOuterNormal());
206 5229970 omega *= extrusionFactor;
207
208 5229970 return omega;
209 }
210
211 5229970 static ReducedComponentMatrix setupMSMatrix_(const Problem& problem,
212 const Element& element,
213 const FVElementGeometry& fvGeometry,
214 const VolumeVariables& volVars,
215 const SubControlVolume& scv,
216 const int phaseIdx)
217 {
218
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 859122 times.
5229970 ReducedComponentMatrix reducedDiffusionMatrix(0.0);
219
220 //this is to not divide by 0 if the saturation in 0 and the effectiveDiffusionCoefficient becomes zero due to that
221
3/6
✗ Branch 0 not taken.
✓ Branch 1 taken 859122 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 859122 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 859122 times.
6948214 if(Dune::FloatCmp::eq<Scalar>(volVars.saturation(phaseIdx), 0))
222 return reducedDiffusionMatrix;
223
224 5229970 const auto avgMolarMass = volVars.averageMolarMass(phaseIdx);
225
226
2/2
✓ Branch 0 taken 8848018 times.
✓ Branch 1 taken 5229970 times.
14077988 for (int compIIdx = 0; compIIdx < numComponents-1; compIIdx++)
227 {
228 8848018 const auto xi = volVars.moleFraction(phaseIdx, compIIdx);
229
230 8848018 const auto Mn = FluidSystem::molarMass(numComponents-1);
231 8848018 Scalar tin = volVars.effectiveDiffusionCoefficient(phaseIdx, compIIdx, numComponents-1);
232
233 // set the entries of the diffusion matrix of the diagonal
234 16084114 reducedDiffusionMatrix[compIIdx][compIIdx] += xi*avgMolarMass/(tin*Mn);
235
236 // now set the rest of the entries (off-diagonal and additional entries for diagonal)
237
2/2
✓ Branch 0 taken 23320210 times.
✓ Branch 1 taken 10459940 times.
33780150 for (int compJIdx = 0; compJIdx < numComponents; compJIdx++)
238 {
239 // we don't want to calculate e.g. water in water diffusion
240
2/2
✓ Branch 0 taken 16084114 times.
✓ Branch 1 taken 8848018 times.
24932132 if (compJIdx == compIIdx)
241 continue;
242
243 16084114 const auto xj = volVars.moleFraction(phaseIdx, compJIdx);
244 16084114 const auto Mi = FluidSystem::molarMass(compIIdx);
245 16084114 const auto Mj = FluidSystem::molarMass(compJIdx);
246 16084114 Scalar tij = volVars.effectiveDiffusionCoefficient(phaseIdx, compIIdx, compJIdx);
247
4/4
✓ Branch 0 taken 7236096 times.
✓ Branch 1 taken 7236096 times.
✓ Branch 2 taken 7236096 times.
✓ Branch 3 taken 7236096 times.
32168228 reducedDiffusionMatrix[compIIdx][compIIdx] += xj*avgMolarMass/(tij*Mi);
248
249
2/2
✓ Branch 0 taken 7236096 times.
✓ Branch 1 taken 7236096 times.
14472192 if (compJIdx < numComponents-1)
250 21708288 reducedDiffusionMatrix[compIIdx][compJIdx] += xi*(avgMolarMass/(tin*Mn) - avgMolarMass/(tij*Mj));
251 }
252 }
253 1611922 return reducedDiffusionMatrix;
254 }
255
256 };
257 } // end namespace
258
259 #endif
260