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 OnePNCModel | ||
10 | * \brief Quantities required by the single-phase, n-component box | ||
11 | * model defined on a vertex. | ||
12 | */ | ||
13 | |||
14 | #ifndef DUMUX_1PNC_VOLUME_VARIABLES_HH | ||
15 | #define DUMUX_1PNC_VOLUME_VARIABLES_HH | ||
16 | |||
17 | #include <dune/common/fvector.hh> | ||
18 | |||
19 | #include <dumux/porousmediumflow/volumevariables.hh> | ||
20 | #include <dumux/porousmediumflow/nonisothermal/volumevariables.hh> | ||
21 | #include <dumux/material/solidstates/updatesolidvolumefractions.hh> | ||
22 | |||
23 | namespace Dumux { | ||
24 | |||
25 | /*! | ||
26 | * \ingroup OnePNCModel | ||
27 | * \brief Contains the quantities which are are constant within a | ||
28 | * finite volume in the one-phase, n-component model. | ||
29 | * | ||
30 | * \note The default value for the phase index given in the fluid property interfaces is not used, | ||
31 | * but is only here to enable calling these functions without handing in a phase index | ||
32 | * (as in a single-phasic context there is only one phase). | ||
33 | */ | ||
34 | template <class Traits> | ||
35 |
26/36✗ Branch 0 not taken.
✓ Branch 1 taken 49769 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 43138 times.
✓ Branch 4 taken 6631 times.
✓ Branch 5 taken 43138 times.
✓ Branch 6 taken 10802 times.
✓ Branch 7 taken 22373 times.
✓ Branch 8 taken 10802 times.
✓ Branch 9 taken 37102 times.
✓ Branch 10 taken 16082 times.
✓ Branch 11 taken 37102 times.
✓ Branch 12 taken 16440 times.
✓ Branch 13 taken 34120 times.
✓ Branch 14 taken 16440 times.
✓ Branch 15 taken 7480 times.
✓ Branch 16 taken 63906 times.
✓ Branch 17 taken 7480 times.
✓ Branch 18 taken 138754 times.
✓ Branch 19 taken 80072 times.
✓ Branch 20 taken 138754 times.
✗ Branch 21 not taken.
✓ Branch 22 taken 218826 times.
✗ Branch 23 not taken.
✓ Branch 25 taken 1088121 times.
✗ Branch 26 not taken.
✓ Branch 28 taken 1050259 times.
✗ Branch 29 not taken.
✓ Branch 31 taken 1050259 times.
✗ Branch 32 not taken.
✗ Branch 33 not taken.
✓ Branch 34 taken 15732 times.
✗ Branch 35 not taken.
✓ Branch 36 taken 15708 times.
✗ Branch 37 not taken.
✓ Branch 38 taken 15708 times.
|
19861833 | class OnePNCVolumeVariables |
36 | : public PorousMediumFlowVolumeVariables<Traits> | ||
37 | , public EnergyVolumeVariables<Traits, OnePNCVolumeVariables<Traits> > | ||
38 | { | ||
39 | using ParentType = PorousMediumFlowVolumeVariables<Traits>; | ||
40 | using EnergyVolVars = EnergyVolumeVariables<Traits, OnePNCVolumeVariables<Traits> >; | ||
41 | using Scalar = typename Traits::PrimaryVariables::value_type; | ||
42 | using PermeabilityType = typename Traits::PermeabilityType; | ||
43 | using EffDiffModel = typename Traits::EffectiveDiffusivityModel; | ||
44 | using Idx = typename Traits::ModelTraits::Indices; | ||
45 | static constexpr int numFluidComps = ParentType::numFluidComponents(); | ||
46 | using DiffusionCoefficients = typename Traits::DiffusionType::DiffusionCoefficientsContainer; | ||
47 | |||
48 | enum | ||
49 | { | ||
50 | // pressure primary variable index | ||
51 | pressureIdx = Idx::pressureIdx | ||
52 | }; | ||
53 | |||
54 | public: | ||
55 | //! Export fluid state type | ||
56 | using FluidState = typename Traits::FluidState; | ||
57 | //! Export fluid system type | ||
58 | using FluidSystem = typename Traits::FluidSystem; | ||
59 | //! Export indices | ||
60 | using Indices = typename Traits::ModelTraits::Indices; | ||
61 | //! Export type of solid state | ||
62 | using SolidState = typename Traits::SolidState; | ||
63 | //! Export type of solid system | ||
64 | using SolidSystem = typename Traits::SolidSystem; | ||
65 | |||
66 | //! Returns whether moles or masses are balanced | ||
67 | static constexpr bool useMoles() { return Traits::ModelTraits::useMoles(); } | ||
68 | |||
69 | /*! | ||
70 | * \brief Updates all quantities for a given control volume. | ||
71 | * | ||
72 | * \param elemSol A vector containing all primary variables connected to the element | ||
73 | * \param problem The object specifying the problem which ought to | ||
74 | * be simulated | ||
75 | * \param element An element which contains part of the control volume | ||
76 | * \param scv The sub-control volume | ||
77 | */ | ||
78 | template<class ElemSol, class Problem, class Element, class Scv> | ||
79 | 24121955 | void update(const ElemSol &elemSol, | |
80 | const Problem &problem, | ||
81 | const Element &element, | ||
82 | const Scv &scv) | ||
83 | { | ||
84 | 24121955 | ParentType::update(elemSol, problem, element, scv); | |
85 | |||
86 | 24121955 | completeFluidState(elemSol, problem, element, scv, fluidState_, solidState_); | |
87 | |||
88 | // calculate the remaining quantities | ||
89 | 24121955 | updateSolidVolumeFractions(elemSol, problem, element, scv, solidState_, numFluidComps); | |
90 | 24121955 | EnergyVolVars::updateSolidEnergyParams(elemSol, problem, element, scv, solidState_); | |
91 | 48219078 | permeability_ = problem.spatialParams().permeability(element, scv, elemSol); | |
92 | 24121955 | EnergyVolVars::updateEffectiveThermalConductivity(); | |
93 | |||
94 | 24121955 | auto getEffectiveDiffusionCoefficient = [&](int phaseIdx, int compIIdx, int compJIdx) | |
95 | { | ||
96 | 18336761 | return EffDiffModel::effectiveDiffusionCoefficient(*this, phaseIdx, compIIdx, compJIdx); | |
97 | }; | ||
98 | |||
99 | 24121955 | effectiveDiffCoeff_.update(getEffectiveDiffusionCoefficient); | |
100 | 24121955 | } | |
101 | |||
102 | /*! | ||
103 | * \brief Sets complete fluid state. | ||
104 | * | ||
105 | * \param elemSol A vector containing all primary variables connected to the element | ||
106 | * \param problem The object specifying the problem which ought to | ||
107 | * be simulated | ||
108 | * \param element An element which contains part of the control volume | ||
109 | * \param scv The sub-control volume | ||
110 | * \param fluidState A container with the current (physical) state of the fluid | ||
111 | * \param solidState A container with the current (physical) state of the solid | ||
112 | */ | ||
113 | template<class ElemSol, class Problem, class Element, class Scv> | ||
114 | 24121955 | void completeFluidState(const ElemSol &elemSol, | |
115 | const Problem& problem, | ||
116 | const Element& element, | ||
117 | const Scv &scv, | ||
118 | FluidState& fluidState, | ||
119 | SolidState& solidState) | ||
120 | { | ||
121 | 24121955 | EnergyVolVars::updateTemperature(elemSol, problem, element, scv, fluidState, solidState); | |
122 | 24121955 | fluidState.setSaturation(0, 1.0); | |
123 | |||
124 | 24121955 | const auto& priVars = elemSol[scv.localDofIndex()]; | |
125 | 48243910 | fluidState.setPressure(0, priVars[pressureIdx]); | |
126 | |||
127 | // Set fluid state mole fractions | ||
128 | if (useMoles()) | ||
129 | { | ||
130 | 23953155 | Scalar sumMoleFracNotMainComp = 0; | |
131 |
2/2✓ Branch 0 taken 17573657 times.
✓ Branch 1 taken 16979353 times.
|
48500614 | for (int compIdx = 1; compIdx < numFluidComps; ++compIdx) |
132 | { | ||
133 | 49094918 | fluidState.setMoleFraction(0, compIdx, priVars[compIdx]); | |
134 | 49094918 | sumMoleFracNotMainComp += priVars[compIdx]; | |
135 | } | ||
136 | 23953155 | fluidState.setMoleFraction(0, 0, 1.0 - sumMoleFracNotMainComp); | |
137 | } | ||
138 | else | ||
139 | { | ||
140 | // for mass fractions we only have to set numComponents-1 mass fractions | ||
141 | // the fluid state will internally compute the remaining mass fraction | ||
142 |
2/2✓ Branch 0 taken 168800 times.
✓ Branch 1 taken 168800 times.
|
337600 | for (int compIdx = 1; compIdx < numFluidComps; ++compIdx) |
143 | 337600 | fluidState.setMassFraction(0, compIdx, priVars[compIdx]); | |
144 | } | ||
145 | |||
146 | typename FluidSystem::ParameterCache paramCache; | ||
147 | 24121955 | paramCache.updateAll(fluidState); | |
148 | |||
149 | 24716259 | Scalar rho = FluidSystem::density(fluidState, paramCache, 0); | |
150 | 24895273 | Scalar rhoMolar = FluidSystem::molarDensity(fluidState, paramCache, 0); | |
151 | 24895273 | Scalar mu = FluidSystem::viscosity(fluidState, paramCache, 0); | |
152 | |||
153 | 33955905 | fluidState.setDensity(0, rho); | |
154 | 33955905 | fluidState.setMolarDensity(0, rhoMolar); | |
155 | 33955905 | fluidState.setViscosity(0, mu); | |
156 | |||
157 | // compute and set the enthalpy | ||
158 | 24121955 | Scalar h = EnergyVolVars::enthalpy(fluidState, paramCache, 0); | |
159 | 48243910 | fluidState.setEnthalpy(0, h); | |
160 | 24121955 | } | |
161 | |||
162 | /*! | ||
163 | * \brief Returns the fluid configuration at the given primary | ||
164 | * variables. | ||
165 | */ | ||
166 | const FluidState &fluidState() const | ||
167 | 123767974 | { return fluidState_; } | |
168 | |||
169 | /*! | ||
170 | * \brief Returns the phase state for the control volume. | ||
171 | */ | ||
172 | const SolidState &solidState() const | ||
173 |
2/4✓ Branch 1 taken 72 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 72 times.
✗ Branch 5 not taken.
|
69859414 | { return solidState_; } |
174 | |||
175 | /*! | ||
176 | * \brief Returns the average molar mass \f$\mathrm{[kg/mol]}\f$ of the fluid phase. | ||
177 | * | ||
178 | * \param phaseIdx The phase index | ||
179 | */ | ||
180 | ✗ | Scalar averageMolarMass(int phaseIdx = 0) const | |
181 | 11343104 | { return fluidState_.averageMolarMass(0); } | |
182 | |||
183 | /*! | ||
184 | * \brief Returns density \f$\mathrm{[kg/m^3]}\f$ the of the fluid phase. | ||
185 | * | ||
186 | * \note the phase index passed to this function is for compatibility reasons | ||
187 | * with multi-phasic models. | ||
188 | */ | ||
189 | ✗ | Scalar density(int phaseIdx = 0) const | |
190 | { | ||
191 |
12/16✓ Branch 0 taken 2553460 times.
✓ Branch 1 taken 118670 times.
✓ Branch 2 taken 2553460 times.
✓ Branch 3 taken 227069 times.
✓ Branch 4 taken 2553508 times.
✓ Branch 5 taken 159894 times.
✓ Branch 6 taken 2553547 times.
✓ Branch 7 taken 51553 times.
✓ Branch 8 taken 72468 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 72487 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 72444 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 72444 times.
✗ Branch 15 not taken.
|
678363115 | return fluidState_.density(0); |
192 | } | ||
193 | |||
194 | /*! | ||
195 | * \brief Returns molar density \f$\mathrm{[mol/m^3]}\f$ the of the fluid phase. | ||
196 | * | ||
197 | * \note the phase index passed to this function is for compatibility reasons | ||
198 | * with multi-phasic models. | ||
199 | */ | ||
200 | ✗ | Scalar molarDensity(int phaseIdx = 0) const | |
201 | { | ||
202 |
10/24✗ Branch 0 not taken.
✓ Branch 1 taken 87963764 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 75960924 times.
✓ Branch 4 taken 227175 times.
✓ Branch 5 taken 31585385 times.
✓ Branch 6 taken 227175 times.
✓ Branch 7 taken 31585385 times.
✓ Branch 8 taken 43616 times.
✓ Branch 9 taken 76978264 times.
✓ Branch 10 taken 43616 times.
✓ Branch 11 taken 76978264 times.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
|
569234315 | return fluidState_.molarDensity(0); |
203 | } | ||
204 | |||
205 | /*! | ||
206 | * \brief Returns the saturation. | ||
207 | * | ||
208 | * This method is here for compatibility reasons with other models. The saturation | ||
209 | * is always 1.0 in a one-phasic context. | ||
210 | */ | ||
211 | ✗ | Scalar saturation(int phaseIdx = 0) const | |
212 | ✗ | { return 1.0; } | |
213 | |||
214 | /*! | ||
215 | * \brief Returns the mole fraction \f$\mathrm{[mol/mol]}\f$ of a component in the phase. | ||
216 | * | ||
217 | * \param phaseIdx The index of the fluid phase | ||
218 | * \param compIdx The index of the component | ||
219 | * | ||
220 | * \note The phase index passed to this function is for compatibility reasons | ||
221 | * with multi-phasic models. | ||
222 | */ | ||
223 | ✗ | Scalar moleFraction(int phaseIdx, int compIdx) const | |
224 | { | ||
225 | // make sure this is only called with admissible indices | ||
226 |
10/34✗ Branch 0 not taken.
✓ Branch 1 taken 87963764 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 31446904 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 76920600 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 10700 times.
✗ Branch 15 not taken.
✓ Branch 16 taken 745600 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 36000 times.
✗ Branch 27 not taken.
✓ Branch 28 taken 61600 times.
✗ Branch 30 not taken.
✗ Branch 31 not taken.
✗ Branch 33 not taken.
✗ Branch 34 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 224 times.
✗ Branch 39 not taken.
✓ Branch 40 taken 12432 times.
✗ Branch 42 not taken.
✓ Branch 43 taken 69888 times.
✗ Branch 45 not taken.
✗ Branch 46 not taken.
✗ Branch 48 not taken.
✗ Branch 49 not taken.
|
197267712 | assert(compIdx < numFluidComps); |
227 |
12/12✓ Branch 0 taken 199308 times.
✓ Branch 1 taken 85828 times.
✓ Branch 2 taken 199308 times.
✓ Branch 3 taken 85828 times.
✓ Branch 4 taken 119393 times.
✓ Branch 5 taken 182439 times.
✓ Branch 6 taken 119393 times.
✓ Branch 7 taken 182439 times.
✓ Branch 16 taken 20506 times.
✓ Branch 17 taken 23494 times.
✓ Branch 18 taken 20506 times.
✓ Branch 19 taken 23494 times.
|
262501644 | return fluidState_.moleFraction(0, compIdx); |
228 | } | ||
229 | |||
230 | /*! | ||
231 | * \brief Returns the mass fraction of a component in the phase. | ||
232 | * | ||
233 | * \param phaseIdx The index of the fluid phase | ||
234 | * \param compIdx The index of the component | ||
235 | * | ||
236 | * \note The phase index passed to this function is for compatibility reasons | ||
237 | * with multi-phasic models. | ||
238 | */ | ||
239 | ✗ | Scalar massFraction(int phaseIdx, int compIdx) const | |
240 | { | ||
241 | // make sure this is only called with admissible indices | ||
242 | ✗ | assert(compIdx < numFluidComps); | |
243 |
1/4✗ Branch 0 not taken.
✗ Branch 1 not taken.
✓ Branch 3 taken 17380 times.
✗ Branch 4 not taken.
|
17380 | return fluidState_.massFraction(0, compIdx); |
244 | } | ||
245 | |||
246 | /*! | ||
247 | * \brief Returns the effective pressure \f$\mathrm{[Pa]}\f$ of a given phase within | ||
248 | * the control volume. | ||
249 | * | ||
250 | * \param phaseIdx The phase index | ||
251 | * | ||
252 | * \note The phase index passed to this function is for compatibility reasons | ||
253 | * with multi-phasic models. | ||
254 | */ | ||
255 | ✗ | Scalar pressure(int phaseIdx = 0) const | |
256 | { | ||
257 |
4/8✗ Branch 0 not taken.
✓ Branch 1 taken 503998 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 503998 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 328526 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 328526 times.
|
301969203 | return fluidState_.pressure(0); |
258 | } | ||
259 | |||
260 | /*! | ||
261 | * \brief Returns the temperature \f$\mathrm{[K]}\f$ inside the sub-control volume. | ||
262 | * | ||
263 | * Note that we assume thermodynamic equilibrium, i.e. the | ||
264 | * temperature of the rock matrix and of all fluid phases are | ||
265 | * identical. | ||
266 | */ | ||
267 | Scalar temperature() const | ||
268 |
10/20✓ Branch 0 taken 3249336 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 3249336 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 22464 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 22464 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 36288 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 36288 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 36288 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 36288 times.
✗ Branch 15 not taken.
✓ Branch 16 taken 36288 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 36288 times.
✗ Branch 19 not taken.
|
182308528 | { return fluidState_.temperature(); } |
269 | |||
270 | /*! | ||
271 | * \brief Returns the mobility \f$\mathrm{[1/(Pa s)]}\f$. | ||
272 | * | ||
273 | * The term mobility is usually not employed in the one phase context. | ||
274 | * The method is here for compatibility reasons with other models. | ||
275 | * | ||
276 | * \note The phase index passed to this function is for compatibility reasons | ||
277 | * with multi-phasic models. | ||
278 | */ | ||
279 | ✗ | Scalar mobility(int phaseIdx = 0) const | |
280 | { | ||
281 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 175472 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 175472 times.
|
205003514 | return 1.0/fluidState_.viscosity(0); |
282 | } | ||
283 | |||
284 | /*! | ||
285 | * \brief Returns the dynamic viscosity \f$\mathrm{[Pa s]}\f$ of the fluid within the | ||
286 | * control volume. | ||
287 | * | ||
288 | * \note The phase index passed to this function is for compatibility reasons | ||
289 | * with multi-phasic models. | ||
290 | */ | ||
291 | ✗ | Scalar viscosity(int phaseIdx = 0) const | |
292 | { | ||
293 | 30363344 | return fluidState_.viscosity(0); | |
294 | } | ||
295 | |||
296 | /*! | ||
297 | * \brief Returns the average porosity \f$\mathrm{[-]}\f$ within the control volume. | ||
298 | */ | ||
299 | Scalar porosity() const | ||
300 |
1/2✓ Branch 2 taken 1216 times.
✗ Branch 3 not taken.
|
273168092 | { return solidState_.porosity(); } |
301 | |||
302 | /*! | ||
303 | * \brief Returns the binary diffusion coefficients for a phase in \f$[m^2/s]\f$. | ||
304 | */ | ||
305 | ✗ | Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const | |
306 | { | ||
307 | typename FluidSystem::ParameterCache paramCache; | ||
308 | 18737805 | paramCache.updatePhase(fluidState_, phaseIdx); | |
309 |
1/4✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 290029 times.
✗ Branch 4 not taken.
|
18472608 | return FluidSystem::binaryDiffusionCoefficient(fluidState_, paramCache, phaseIdx, compIIdx, compJIdx); |
310 | } | ||
311 | |||
312 | /*! | ||
313 | * \brief Returns the effective diffusion coefficients for a phase in \f$[m^2/s]\f$. | ||
314 | */ | ||
315 | Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const | ||
316 | 92740600 | { return effectiveDiffCoeff_(phaseIdx, compIIdx, compJIdx); } | |
317 | |||
318 | /*! | ||
319 | * \brief Returns the molarity of a component in the phase. | ||
320 | * | ||
321 | * \param compIdx The index of the component | ||
322 | */ | ||
323 | Scalar molarity(int compIdx) const // [moles/m^3] | ||
324 | { | ||
325 | assert(compIdx < numFluidComps); | ||
326 | return fluidState_.molarity(0, compIdx); | ||
327 | } | ||
328 | |||
329 | /*! | ||
330 | * \brief Returns the mass fraction of a component in the phase. | ||
331 | * | ||
332 | * \param compIdx The index of the component | ||
333 | */ | ||
334 | Scalar massFraction(int compIdx) const | ||
335 | { | ||
336 | assert(compIdx < numFluidComps); | ||
337 | return fluidState_.massFraction(0, compIdx); | ||
338 | } | ||
339 | |||
340 | /*! | ||
341 | * \brief Returns the permeability within the control volume in \f$[m^2]\f$. | ||
342 | */ | ||
343 | const PermeabilityType& permeability() const | ||
344 | ✗ | { return permeability_; } | |
345 | |||
346 | protected: | ||
347 | FluidState fluidState_; | ||
348 | SolidState solidState_; | ||
349 | |||
350 | private: | ||
351 | PermeabilityType permeability_; | ||
352 | |||
353 | // Effective diffusion coefficients for the phases | ||
354 | DiffusionCoefficients effectiveDiffCoeff_; | ||
355 | }; | ||
356 | |||
357 | } // end namespace Dumux | ||
358 | |||
359 | #endif | ||
360 |