GCC Code Coverage Report


Directory: ../../../builds/dumux-repositories/
File: /builds/dumux-repositories/dumux/dumux/material/fluidsystems/3pimmiscible.hh
Date: 2024-09-21 20:52:54
Exec Total Coverage
Lines: 74 90 82.2%
Functions: 15 17 88.2%
Branches: 61 251 24.3%

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 FluidSystems
10 * \copybrief Dumux::FluidSystems::ThreePImmiscible
11 */
12 #ifndef DUMUX_3P_IMMISCIBLE_FLUID_SYSTEM_HH
13 #define DUMUX_3P_IMMISCIBLE_FLUID_SYSTEM_HH
14
15 #include <cassert>
16 #include <limits>
17 #include <iostream>
18
19 #include <dune/common/exceptions.hh>
20
21 #include <dumux/material/fluidsystems/1pliquid.hh>
22 #include <dumux/material/fluidsystems/1pgas.hh>
23 #include <dumux/material/fluidstates/immiscible.hh>
24 #include <dumux/material/components/base.hh>
25 #include <dumux/io/name.hh>
26
27 namespace Dumux {
28 namespace FluidSystems {
29
30 /*!
31 * \ingroup FluidSystems
32 * \brief A fluid system for three-phase models assuming immiscibility and
33 * thermodynamic equilibrium
34 *
35 * The fluid phases are completely specified by means of their
36 * constituting components.
37 * The wetting and the nonwetting phase can be defined individually
38 * via FluidSystem::OnePLiquid<Scalar, Component>. The gas phase can be defined via
39 * FluidSystems::OnePGas<Scalar, Component>
40 * These phases consist of one pure component.
41 *
42 * \tparam Scalar the scalar type
43 * \tparam WettingFluid the wetting phase fluid system (use FluidSystem::OnePLiquid<Scalar, Component>)
44 * \tparam NonwettingFluid the wetting phase fluid system (use FluidSystem::OnePLiquid<Scalar, Component>)
45 * \tparam Gas the gas phase fluid system (use FluidSystem::OnePGas<Scalar, Component>)
46 */
47 template <class Scalar, class WettingFluid, class NonwettingFluid, class Gas>
48 class ThreePImmiscible
49 : public Base<Scalar, ThreePImmiscible<Scalar, WettingFluid, NonwettingFluid, Gas> >
50 {
51 static_assert((WettingFluid::numPhases == 1), "WettingFluid has more than one phase");
52 static_assert((NonwettingFluid::numPhases == 1), "NonwettingFluid has more than one phase");
53 static_assert((Gas::numPhases == 1), "Gas has more than one phase");
54 static_assert((WettingFluid::numComponents == 1), "WettingFluid has more than one component");
55 static_assert((NonwettingFluid::numComponents == 1), "NonwettingFluid has more than one component");
56 static_assert((Gas::numComponents == 1), "Gas has more than one component");
57
58 using ThisType = ThreePImmiscible<Scalar, WettingFluid, NonwettingFluid, Gas>;
59
60 public:
61 /****************************************
62 * Fluid phase related static parameters
63 ****************************************/
64
65 //! Number of phases in the fluid system
66 static constexpr int numPhases = 3;
67
68 //! Index of the wetting phase
69 static constexpr int wPhaseIdx = 0;
70 //! Index of the nonwetting phase
71 static constexpr int nPhaseIdx = 1;
72 //! Index of the gas phase
73 static constexpr int gPhaseIdx = 2;
74
75 /*!
76 * \brief Return the human readable name of a fluid phase
77 * \param phaseIdx The index of the fluid phase to consider
78 */
79 21 static std::string phaseName(int phaseIdx)
80 {
81
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 21 times.
21 assert(0 <= phaseIdx && phaseIdx < numPhases);
82
3/4
✓ Branch 0 taken 7 times.
✓ Branch 1 taken 7 times.
✓ Branch 2 taken 7 times.
✗ Branch 3 not taken.
21 switch (phaseIdx)
83 {
84 case wPhaseIdx: return Components::IsAqueous<typename WettingFluid::Component>::value
85 7 ? IOName::aqueousPhase() : IOName::naplPhase();
86 case nPhaseIdx: return Components::IsAqueous<typename NonwettingFluid::Component>::value
87 7 ? IOName::aqueousPhase() : IOName::naplPhase();
88 7 case gPhaseIdx: return IOName::gaseousPhase();
89 }
90 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
91 }
92
93 /*!
94 * \brief Returns whether the fluids are miscible
95 */
96 static constexpr bool isMiscible()
97 { return false; }
98
99 /*!
100 * \brief Return whether a phase is gaseous
101 * \param phaseIdx The index of the fluid phase to consider
102 */
103 static constexpr bool isGas(int phaseIdx)
104 {
105 assert(0 <= phaseIdx && phaseIdx < numPhases);
106
107 switch (phaseIdx)
108 {
109 case wPhaseIdx: return WettingFluid::isGas(); break;
110 case nPhaseIdx: return NonwettingFluid::isGas(); break;
111 case gPhaseIdx: return Gas::isGas(); break;
112 default: return false; // TODO: constexpr-compatible throw
113 }
114 }
115
116 /*!
117 * \brief Returns true if and only if a fluid phase is assumed to
118 * be an ideal mixture.
119 * \param phaseIdx The index of the fluid phase to consider
120 *
121 * We define an ideal mixture as a fluid phase where the fugacity
122 * coefficients of all components times the pressure of the phase
123 * are independent on the fluid composition. This assumption is true
124 * if immiscibility is assumed. If you are unsure what
125 * this function should return, it is safe to return false. The
126 * only damage done will be (slightly) increased computation times
127 * in some cases.
128 */
129 static constexpr bool isIdealMixture(int phaseIdx)
130 { return true; }
131
132 /*!
133 * \brief Returns true if and only if a fluid phase is assumed to
134 * be compressible.
135 *
136 * Compressible means. that the partial derivative of the density
137 * to the fluid pressure is always larger than zero.
138 *
139 * \param phaseIdx The index of the fluid phase to consider
140 */
141 static constexpr bool isCompressible(int phaseIdx)
142 {
143 assert(0 <= phaseIdx && phaseIdx < numPhases);
144
145 // let the fluids decide
146 switch(phaseIdx)
147 {
148 case wPhaseIdx: return WettingFluid::isCompressible(); break;
149 case nPhaseIdx: return NonwettingFluid::isCompressible(); break;
150 case gPhaseIdx: return Gas::isCompressible(); break;
151 default: return false; // TODO: constexpr-compatible throw
152 }
153 }
154
155 /*!
156 * \brief Returns true if and only if a fluid phase is assumed to
157 * be an ideal gas.
158 *
159 * \param phaseIdx The index of the fluid phase to consider
160 */
161 static constexpr bool isIdealGas(int phaseIdx)
162 {
163 assert(0 <= phaseIdx && phaseIdx < numPhases);
164
165 // let the fluids decide
166 switch(phaseIdx)
167 {
168 case wPhaseIdx: return WettingFluid::isIdealGas(); break;
169 case nPhaseIdx: return NonwettingFluid::isIdealGas(); break;
170 case gPhaseIdx: return Gas::isIdealGas(); break;
171 default: return false; // TODO: constexpr-compatible throw
172 }
173 }
174
175 /****************************************
176 * Component related static parameters
177 ****************************************/
178
179 //! Number of components in the fluid system
180 static constexpr int numComponents = 3;//WettingFluid::numComponents + NonwettingFluid::numComponents + Gas::numComponents; // TODO: 3??
181
182 //! Index of the wetting phase's component
183 static constexpr int wCompIdx = 0;
184 //! Index of the nonwetting phase's component
185 static constexpr int nCompIdx = 1;
186 //! Index of the gas phase's component
187 static constexpr int gCompIdx = 2; // TODO: correct??
188
189 /*!
190 * \brief Return the human readable name of a component
191 *
192 * \param compIdx index of the component
193 */
194 3 static std::string componentName(int compIdx)
195 {
196
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 3 times.
3 assert(0 <= compIdx && compIdx < numComponents);
197
198
3/4
✓ Branch 0 taken 1 times.
✓ Branch 1 taken 1 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
3 switch(compIdx)
199 {
200 1 case wCompIdx: return WettingFluid::name(); break;
201 1 case nCompIdx: return NonwettingFluid::name(); break;
202 1 case gCompIdx: return Gas::name(); break;
203 default: DUNE_THROW(Dune::InvalidStateException, "Invalid component index");
204 }
205 }
206
207 /*!
208 * \brief Return the molar mass of a component in \f$\mathrm{[kg/mol]}\f$.
209 * \param compIdx index of the component
210 */
211 30 static Scalar molarMass(int compIdx)
212 {
213
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 30 times.
30 assert(0 <= compIdx && compIdx < numComponents);
214
215
2/3
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 20 times.
30 switch(compIdx)
216 {
217 case wCompIdx: return WettingFluid::molarMass(); break;
218 case nCompIdx: return NonwettingFluid::molarMass(); break;
219 10 case gCompIdx: return Gas::molarMass(); break;
220 default: DUNE_THROW(Dune::InvalidStateException, "Invalid component index");
221 }
222 }
223
224 /*!
225 * \brief Critical temperature of a component \f$\mathrm{[K]}\f$.
226 * \param compIdx index of the component
227 */
228 static Scalar criticalTemperature(int compIdx)
229 {
230 assert(0 <= compIdx && compIdx < numComponents);
231
232 switch(compIdx)
233 {
234 case wCompIdx: return WettingFluid::criticalTemperature(); break;
235 case nCompIdx: return NonwettingFluid::criticalTemperature(); break;
236 case gCompIdx: return Gas::criticalTemperature(); break;
237 default: DUNE_THROW(Dune::InvalidStateException, "Invalid component index");
238 }
239 }
240
241 /*!
242 * \brief Critical pressure of a component \f$\mathrm{[Pa]}\f$.
243 * \param compIdx index of the component
244 */
245 static Scalar criticalPressure(int compIdx)
246 {
247 assert(0 <= compIdx && compIdx < numComponents);
248
249 switch(compIdx)
250 {
251 case wCompIdx: return WettingFluid::criticalPressure(); break;
252 case nCompIdx: return NonwettingFluid::criticalPressure(); break;
253 case gCompIdx: return Gas::criticalPressure(); break;
254 default: DUNE_THROW(Dune::InvalidStateException, "Invalid component index");
255 }
256 }
257
258 /*!
259 * \brief The acentric factor of a component \f$\mathrm{[-]}\f$.
260 * \param compIdx index of the component
261 */
262 static Scalar acentricFactor(int compIdx)
263 {
264 assert(0 <= compIdx && compIdx < numComponents);
265
266 switch(compIdx)
267 {
268 case wCompIdx: return WettingFluid::Component::acentricFactor(); break;
269 case nCompIdx: return NonwettingFluid::Component::acentricFactor(); break;
270 case gCompIdx: return Gas::Component::acentricFactor(); break;
271 default: DUNE_THROW(Dune::InvalidStateException, "Invalid component index");
272 }
273 }
274
275 /****************************************
276 * thermodynamic relations
277 ****************************************/
278
279 /*!
280 * \brief Initialize the fluid system's static parameters
281 */
282 static constexpr void init()
283 {
284 // two gaseous phases at once do not make sense physically!
285 // (But two liquids are fine)
286 static_assert(!WettingFluid::isGas() && !NonwettingFluid::isGas() && Gas::isGas(), "There can only be one gaseous phase!");
287 }
288
289 /*!
290 * \brief Initialize the fluid system's static parameters using
291 * problem specific temperature and pressure ranges
292 *
293 * \param tempMin The minimum temperature used for tabulation of water \f$\mathrm{[K]}\f$
294 * \param tempMax The maximum temperature used for tabulation of water \f$\mathrm{[K]}\f$
295 * \param nTemp The number of ticks on the temperature axis of the table of water
296 * \param pressMin The minimum pressure used for tabulation of water \f$\mathrm{[Pa]}\f$
297 * \param pressMax The maximum pressure used for tabulation of water \f$\mathrm{[Pa]}\f$
298 * \param nPress The number of ticks on the pressure axis of the table of water
299 */
300 2 static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp,
301 Scalar pressMin, Scalar pressMax, unsigned nPress)
302 {
303 // two gaseous phases at once do not make sense physically!
304 static_assert(!WettingFluid::isGas() && !NonwettingFluid::isGas() && Gas::isGas(), "There can only be one gaseous phase!");
305
306 if (WettingFluid::Component::isTabulated)
307 {
308 2 std::cout << "Initializing tables for the wetting fluid properties ("
309 2 << nTemp*nPress
310 2 << " entries).\n";
311
312 2 WettingFluid::Component::init(tempMin, tempMax, nTemp,
313 pressMin, pressMax, nPress);
314
315 }
316
317 if (NonwettingFluid::Component::isTabulated)
318 {
319 std::cout << "Initializing tables for the nonwetting fluid properties ("
320 << nTemp*nPress
321 << " entries).\n";
322
323 NonwettingFluid::Component::init(tempMin, tempMax, nTemp,
324 pressMin, pressMax, nPress);
325
326 }
327
328 if (Gas::Component::isTabulated)
329 {
330 std::cout << "Initializing tables for the gas fluid properties ("
331 << nTemp*nPress
332 << " entries).\n";
333
334 Gas::Component::init(tempMin, tempMax, nTemp,
335 pressMin, pressMax, nPress);
336
337 }
338 2 }
339
340 using Base<Scalar, ThisType>::density;
341 //! \copydoc Base<Scalar,ThisType>::density(const FluidState&,int)
342 template <class FluidState>
343 758703 static Scalar density(const FluidState &fluidState,
344 int phaseIdx)
345 {
346
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 758703 times.
758703 assert(0 <= phaseIdx && phaseIdx < numPhases);
347
348
3/4
✓ Branch 0 taken 252900 times.
✓ Branch 1 taken 252900 times.
✓ Branch 2 taken 252900 times.
✗ Branch 3 not taken.
758703 Scalar temperature = fluidState.temperature(phaseIdx);
349
3/4
✓ Branch 0 taken 252900 times.
✓ Branch 1 taken 252900 times.
✓ Branch 2 taken 252900 times.
✗ Branch 3 not taken.
758703 Scalar pressure = fluidState.pressure(phaseIdx);
350
351
3/4
✓ Branch 0 taken 252901 times.
✓ Branch 1 taken 252901 times.
✓ Branch 2 taken 252901 times.
✗ Branch 3 not taken.
758703 switch(phaseIdx)
352 {
353 505802 case wPhaseIdx: return WettingFluid::density(temperature, pressure); break;
354 252902 case nPhaseIdx: return NonwettingFluid::density(temperature, pressure); break;
355 505802 case gPhaseIdx: return Gas::density(temperature, pressure); break;
356 default: DUNE_THROW(Dune::InvalidStateException, "Invalid phase index");
357 }
358 }
359
360 using Base<Scalar, ThisType>::molarDensity;
361 /*!
362 * \brief The molar density \f$\rho_{mol,\alpha}\f$
363 * of a fluid phase \f$\alpha\f$ in \f$\mathrm{[mol/m^3]}\f$
364 *
365 * The molar density is defined by the
366 * mass density \f$\rho_\alpha\f$ and the component molar mass \f$M_\alpha\f$:
367 *
368 * \f[\rho_{mol,\alpha} = \frac{\rho_\alpha}{M_\alpha} \;.\f]
369 */
370 template <class FluidState>
371 3 static Scalar molarDensity(const FluidState &fluidState,
372 int phaseIdx)
373 {
374
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 3 times.
3 assert(0 <= phaseIdx && phaseIdx < numPhases);
375
376 3 Scalar temperature = fluidState.temperature(phaseIdx);
377 3 Scalar pressure = fluidState.pressure(phaseIdx);
378
379
3/4
✓ Branch 0 taken 1 times.
✓ Branch 1 taken 1 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
3 switch(phaseIdx)
380 {
381 1 case wPhaseIdx: return WettingFluid::molarDensity(temperature, pressure);
382 1 case nPhaseIdx: return NonwettingFluid::molarDensity(temperature, pressure);
383 2 case gPhaseIdx: return Gas::molarDensity(temperature, pressure);
384 default: DUNE_THROW(Dune::InvalidStateException, "Invalid phase index");
385 }
386 }
387
388 using Base<Scalar, ThisType>::viscosity;
389 /*!
390 * \brief Return the viscosity of a phase \f$\mathrm{[Pa*s]}\f$.
391 * \param fluidState The fluid state of the two-phase model
392 * \param phaseIdx Index of the fluid phase
393 */
394 template <class FluidState>
395 758703 static Scalar viscosity(const FluidState &fluidState,
396 int phaseIdx)
397 {
398
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 758703 times.
758703 assert(0 <= phaseIdx && phaseIdx < numPhases);
399
400
3/4
✓ Branch 0 taken 252900 times.
✓ Branch 1 taken 252900 times.
✓ Branch 2 taken 252900 times.
✗ Branch 3 not taken.
758703 Scalar temperature = fluidState.temperature(phaseIdx);
401
3/4
✓ Branch 0 taken 252900 times.
✓ Branch 1 taken 252900 times.
✓ Branch 2 taken 252900 times.
✗ Branch 3 not taken.
758703 Scalar pressure = fluidState.pressure(phaseIdx);
402
403
3/4
✓ Branch 0 taken 252901 times.
✓ Branch 1 taken 252901 times.
✓ Branch 2 taken 252901 times.
✗ Branch 3 not taken.
758703 switch(phaseIdx)
404 {
405 505802 case wPhaseIdx: return WettingFluid::viscosity(temperature, pressure); break;
406 505802 case nPhaseIdx: return NonwettingFluid::viscosity(temperature, pressure); break;
407 505802 case gPhaseIdx: return Gas::viscosity(temperature, pressure); break;
408 default: DUNE_THROW(Dune::InvalidStateException, "Invalid phase index");
409 }
410 }
411
412 using Base<Scalar, ThisType>::fugacityCoefficient;
413 /*!
414 * \brief Calculate the fugacity coefficient \f$\mathrm{[-]}\f$ of an individual
415 * component in a fluid phase
416 *
417 * The fugacity coefficient \f$\mathrm{\phi^\kappa_\alpha}\f$ is connected to the
418 * fugacity \f$\mathrm{f^\kappa_\alpha}\f$ and the component's mole
419 * fraction \f$\mathrm{x^\kappa_\alpha}\f$ by means of the relation
420 *
421 * \f[
422 f^\kappa_\alpha = \phi^\kappa_\alpha\;x^\kappa_\alpha\;p_\alpha
423 * \f]
424 *
425 * \param fluidState The fluid state of the two-phase model
426 * \param phaseIdx Index of the fluid phase
427 * \param compIdx index of the component
428 */
429 template <class FluidState>
430 static Scalar fugacityCoefficient(const FluidState &fluidState,
431 int phaseIdx,
432 int compIdx)
433 {
434 assert(0 <= phaseIdx && phaseIdx < numPhases);
435 assert(0 <= compIdx && compIdx < numComponents);
436
437 if (phaseIdx == compIdx)
438 // We could calculate the real fugacity coefficient of
439 // the component in the fluid. Probably that's not worth
440 // the effort, since the fugacity coefficient of the other
441 // component is infinite anyway...
442 return 1.0;
443 return std::numeric_limits<Scalar>::infinity();
444 }
445
446 using Base<Scalar, ThisType>::diffusionCoefficient;
447 /*!
448 * \brief Calculate the binary molecular diffusion coefficient for
449 * a component in a fluid phase \f$\mathrm{[mol^2 * s / (kg*m^3)]}\f$
450 * \param fluidState The fluid state of the two-phase model
451 * \param phaseIdx Index of the fluid phase
452 * \param compIdx index of the component
453 *
454 * Molecular diffusion of a component \f$\mathrm{\kappa}\f$ is caused by a
455 * gradient of the chemical potential and follows the law
456 *
457 * \f[ J = - D \nabla \mu_\kappa \f]
458 *
459 * where \f$\mathrm{\mu_\kappa]}\f$ is the component's chemical potential,
460 * \f$\mathrm{D}\f$ is the diffusion coefficient and \f$\mathrm{J}\f$ is the
461 * diffusive flux. \f$\mathrm{\mu_\kappa}\f$ is connected to the component's
462 * fugacity \f$\mathrm{f_\kappa}\f$ by the relation
463 *
464 * \f[ \mu_\kappa = R T_\alpha \mathrm{ln} \frac{f_\kappa}{p_\alpha} \f]
465 *
466 * where \f$\mathrm{p_\alpha}\f$ and \f$\mathrm{T_\alpha}\f$ are the fluid phase'
467 * pressure and temperature.
468 */
469 template <class FluidState>
470 9 static Scalar diffusionCoefficient(const FluidState &fluidState,
471 int phaseIdx,
472 int compIdx)
473 {
474
7/16
✓ Branch 2 taken 9 times.
✗ Branch 3 not taken.
✓ Branch 11 taken 9 times.
✗ Branch 12 not taken.
✓ Branch 15 taken 9 times.
✗ Branch 16 not taken.
✓ Branch 18 taken 9 times.
✗ Branch 19 not taken.
✓ Branch 21 taken 9 times.
✗ Branch 22 not taken.
✓ Branch 23 taken 9 times.
✗ Branch 24 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 9 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
99 DUNE_THROW(Dune::InvalidStateException,
475 "Diffusion coefficients of components are meaningless if"
476 " immiscibility is assumed");
477 }
478
479 using Base<Scalar, ThisType>::binaryDiffusionCoefficient;
480 /*!
481 * \brief Given a phase's composition, temperature and pressure,
482 * return the binary diffusion coefficient \f$\mathrm{[m^2/s]}\f$ for components
483 * \f$\mathrm{i}\f$ and \f$\mathrm{j}\f$ in this phase.
484 * \param fluidState The fluid state of the two-phase model
485 * \param phaseIdx Index of the fluid phase
486 * \param compIIdx index of the component i
487 * \param compJIdx index of the component j
488 */
489 template <class FluidState>
490 27 static Scalar binaryDiffusionCoefficient(const FluidState &fluidState,
491 int phaseIdx,
492 int compIIdx,
493 int compJIdx)
494
495 {
496
7/16
✓ Branch 2 taken 27 times.
✗ Branch 3 not taken.
✓ Branch 11 taken 27 times.
✗ Branch 12 not taken.
✓ Branch 15 taken 27 times.
✗ Branch 16 not taken.
✓ Branch 18 taken 27 times.
✗ Branch 19 not taken.
✓ Branch 21 taken 27 times.
✗ Branch 22 not taken.
✓ Branch 23 taken 27 times.
✗ Branch 24 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 27 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
297 DUNE_THROW(Dune::InvalidStateException,
497 "Binary diffusion coefficients of components are meaningless if"
498 " immiscibility is assumed");
499 }
500
501 using Base<Scalar, ThisType>::enthalpy;
502 /*!
503 * \brief Return the specific enthalpy of a fluid phase \f$\mathrm{[J/kg]}\f$.
504 * \param fluidState The fluid state of the two-phase model
505 * \param phaseIdx Index of the fluid phase
506 */
507 template <class FluidState>
508 3 static Scalar enthalpy(const FluidState &fluidState,
509 int phaseIdx)
510 {
511
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 3 times.
3 assert(0 <= phaseIdx && phaseIdx < numPhases);
512
513 3 Scalar temperature = fluidState.temperature(phaseIdx);
514 3 Scalar pressure = fluidState.pressure(phaseIdx);
515
516
3/4
✓ Branch 0 taken 1 times.
✓ Branch 1 taken 1 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
3 switch(phaseIdx)
517 {
518 2 case wPhaseIdx: return WettingFluid::enthalpy(temperature, pressure); break;
519 2 case nPhaseIdx: return NonwettingFluid::enthalpy(temperature, pressure); break;
520 1 case gPhaseIdx: return Gas::enthalpy(temperature, pressure); break;
521 default: DUNE_THROW(Dune::InvalidStateException, "Invalid phase index");
522 }
523 }
524
525 using Base<Scalar, ThisType>::thermalConductivity;
526 /*!
527 * \brief Thermal conductivity of a fluid phase \f$\mathrm{[W/(m K)]}\f$.
528 * \param fluidState The fluid state of the two-phase model
529 * \param phaseIdx Index of the fluid phase
530 */
531 template <class FluidState>
532 3 static Scalar thermalConductivity(const FluidState &fluidState,
533 int phaseIdx)
534 {
535
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 3 times.
3 assert(0 <= phaseIdx && phaseIdx < numPhases);
536
537 3 Scalar temperature = fluidState.temperature(phaseIdx);
538 3 Scalar pressure = fluidState.pressure(phaseIdx);
539
540
3/4
✓ Branch 0 taken 1 times.
✓ Branch 1 taken 1 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
3 switch(phaseIdx)
541 {
542 2 case wPhaseIdx: return WettingFluid::thermalConductivity(temperature, pressure); break;
543 2 case nPhaseIdx: return NonwettingFluid::thermalConductivity(temperature, pressure); break;
544 2 case gPhaseIdx: return Gas::thermalConductivity(temperature, pressure); break;
545 default: DUNE_THROW(Dune::InvalidStateException, "Invalid phase index");
546 }
547 }
548
549 using Base<Scalar, ThisType>::heatCapacity;
550 /*!
551 * @copybrief Base::thermalConductivity
552 *
553 * Additional comments:
554 *
555 * Specific isobaric heat capacity of a fluid phase.
556 * \f$\mathrm{[J/(kg*K)]}\f$.
557 *
558 * \param fluidState The fluid state of the two-phase model
559 * \param phaseIdx for which phase to give back the heat capacity
560 */
561 template <class FluidState>
562 3 static Scalar heatCapacity(const FluidState &fluidState,
563 int phaseIdx)
564 {
565
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 3 times.
3 assert(0 <= phaseIdx && phaseIdx < numPhases);
566
567 3 Scalar temperature = fluidState.temperature(phaseIdx);
568 3 Scalar pressure = fluidState.pressure(phaseIdx);
569
570
3/4
✓ Branch 0 taken 1 times.
✓ Branch 1 taken 1 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
3 switch(phaseIdx)
571 {
572 2 case wPhaseIdx: return WettingFluid::heatCapacity(temperature, pressure); break;
573 2 case nPhaseIdx: return NonwettingFluid::heatCapacity(temperature, pressure); break;
574 1 case gPhaseIdx: return Gas::heatCapacity(temperature, pressure); break;
575 default: DUNE_THROW(Dune::InvalidStateException, "Invalid phase index");
576 }
577 }
578 };
579
580 } // end namespace FluidSystems
581 } // end namespace Dumux
582
583 #endif
584