GCC Code Coverage Report


Directory: ../../../builds/dumux-repositories/
File: /builds/dumux-repositories/dumux/dumux/flux/staggered/freeflow/maxwellstefanslaw.hh
Date: 2024-09-21 20:52:54
Exec Total Coverage
Lines: 76 76 100.0%
Functions: 4 4 100.0%
Branches: 30 32 93.8%

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 StaggeredFlux
10 * \brief Specialization of Maxwell Stefan's Law for the Staggered method.
11 */
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
17 #include <dumux/common/properties.hh>
18 #include <dumux/common/parameters.hh>
19 #include <dumux/discretization/method.hh>
20 #include <dumux/discretization/extrusion.hh>
21
22 #include <dumux/flux/fluxvariablescaching.hh>
23 #include <dumux/flux/referencesystemformulation.hh>
24 #include <dumux/flux/maxwellstefandiffusioncoefficients.hh>
25
26 namespace Dumux {
27
28 // forward declaration
29 template <class TypeTag, class DiscretizationMethod, ReferenceSystemFormulation referenceSystem>
30 class MaxwellStefansLawImplementation;
31
32 /*!
33 * \ingroup StaggeredFlux
34 * \brief Specialization of Maxwell Stefan's Law for the Staggered method.
35 */
36 template <class TypeTag, ReferenceSystemFormulation referenceSystem>
37 class MaxwellStefansLawImplementation<TypeTag, DiscretizationMethods::Staggered, referenceSystem>
38 {
39 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
40 using Problem = GetPropType<TypeTag, Properties::Problem>;
41 using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>;
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;
47 using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>;
48 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
49 using Indices = typename GetPropType<TypeTag, Properties::ModelTraits>::Indices;
50 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
51
52 using ModelTraits = GetPropType<TypeTag, Properties::ModelTraits>;
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
65 public:
66 using DiscretizationMethod = DiscretizationMethods::Staggered;
67 // state the discretization method this implementation belongs to
68 static constexpr DiscretizationMethod discMethod{};
69
70 //return the reference system
71 static constexpr ReferenceSystemFormulation referenceSystemFormulation()
72 { return referenceSystem; }
73
74 //! state the type for the corresponding cache and its filler
75 //! We don't cache anything for this law
76 using Cache = FluxVariablesCaching::EmptyDiffusionCache;
77 using CacheFiller = FluxVariablesCaching::EmptyCacheFiller;
78
79 using DiffusionCoefficientsContainer = MaxwellStefanDiffusionCoefficients<Scalar, numPhases, numComponents>;
80
81 /*!
82 * \brief Returns the diffusive fluxes of all components within
83 * a fluid phase across the given sub-control volume face.
84 * The computed fluxes are given in kg/s.
85 */
86 template<class ElementVolumeVariables>
87 2574720 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 2574720 CellCenterPrimaryVariables componentFlux(0.0);
96 2574720 ReducedComponentVector moleFracInside(0.0);
97 2574720 ReducedComponentVector moleFracOutside(0.0);
98 2574720 ReducedComponentVector reducedFlux(0.0);
99 2574720 ReducedComponentMatrix reducedDiffusionMatrixInside(0.0);
100 2574720 ReducedComponentMatrix reducedDiffusionMatrixOutside(0.0);
101
102 // get inside/outside volume variables
103 5149440 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
104 5149440 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
105 2574720 const auto rhoInside = insideVolVars.density();
106 2574720 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
2/2
✓ Branch 0 taken 6762584 times.
✓ Branch 1 taken 2546784 times.
9309368 for(int compIdx = 0; compIdx < numComponents; ++compIdx)
110 {
111
2/2
✓ Branch 0 taken 4187864 times.
✓ Branch 1 taken 2574720 times.
6762584 if(compIdx == FluidSystem::getMainComponent(0))
112 continue;
113
114 // get equation index
115 4187864 const auto eqIdx = Indices::conti0EqIdx + compIdx;
116
2/2
✓ Branch 0 taken 174120 times.
✓ Branch 1 taken 4013744 times.
4187864 if(scvf.boundary())
117 {
118 174120 const auto bcTypes = problem.boundaryTypes(element, scvf);
119
6/8
✓ Branch 0 taken 146184 times.
✓ Branch 1 taken 27936 times.
✓ Branch 2 taken 146184 times.
✓ Branch 3 taken 27936 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 146184 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 146184 times.
348240 if((bcTypes.isOutflow(eqIdx)) || (bcTypes.isSymmetry()))
120 27936 return componentFlux;
121 }
122 }
123
124 //calculate the mole fraction vectors
125
2/2
✓ Branch 0 taken 4159928 times.
✓ Branch 1 taken 2546784 times.
6706712 for (int compIdx = 0; compIdx < numComponents-1; compIdx++)
126 {
127 //calculate x_inside
128 4159928 const auto xInside = insideVolVars.moleFraction(compIdx);
129 //calculate outside molefraction with the respective transmissibility
130 4159928 const auto xOutside = outsideVolVars.moleFraction(compIdx);
131
132 4159928 moleFracInside[compIdx] = xInside;
133 7386216 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 2546784 reducedDiffusionMatrixInside = setupMSMatrix_(problem, fvGeometry, insideVolVars, scvf);
138
139 2546784 reducedDiffusionMatrixOutside = setupMSMatrix_(problem, fvGeometry, outsideVolVars, scvf);
140
141 2546784 const auto insideScvIdx = scvf.insideScvIdx();
142 2546784 const auto& insideScv = fvGeometry.scv(insideScvIdx);
143 2546784 const auto outsideScvIdx = scvf.outsideScvIdx();
144 2546784 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
145
146 10187136 const Scalar insideDistance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
147
148 5093568 const Scalar omegai = calculateOmega_(insideDistance, insideVolVars.extrusionFactor());
149
150 //if on boundary
151 2546784 if (scvf.boundary())
152 {
153 89712 moleFracOutside -= moleFracInside;
154 89712 reducedDiffusionMatrixInside.solve(reducedFlux, moleFracOutside);
155 89712 reducedFlux *= omegai*rhoInside;
156 }
157 else
158 {
159 9828288 const Scalar outsideDistance = (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm();
160 4914144 const Scalar omegaj = calculateOmega_(outsideDistance, outsideVolVars.extrusionFactor());
161
162 4914144 reducedDiffusionMatrixInside.invert();
163 2457072 reducedDiffusionMatrixOutside.invert();
164 2457072 reducedDiffusionMatrixInside *= omegai*rhoInside;
165 2457072 reducedDiffusionMatrixOutside *= omegaj*rhoOutside;
166
167 //in the helpervector we store the values for x*
168 2457072 ReducedComponentVector helperVector(0.0);
169 2457072 ReducedComponentVector gradientVectori(0.0);
170 2457072 ReducedComponentVector gradientVectorj(0.0);
171
172 900400 reducedDiffusionMatrixInside.mv(moleFracInside, gradientVectori);
173 2457072 reducedDiffusionMatrixOutside.mv(moleFracOutside, gradientVectorj);
174
175 2457072 auto gradientVectorij = (gradientVectori + gradientVectorj);
176
177 //add the two matrixes to each other
178 3357472 reducedDiffusionMatrixOutside += reducedDiffusionMatrixInside;
179
180 2457072 reducedDiffusionMatrixOutside.solve(helperVector, gradientVectorij);
181
182 //Bi^-1 omegai rho(x*-xi)
183 helperVector -=moleFracInside;
184 1556672 reducedDiffusionMatrixInside.mv(helperVector, reducedFlux);
185 }
186
187 5093568 reducedFlux *= -Extrusion::area(fvGeometry, scvf);
188
189
2/2
✓ Branch 0 taken 4159928 times.
✓ Branch 1 taken 2546784 times.
6706712 for (int compIdx = 0; compIdx < numComponents-1; compIdx++)
190 {
191 7386216 componentFlux[compIdx] = reducedFlux[compIdx];
192 11546144 componentFlux[numComponents-1] -= reducedFlux[compIdx];
193 }
194
195 933640 return componentFlux ;
196 }
197
198 private:
199 static Scalar calculateOmega_(const Scalar distance,
200 const Scalar extrusionFactor)
201 {
202 5003856 Scalar omega = 1/distance;
203
2/2
✓ Branch 0 taken 89712 times.
✓ Branch 1 taken 2457072 times.
5003856 omega *= extrusionFactor;
204
205 return omega;
206 }
207
208 5093568 static ReducedComponentMatrix setupMSMatrix_(const Problem& problem,
209 const FVElementGeometry& fvGeometry,
210 const VolumeVariables& volVars,
211 const SubControlVolumeFace& scvf)
212 {
213 5093568 ReducedComponentMatrix reducedDiffusionMatrix(0.0);
214
215
2/2
✓ Branch 0 taken 8319856 times.
✓ Branch 1 taken 5093568 times.
13413424 for (int compIIdx = 0; compIIdx < numComponents-1; compIIdx++)
216 {
217 8319856 const auto xi = volVars.moleFraction(compIIdx);
218 8319856 const auto avgMolarMass = volVars.averageMolarMass(0);
219 8319856 const auto Mn = FluidSystem::molarMass(numComponents-1);
220 16639712 const Scalar tin = getEffectiveDiffusionCoefficient_(volVars, compIIdx, numComponents-1);
221
222 // set the entries of the diffusion matrix of the diagonal
223 10187136 reducedDiffusionMatrix[compIIdx][compIIdx] += xi*avgMolarMass/(tin*Mn);
224
225
2/2
✓ Branch 0 taken 21225008 times.
✓ Branch 1 taken 10187136 times.
31412144 for (int compJIdx = 0; compJIdx < numComponents; compJIdx++)
226 {
227 // we don't want to calculate e.g. water in water diffusion
228
2/2
✓ Branch 0 taken 14772432 times.
✓ Branch 1 taken 8319856 times.
23092288 if (compJIdx == compIIdx)
229 continue;
230
231 14772432 const auto xj = volVars.moleFraction(compJIdx);
232 14772432 const auto Mi = FluidSystem::molarMass(compIIdx);
233 14772432 const auto Mj = FluidSystem::molarMass(compJIdx);
234 29544864 const Scalar tij = getEffectiveDiffusionCoefficient_(volVars, compIIdx, compJIdx);
235
2/2
✓ Branch 0 taken 6452576 times.
✓ Branch 1 taken 6452576 times.
14772432 reducedDiffusionMatrix[compIIdx][compIIdx] += xj*avgMolarMass/(tij*Mi);
236
2/2
✓ Branch 0 taken 6452576 times.
✓ Branch 1 taken 6452576 times.
12905152 if (compJIdx < numComponents-1)
237 19357728 reducedDiffusionMatrix[compIIdx][compJIdx] += xi*(avgMolarMass/(tin*Mn) - avgMolarMass/(tij*Mj));
238 }
239 }
240 5093568 return reducedDiffusionMatrix;
241 }
242
243 static Scalar getEffectiveDiffusionCoefficient_(const VolumeVariables& volVars, const int phaseIdx, const int compIdx)
244 {
245
2/2
✓ Branch 2 taken 6452576 times.
✓ Branch 3 taken 6452576 times.
23092288 return volVars.effectiveDiffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
246 }
247 };
248 } // end namespace
249
250 #endif
251