GCC Code Coverage Report


Directory: ../../../builds/dumux-repositories/
File: /builds/dumux-repositories/dumux/dumux/porousmediumflow/nonequilibrium/localresidual.hh
Date: 2024-09-21 20:52:54
Exec Total Coverage
Lines: 63 64 98.4%
Functions: 11 15 73.3%
Branches: 33 59 55.9%

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 NonEquilibriumModel
10 * \brief The local residual for the kinetic mass transfer module of
11 * the compositional multi-phase model.
12 */
13
14 #ifndef DUMUX_NONEQUILIBRIUM_LOCAL_RESIDUAL_HH
15 #define DUMUX_NONEQUILIBRIUM_LOCAL_RESIDUAL_HH
16
17 #include <cmath>
18 #include <dumux/common/properties.hh>
19 #include <dumux/common/numeqvector.hh>
20 #include <dumux/porousmediumflow/nonequilibrium/thermal/localresidual.hh>
21
22 namespace Dumux {
23
24 // forward declaration
25 template<class TypeTag, bool enableChemicalNonEquilibrium>
26 class NonEquilibriumLocalResidualImplementation;
27
28 template <class TypeTag>
29 using NonEquilibriumLocalResidual = NonEquilibriumLocalResidualImplementation<TypeTag, GetPropType<TypeTag, Properties::ModelTraits>::enableChemicalNonEquilibrium()>;
30
31 /*!
32 * \ingroup NonEquilibriumModel
33 * \brief The local residual for a model without chemical non-equilibrium
34 * but potentially with thermal non-equilibrium
35 */
36 template<class TypeTag>
37 class NonEquilibriumLocalResidualImplementation<TypeTag, false>: public GetPropType<TypeTag, Properties::EquilibriumLocalResidual>
38 {
39 using ParentType = GetPropType<TypeTag, Properties::EquilibriumLocalResidual>;
40 using Problem = GetPropType<TypeTag, Properties::Problem>;
41 using NumEqVector = Dumux::NumEqVector<GetPropType<TypeTag, Properties::PrimaryVariables>>;
42 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
43 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
44 using GridView = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView;
45 using Element = typename GridView::template Codim<0>::Entity;
46 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
47 using EnergyLocalResidual = GetPropType<TypeTag, Properties::EnergyLocalResidual>;
48 using ModelTraits = GetPropType<TypeTag, Properties::ModelTraits>;
49
50 public:
51
2/4
✓ Branch 1 taken 67600 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 67600 times.
✗ Branch 5 not taken.
135200 using ParentType::ParentType;
52
53 /*!
54 * \brief Calculates the source term of the equation.
55 *
56 * \param problem The source term
57 * \param element An element which contains part of the control volume
58 * \param fvGeometry The finite-volume geometry
59 * \param elemVolVars The volume variables of the current element
60 * \param scv The sub-control volume over which we integrate the source term
61 */
62 2505200 NumEqVector computeSource(const Problem& problem,
63 const Element& element,
64 const FVElementGeometry& fvGeometry,
65 const ElementVolumeVariables& elemVolVars,
66 const SubControlVolume &scv) const
67 {
68 2505200 NumEqVector source(0.0);
69
70 // add contributions from volume flux sources
71 2505200 source += problem.source(element, fvGeometry, elemVolVars, scv);
72
73 // add contribution from possible point sources
74 2505200 source += problem.scvPointSources(element, fvGeometry, elemVolVars, scv);
75
76 // Call the (kinetic) Energy module, for the source term.
77 // it has to be called from here, because the mass transferred has to be known.
78 if constexpr(ModelTraits::enableThermalNonEquilibrium())
79 {
80 2505200 EnergyLocalResidual::computeSourceEnergy(source,
81 element,
82 fvGeometry,
83 elemVolVars,
84 scv);
85 }
86
87 2505200 return source;
88 }
89
90 };
91
92 /*!
93 * \brief The local residual for a model assuming chemical non-equilibrium
94 * and potentially thermal non-equilibrium
95 */
96 template<class TypeTag>
97 class NonEquilibriumLocalResidualImplementation<TypeTag, true>: public GetPropType<TypeTag, Properties::EquilibriumLocalResidual>
98 {
99 using ParentType = GetPropType<TypeTag, Properties::EquilibriumLocalResidual>;
100 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
101 using Problem = GetPropType<TypeTag, Properties::Problem>;
102 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
103 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
104 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
105 using NumEqVector = Dumux::NumEqVector<GetPropType<TypeTag, Properties::PrimaryVariables>>;
106 using FluxVariables = GetPropType<TypeTag, Properties::FluxVariables>;
107 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
108 using GridView = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView;
109 using Element = typename GridView::template Codim<0>::Entity;
110 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
111 using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>;
112 using EnergyLocalResidual = GetPropType<TypeTag, Properties::EnergyLocalResidual>;
113 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
114
115 using ModelTraits = GetPropType<TypeTag, Properties::ModelTraits>;
116 using Indices = typename ModelTraits::Indices;
117
118 static constexpr int numPhases = ModelTraits::numFluidPhases();
119 static constexpr int numComponents = ModelTraits::numFluidComponents();
120
121 static constexpr auto conti0EqIdx = Indices::conti0EqIdx;
122 static constexpr auto comp0Idx = FluidSystem::comp0Idx;
123 static constexpr auto phase0Idx = FluidSystem::phase0Idx;
124 static constexpr auto phase1Idx = FluidSystem::phase1Idx;
125
126 static_assert(numPhases > 1,
127 "chemical non-equlibrium only makes sense for multiple phases");
128 static_assert(numPhases == numComponents,
129 "currently chemical non-equilibrium is only available when numPhases equals numComponents");
130 static_assert(ModelTraits::useMoles(),
131 "chemical nonequilibrium can only be calculated based on mole fractions not mass fractions");
132
133 public:
134
2/4
✓ Branch 1 taken 29760 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 29760 times.
✗ Branch 5 not taken.
59520 using ParentType::ParentType;
135 /*!
136 * \brief Calculates the storage for all mass balance equations.
137 *
138 * \param problem The object specifying the problem which ought to be simulated
139 * \param scv The sub-control volume
140 * \param volVars The volume variables
141 */
142 5981760 NumEqVector computeStorage(const Problem& problem,
143 const SubControlVolume& scv,
144 const VolumeVariables& volVars) const
145 {
146 5981760 NumEqVector storage(0.0);
147
148 // compute storage term of all components within all phases
149
2/2
✓ Branch 0 taken 11963520 times.
✓ Branch 1 taken 5981760 times.
17945280 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
150 {
151
2/2
✓ Branch 0 taken 23927040 times.
✓ Branch 1 taken 11963520 times.
35890560 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
152 {
153 23927040 auto eqIdx = conti0EqIdx + phaseIdx*numComponents + compIdx;
154 23927040 storage[eqIdx] += volVars.porosity()
155 47854080 * volVars.saturation(phaseIdx)
156 23927040 * volVars.molarDensity(phaseIdx)
157 47854080 * volVars.moleFraction(phaseIdx, compIdx);
158 }
159 //! The energy storage in the fluid phase with index phaseIdx
160 23927040 EnergyLocalResidual::fluidPhaseStorage(storage, problem, scv, volVars, phaseIdx);
161 }
162 //! The energy storage in the solid matrix
163 7239360 EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars);
164 5981760 return storage;
165 }
166
167 /*!
168 * \brief Calculates the flux for all mass balance equations.
169 *
170 * \param problem The object specifying the problem which ought to be simulated
171 * \param element An element which contains part of the control volume
172 * \param fvGeometry The finite-volume geometry
173 * \param elemVolVars The volume variables of the current element
174 * \param scvf The sub control volume face to compute the flux on
175 * \param elemFluxVarsCache The cache related to flux computation
176 */
177 3254040 NumEqVector computeFlux(const Problem& problem,
178 const Element& element,
179 const FVElementGeometry& fvGeometry,
180 const ElementVolumeVariables& elemVolVars,
181 const SubControlVolumeFace& scvf,
182 const ElementFluxVariablesCache& elemFluxVarsCache) const
183 {
184 3254040 FluxVariables fluxVars;
185 3254040 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
186 3254040 NumEqVector flux(0.0);
187
188 // advective fluxes
189
2/2
✓ Branch 0 taken 6508080 times.
✓ Branch 1 taken 3254040 times.
9762120 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
190 {
191 6508080 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
192
3/3
✓ Branch 0 taken 3567840 times.
✓ Branch 1 taken 11232240 times.
✓ Branch 2 taken 4724160 times.
19524240 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
193 {
194 // get equation index
195 13016160 const auto eqIdx = conti0EqIdx + phaseIdx*numComponents + compIdx;
196
197 // the physical quantities for which we perform upwinding
198 72216480 const auto upwindTerm = [phaseIdx, compIdx] (const auto& volVars)
199 104129280 { return volVars.molarDensity(phaseIdx)*volVars.moleFraction(phaseIdx, compIdx)*volVars.mobility(phaseIdx); };
200
201
2/2
✓ Branch 1 taken 6508080 times.
✓ Branch 2 taken 6508080 times.
13016160 flux[eqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
202
203 // do not add diffusive flux of main component, as that is not done in master as well
204
2/2
✓ Branch 0 taken 6508080 times.
✓ Branch 1 taken 6508080 times.
13016160 if (compIdx == phaseIdx)
205 6508080 continue;
206 //check for the reference system and adapt units of the diffusive flux accordingly.
207 if (FluxVariables::MolecularDiffusionType::referenceSystemFormulation() == ReferenceSystemFormulation::massAveraged)
208 22464480 flux[eqIdx] += diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx);
209 else
210 flux[eqIdx] += diffusiveFluxes[compIdx];
211 }
212 //! Add advective phase energy fluxes. For isothermal model the contribution is zero.
213 6508080 EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, phaseIdx);
214 }
215
216 //! Add diffusive energy fluxes. For isothermal model the contribution is zero.
217 3254040 EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
218
219 3254040 return flux;
220 }
221
222 /*!
223 * \brief Calculates the source term of the equation.
224 *
225 * \param problem The source term
226 * \param element An element which contains part of the control volume
227 * \param fvGeometry The finite-volume geometry
228 * \param elemVolVars The volume variables of the current element
229 * \param scv The sub-control volume over which we integrate the source term
230 */
231 2990880 NumEqVector computeSource(const Problem& problem,
232 const Element& element,
233 const FVElementGeometry& fvGeometry,
234 const ElementVolumeVariables& elemVolVars,
235 const SubControlVolume &scv) const
236 {
237 2990880 NumEqVector source(0.0);
238 // In the case of a kinetic consideration, mass transfer
239 // between phases is realized via source terms there is a
240 // balance equation for each component in each phase
241 2990880 const auto& volVars = elemVolVars[scv];
242 2990880 std::array<std::array<Scalar, numComponents>, numPhases> componentIntoPhaseMassTransfer = {{{0.0},{0.0}}};
243
244 //get characteristic length
245 2990880 const Scalar characteristicLength = volVars.characteristicLength() ;
246 2990880 const Scalar factorMassTransfer = volVars.factorMassTransfer() ;
247
248 2990880 const Scalar awn = volVars.interfacialArea(phase0Idx, phase1Idx);
249
250
2/2
✓ Branch 0 taken 5981760 times.
✓ Branch 1 taken 2990880 times.
8972640 for(int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
251 {
252 5981760 const Scalar sherwoodNumber = volVars.sherwoodNumber(phaseIdx);
253
254
2/2
✓ Branch 0 taken 11963520 times.
✓ Branch 1 taken 5981760 times.
17945280 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
255 {
256 if (compIdx <= numPhases)
257 {
258 //we only have to do the source calculation for one component going into the other phase as for each component: n_wn -> n_n = - n_n -> n_wn
259
2/2
✓ Branch 0 taken 5981760 times.
✓ Branch 1 taken 5981760 times.
11963520 if (phaseIdx == compIdx)
260 continue;
261
262 5981760 const Scalar xNonEquil = volVars.moleFraction(phaseIdx, compIdx);
263
264 //additionally get equilibrium values from volume variables
265 5981760 const Scalar xEquil = volVars.xEquil(phaseIdx, compIdx);
266 //get the diffusion coefficient
267 5981760 const Scalar diffCoeff = volVars.diffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
268
269 //now compute the flux
270 5981760 const Scalar compFluxIntoOtherPhase = factorMassTransfer * (xEquil-xNonEquil)/characteristicLength * awn * volVars.molarDensity(phaseIdx) * diffCoeff * sherwoodNumber;
271
272 11963520 componentIntoPhaseMassTransfer[phaseIdx][compIdx] += compFluxIntoOtherPhase;
273 17945280 componentIntoPhaseMassTransfer[compIdx][compIdx] -= compFluxIntoOtherPhase;
274 }
275 }
276 }
277
278 // Actually add the mass transfer to the sources which might
279 // exist externally
280
2/2
✓ Branch 0 taken 5981760 times.
✓ Branch 1 taken 2990880 times.
8972640 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
281 {
282
2/2
✓ Branch 0 taken 11963520 times.
✓ Branch 1 taken 5981760 times.
17945280 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
283 {
284 11963520 const unsigned int eqIdx = conti0EqIdx + compIdx + phaseIdx*numComponents;
285
3/6
✗ Branch 0 not taken.
✓ Branch 1 taken 11963520 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 11963520 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 11963520 times.
35890560 source[eqIdx] += componentIntoPhaseMassTransfer[phaseIdx][compIdx];
286
287 using std::isfinite;
288
3/6
✗ Branch 0 not taken.
✓ Branch 1 taken 11963520 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 11963520 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 11963520 times.
35890560 if (!isfinite(source[eqIdx]))
289 DUNE_THROW(NumericalProblem, "Calculated non-finite source");
290 }
291 }
292
293 if constexpr (ModelTraits::enableThermalNonEquilibrium())
294 {
295 // Call the (kinetic) Energy module, for the source term.
296 // it has to be called from here, because the mass transferred has to be known.
297 2362080 EnergyLocalResidual::computeSourceEnergy(source,
298 element,
299 fvGeometry,
300 elemVolVars,
301 scv);
302 }
303
304 // add contributions from volume flux sources
305 2990880 source += problem.source(element, fvGeometry, elemVolVars, scv);
306
307 // add contribution from possible point sources
308 2990880 source += problem.scvPointSources(element, fvGeometry, elemVolVars, scv);
309
310 2990880 return source;
311 }
312 };
313
314 } // end namespace Dumux
315
316 #endif
317