3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
h2oair.hh
Go to the documentation of this file.
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
19
25#ifndef DUMUX_H2O_AIR_SYSTEM_HH
26#define DUMUX_H2O_AIR_SYSTEM_HH
27
28#include <cassert>
29#include <iomanip>
30
33
38
40
41#include <dumux/io/name.hh>
42
43namespace Dumux {
44namespace FluidSystems {
49template<bool fastButSimplifiedRelations = false>
51{
52 static constexpr bool useH2ODensityAsLiquidMixtureDensity() { return fastButSimplifiedRelations; }
53 static constexpr bool useIdealGasDensity() { return fastButSimplifiedRelations; }
54 static constexpr bool useAirViscosityAsGasMixtureViscosity() { return fastButSimplifiedRelations; }
55};
56
68template <class Scalar,
70 class Policy = H2OAirDefaultPolicy<>,
71 bool useKelvinVaporPressure = false>
72class H2OAir
73: public Base<Scalar, H2OAir<Scalar, H2Otype, Policy> >
74{
78
79public:
80 using H2O = H2Otype;
82
83 static constexpr int numPhases = 2;
84 static constexpr int numComponents = 2;
85
86 static constexpr int liquidPhaseIdx = 0;
87 static constexpr int gasPhaseIdx = 1;
88 static constexpr int phase0Idx = liquidPhaseIdx;
89 static constexpr int phase1Idx = gasPhaseIdx;
90
91 static constexpr int H2OIdx = 0;
92 static constexpr int AirIdx = 1;
93 static constexpr int comp0Idx = H2OIdx;
94 static constexpr int comp1Idx = AirIdx;
95 static constexpr int liquidCompIdx = H2OIdx;
96 static constexpr int gasCompIdx = AirIdx;
97
103 static std::string phaseName(int phaseIdx)
104 {
105 assert(0 <= phaseIdx && phaseIdx < numPhases);
106 switch (phaseIdx)
107 {
108 case liquidPhaseIdx: return IOName::liquidPhase();
109 case gasPhaseIdx: return IOName::gaseousPhase();
110 }
111 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
112 }
113
117 static constexpr bool isMiscible()
118 { return true; }
119
125 static constexpr bool isGas(int phaseIdx)
126 {
127 assert(0 <= phaseIdx && phaseIdx < numPhases);
128 return phaseIdx == gasPhaseIdx;
129 }
130
145 static constexpr bool isIdealMixture(int phaseIdx)
146 {
147 assert(0 <= phaseIdx && phaseIdx < numPhases);
148 // we assume Henry's and Raoult's laws for the water phase and
149 // and no interaction between gas molecules of different
150 // components, so all phases are ideal mixtures!
151 return true;
152 }
153
163 static constexpr bool isCompressible(int phaseIdx)
164 {
165 assert(0 <= phaseIdx && phaseIdx < numPhases);
166 // ideal gases are always compressible
167 if (phaseIdx == gasPhaseIdx)
168 return true;
169 // the water component decides for the liquid phase...
170 return H2O::liquidIsCompressible();
171 }
172
179 static constexpr bool isIdealGas(int phaseIdx)
180 {
181 assert(0 <= phaseIdx && phaseIdx < numPhases);
182
183 // let the fluids decide
184 if (phaseIdx == gasPhaseIdx)
185 return H2O::gasIsIdeal() && Air::gasIsIdeal();
186 return false; // not a gas
187 }
188
189 /****************************************
190 * Component related static parameters
191 ****************************************/
197 static std::string componentName(int compIdx)
198 {
199 switch (compIdx)
200 {
201 case H2OIdx: return H2O::name();
202 case AirIdx: return Air::name();
203 }
204 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx);
205 }
206
212 static Scalar molarMass(int compIdx)
213 {
214 switch (compIdx)
215 {
216 case H2OIdx: return H2O::molarMass();
217 case AirIdx: return Air::molarMass();
218 }
219 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx);
220 }
221
227 static Scalar criticalTemperature(int compIdx)
228 {
229 static const Scalar TCrit[] = {
230 H2O::criticalTemperature(),
232 };
233
234 assert(0 <= compIdx && compIdx < numComponents);
235 return TCrit[compIdx];
236 }
237
243 static Scalar criticalPressure(int compIdx)
244 {
245 static const Scalar pCrit[] = {
246 H2O::criticalPressure(),
248 };
249
250 assert(0 <= compIdx && compIdx < numComponents);
251 return pCrit[compIdx];
252 }
253
260 template <class FluidState>
261 static Scalar vaporPressure(const FluidState &fluidState, int compIdx)
262 {
263 if (compIdx == H2OIdx)
264 {
265 const auto t = fluidState.temperature(H2OIdx);
266 if (!useKelvinVaporPressure)
267 return H2O::vaporPressure(t);
268 else
269 {
270 const auto pc = (fluidState.wettingPhase() == (int) H2OIdx)
271 ? fluidState.pressure(AirIdx)-fluidState.pressure(H2OIdx)
272 : fluidState.pressure(H2OIdx)-fluidState.pressure(AirIdx);
273 return H2O::vaporPressure(t)*exp( -pc * molarMass(H2OIdx)
274 / density(fluidState, H2OIdx)
276 }
277 }
278 else if (compIdx == AirIdx)
279 // return Air::vaporPressure(fluidState.temperature(AirIdx));
280 DUNE_THROW(Dune::NotImplemented, "Air::vaporPressure(t)");
281 else
282 DUNE_THROW(Dune::NotImplemented, "Invalid component index " << compIdx);
283 }
284
290 static Scalar criticalMolarVolume(int compIdx)
291 {
292 DUNE_THROW(Dune::NotImplemented,
293 "H2OAirFluidSystem::criticalMolarVolume()");
294 }
295
301 static Scalar acentricFactor(int compIdx)
302 {
303 static const Scalar accFac[] = {
304 H2O::acentricFactor(),
305 Air::acentricFactor()
306 };
307
308 assert(0 <= compIdx && compIdx < numComponents);
309 return accFac[compIdx];
310 }
311
312 /****************************************
313 * thermodynamic relations
314 ****************************************/
315
322 static void init()
323 {
324 init(/*tempMin=*/273.15,
325 /*tempMax=*/623.15,
326 /*numTemp=*/100,
327 /*pMin=*/-10.,
328 /*pMax=*/20e6,
329 /*numP=*/200);
330 }
331
343 static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp,
344 Scalar pressMin, Scalar pressMax, unsigned nPress)
345 {
346 std::cout << "The H2O-air fluid system was configured with the following policy:\n";
347 std::cout << " - use H2O density as liquid mixture density: " << std::boolalpha << Policy::useH2ODensityAsLiquidMixtureDensity() << "\n";
348 std::cout << " - use ideal gas density: " << std::boolalpha << Policy::useIdealGasDensity() << "\n";
349 std::cout << " - use air viscosity as gas mixture viscosity: " << std::boolalpha << Policy::useAirViscosityAsGasMixtureViscosity() << std::endl;
350
351 if (H2O::isTabulated)
352 {
353 H2O::init(tempMin, tempMax, nTemp,
354 pressMin, pressMax, nPress);
355 }
356 }
357
358 using Base::density;
372 template <class FluidState>
373 static Scalar density(const FluidState &fluidState,
374 const int phaseIdx)
375 {
376 assert(0 <= phaseIdx && phaseIdx < numPhases);
377
378 const Scalar T = fluidState.temperature(phaseIdx);
379 const Scalar p = fluidState.pressure(phaseIdx);
380
381 if (phaseIdx == phase0Idx)
382 {
383 if (Policy::useH2ODensityAsLiquidMixtureDensity())
384 // assume pure water
385 return H2O::liquidDensity(T, p);
386 else
387 {
388 // See: Eq. (7) in Class et al. (2002a)
389 // This assumes each gas molecule displaces exactly one
390 // molecule in the liquid.
391 return H2O::liquidMolarDensity(T, p)
392 * (H2O::molarMass()*fluidState.moleFraction(liquidPhaseIdx, H2OIdx)
393 + Air::molarMass()*fluidState.moleFraction(liquidPhaseIdx, AirIdx));
394 }
395 }
396 else if (phaseIdx == gasPhaseIdx)
397 {
398 if (Policy::useIdealGasDensity())
399 // for the gas phase assume an ideal gas
400 {
401 const Scalar averageMolarMass = fluidState.averageMolarMass(gasPhaseIdx);
402 return IdealGas::density(averageMolarMass, T, p);
403 }
404
405 return H2O::gasDensity(T, fluidState.partialPressure(gasPhaseIdx, H2OIdx))
406 + Air::gasDensity(T, fluidState.partialPressure(gasPhaseIdx, AirIdx));
407 }
408 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
409 }
410
411 using Base::molarDensity;
421 template <class FluidState>
422 static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
423 {
424 assert(0 <= phaseIdx && phaseIdx < numPhases);
425
426 const Scalar T = fluidState.temperature(phaseIdx);
427 const Scalar p = fluidState.pressure(phaseIdx);
428
429 if (phaseIdx == phase0Idx)
430 {
431 // assume pure water or that each gas molecule displaces exactly one
432 // molecule in the liquid.
433 return H2O::liquidMolarDensity(T, p);
434 }
435 else if (phaseIdx == phase1Idx)
436 {
437 if (Policy::useIdealGasDensity())
438 // for the gas phase assume an ideal gas
439 { return IdealGas::molarDensity(T, p); }
440
441 return H2O::gasMolarDensity(T, fluidState.partialPressure(phase1Idx, H2OIdx))
442 + Air::gasMolarDensity(T, fluidState.partialPressure(phase1Idx, AirIdx));
443 }
444 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
445 }
446
447 using Base::viscosity;
460 template <class FluidState>
461 static Scalar viscosity(const FluidState &fluidState,
462 int phaseIdx)
463 {
464 assert(0 <= phaseIdx && phaseIdx < numPhases);
465
466 Scalar T = fluidState.temperature(phaseIdx);
467 Scalar p = fluidState.pressure(phaseIdx);
468
469 if (phaseIdx == liquidPhaseIdx)
470 {
471 // assume pure water for the liquid phase
472 return H2O::liquidViscosity(T, p);
473 }
474 else if (phaseIdx == gasPhaseIdx)
475 {
476 if(Policy::useAirViscosityAsGasMixtureViscosity()){
477 return Air::gasViscosity(T, p);
478 }
479 else //using a complicated version of this fluid system
480 {
481 // Wilke method (Reid et al.):
482 Scalar muResult = 0;
483 const Scalar mu[numComponents] = {
486 };
487
488 // molar masses
489 const Scalar M[numComponents] = {
490 H2O::molarMass(),
492 };
493
494 for (int i = 0; i < numComponents; ++i)
495 {
496 Scalar divisor = 0;
497 using std::sqrt;
498 for (int j = 0; j < numComponents; ++j)
499 {
500 // 1 + (mu[i]/mu[j]^1/2 * (M[i]/M[j])^1/4)
501 Scalar phiIJ = 1 + sqrt(mu[i]/mu[j] * sqrt(M[j]/M[i]));
502 phiIJ = phiIJ * phiIJ / sqrt(8*(1 + M[i]/M[j]));
503 divisor += fluidState.moleFraction(phaseIdx, j)*phiIJ;
504 }
505 muResult += fluidState.moleFraction(phaseIdx, i)*mu[i] / divisor;
506 }
507 return muResult;
508 }
509 }
510 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
511 }
512
532 template <class FluidState>
533 static Scalar fugacityCoefficient(const FluidState &fluidState,
534 int phaseIdx,
535 int compIdx)
536 {
537 assert(0 <= phaseIdx && phaseIdx < numPhases);
538 assert(0 <= compIdx && compIdx < numComponents);
539
540 Scalar T = fluidState.temperature(phaseIdx);
541 Scalar p = fluidState.pressure(phaseIdx);
542
543 if (phaseIdx == liquidPhaseIdx) {
544 if (compIdx == H2OIdx)
545 return vaporPressure(fluidState, compIdx)/p;
546 return BinaryCoeff::H2O_Air::henry(T)/p;
547 }
548
549 // for the gas phase, assume an ideal gas when it comes to
550 // fugacity (-> fugacity == partial pressure)
551 return 1.0;
552 }
553
560 template <class FluidState>
561 static Scalar relativeHumidity(const FluidState &fluidState)
562 {
563 return fluidState.partialPressure(gasPhaseIdx, comp0Idx)
564 / H2O::vaporPressure(fluidState.temperature(gasPhaseIdx));
565 }
566
568 template <class FluidState>
569 static Scalar diffusionCoefficient(const FluidState &fluidState,
570 int phaseIdx,
571 int compIdx)
572 {
573 DUNE_THROW(Dune::NotImplemented, "FluidSystems::H2OAir::diffusionCoefficient()");
574 }
575
587 template <class FluidState>
588 static Scalar binaryDiffusionCoefficient(const FluidState &fluidState,
589 int phaseIdx,
590 int compIIdx,
591 int compJIdx)
592 {
593 using std::swap;
594 if (compIIdx > compJIdx)
595 swap(compIIdx, compJIdx);
596
597 Scalar T = fluidState.temperature(phaseIdx);
598 Scalar p = fluidState.pressure(phaseIdx);
599
600 // we are in the liquid phase
601 if (phaseIdx == liquidPhaseIdx)
602 {
603 if (compIIdx == H2OIdx && compJIdx == AirIdx)
605 else
606 DUNE_THROW(Dune::InvalidStateException,
607 "Binary diffusion coefficient of components "
608 << compIIdx << " and " << compJIdx
609 << " in phase " << phaseIdx << " is undefined!\n");
610 }
611
612 // we are in the gas phase
613 else if (phaseIdx == gasPhaseIdx)
614 {
615 if (compIIdx == H2OIdx && compJIdx == AirIdx)
617 else
618 DUNE_THROW(Dune::InvalidStateException,
619 "Binary diffusion coefficient of components "
620 << compIIdx << " and " << compJIdx
621 << " in phase " << phaseIdx << " is undefined!\n");
622 }
623
624 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
625 }
626
627 using Base::enthalpy;
645 template <class FluidState>
646 static Scalar enthalpy(const FluidState &fluidState,
647 int phaseIdx)
648 {
649 const Scalar T = fluidState.temperature(phaseIdx);
650 const Scalar p = fluidState.pressure(phaseIdx);
651
652 if (phaseIdx == liquidPhaseIdx)
653 return H2O::liquidEnthalpy(T, p);
654
655 else if (phaseIdx == gasPhaseIdx)
656 return H2O::gasEnthalpy(T, p)*fluidState.massFraction(gasPhaseIdx, H2OIdx)
657 + Air::gasEnthalpy(T, p)*fluidState.massFraction(gasPhaseIdx, AirIdx);
658
659 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
660 }
661
669 template <class FluidState>
670 static Scalar componentEnthalpy(const FluidState &fluidState,
671 int phaseIdx,
672 int componentIdx)
673 {
674 const Scalar T = fluidState.temperature(phaseIdx);
675 const Scalar p = fluidState.pressure(phaseIdx);
676
677 if (phaseIdx == liquidPhaseIdx)
678 {
679 // the liquid enthalpy is constant
680 return H2O::liquidEnthalpy(T, p);
681 }
682 else if (phaseIdx == gasPhaseIdx)
683 {
684 if (componentIdx == H2OIdx)
685 {
686 return H2O::gasEnthalpy(T, p);
687 }
688 else if (componentIdx == AirIdx)
689 {
690 return Air::gasEnthalpy(T, p);
691 }
692 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx);
693 }
694 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
695 }
696
707 template <class FluidState>
708 static Scalar thermalConductivity(const FluidState &fluidState,
709 int phaseIdx)
710 {
711 assert(0 <= phaseIdx && phaseIdx < numPhases);
712
713 const Scalar temperature = fluidState.temperature(phaseIdx) ;
714 const Scalar pressure = fluidState.pressure(phaseIdx);
715 if (phaseIdx == liquidPhaseIdx)
716 {
717 return H2O::liquidThermalConductivity(temperature, pressure);
718 }
719 else if (phaseIdx == gasPhaseIdx)
720 {
722 }
723 else
724 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
725 }
726
727 using Base::heatCapacity;
738 template <class FluidState>
739 static Scalar heatCapacity(const FluidState &fluidState,
740 int phaseIdx)
741 {
742 const Scalar temperature = fluidState.temperature(phaseIdx);
743 const Scalar pressure = fluidState.pressure(phaseIdx);
744 if (phaseIdx == liquidPhaseIdx)
745 {
746 // influence of air is neglected
747 return H2O::liquidHeatCapacity(temperature, pressure);
748 }
749 else if (phaseIdx == gasPhaseIdx)
750 {
751 return Air::gasHeatCapacity(temperature, pressure) * fluidState.moleFraction(gasPhaseIdx, AirIdx)
752 + H2O::gasHeatCapacity(temperature, pressure) * fluidState.moleFraction(gasPhaseIdx, H2OIdx);
753 }
754 else
755 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
756 }
757};
758
759} // end namespace FluidSystems
760} // end namespace Dumux
761
762#endif
Some exceptions thrown in DuMux
A collection of input/output field names for common physical quantities.
Binary coefficients for water and air.
Tabulates all thermodynamic properties of a given untabulated chemical species.
A simple class for the air fluid properties.
Material properties of pure water .
Relations valid for an ideal gas.
Definition: adapt.hh:29
std::string temperature() noexcept
I/O name of temperature for equilibrium models.
Definition: name.hh:51
std::string gaseousPhase() noexcept
I/O name of gaseous phase.
Definition: name.hh:123
std::string liquidPhase() noexcept
I/O name of liquid phase.
Definition: name.hh:119
std::string pressure(int phaseIdx) noexcept
I/O name of pressure for multiphase systems.
Definition: name.hh:34
Scalar h2oGasViscosityInMixture(Scalar temperature, Scalar pressure)
The dynamic viscosity of steam in a gas mixture.
Definition: h2o.hh:977
static Scalar henry(Scalar temperature)
Henry coefficient for air in liquid water.
Definition: h2o_air.hh:48
static Scalar gasDiffCoeff(Scalar temperature, Scalar pressure)
Binary diffusion coefficient for molecular water and air.
Definition: h2o_air.hh:68
static Scalar liquidDiffCoeff(Scalar temperature, Scalar pressure)
Diffusion coefficient for molecular nitrogen in liquid water.
Definition: h2o_air.hh:101
A class for the air fluid properties.
Definition: air.hh:47
static Scalar gasDensity(Scalar temperature, Scalar pressure)
The density of Air at a given pressure and temperature.
Definition: air.hh:85
static constexpr Scalar molarMass()
The molar mass in of Air.
Definition: air.hh:62
static const Scalar gasHeatCapacity(Scalar temperature, Scalar pressure)
Specific isobaric heat capacity of pure air.
Definition: air.hh:307
static Scalar criticalPressure()
Returns the critical pressure of Air.
Definition: air.hh:74
static Scalar gasViscosity(Scalar temperature, Scalar pressure)
The dynamic viscosity of Air at a given pressure and temperature.
Definition: air.hh:187
static Scalar gasThermalConductivity(Scalar temperature, Scalar pressure)
Thermal conductivity of air.
Definition: air.hh:345
static constexpr bool gasIsIdeal()
Returns true, the gas phase is assumed to be ideal.
Definition: air.hh:109
static Scalar gasEnthalpy(Scalar temperature, Scalar pressure)
Specific enthalpy of Air with 273.15 as basis.
Definition: air.hh:270
static Scalar criticalTemperature()
Returns the critical temperature of Air.
Definition: air.hh:68
static std::string name()
A human readable name for Air.
Definition: air.hh:54
static Scalar gasMolarDensity(Scalar temperature, Scalar pressure)
The molar density of air in , depending on pressure and temperature.
Definition: air.hh:97
Tabulates all thermodynamic properties of a given untabulated chemical species.
Definition: tabulatedcomponent.hh:82
A central place for various physical constants occuring in some equations.
Definition: constants.hh:39
Fluid system base class.
Definition: fluidsystems/base.hh:45
static Scalar density(const FluidState &fluidState, int phaseIdx)
Calculate the density of a fluid phase.
Definition: fluidsystems/base.hh:134
static Scalar thermalConductivity(const FluidState &fluidState, int phaseIdx)
Thermal conductivity of a fluid phase .
Definition: fluidsystems/base.hh:390
static Scalar fugacityCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Calculate the fugacity coefficient of an individual component in a fluid phase.
Definition: fluidsystems/base.hh:197
static Scalar diffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Calculate the binary molecular diffusion coefficient for a component in a fluid phase .
Definition: fluidsystems/base.hh:278
static Scalar binaryDiffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIIdx, int compJIdx)
Given a phase's composition, temperature and pressure, return the binary diffusion coefficient for c...
Definition: fluidsystems/base.hh:326
static Scalar enthalpy(const FluidState &fluidState, int phaseIdx)
Given a phase's composition, temperature, pressure and density, calculate its specific enthalpy .
Definition: fluidsystems/base.hh:363
static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
Calculate the molar density of a fluid phase.
Definition: fluidsystems/base.hh:160
static Scalar viscosity(const FluidState &fluidState, int phaseIdx)
Calculate the dynamic viscosity of a fluid phase .
Definition: fluidsystems/base.hh:236
static Scalar heatCapacity(const FluidState &fluidState, int phaseIdx)
Specific isobaric heat capacity of a fluid phase .
Definition: fluidsystems/base.hh:424
Policy for the H2O-air fluid system.
Definition: h2oair.hh:51
static constexpr bool useAirViscosityAsGasMixtureViscosity()
Definition: h2oair.hh:54
static constexpr bool useIdealGasDensity()
Definition: h2oair.hh:53
static constexpr bool useH2ODensityAsLiquidMixtureDensity()
Definition: h2oair.hh:52
A compositional two-phase fluid system with water and air as components in both, the liquid and the g...
Definition: h2oair.hh:74
static Scalar viscosity(const FluidState &fluidState, int phaseIdx)
Calculate the dynamic viscosity of a fluid phase .
Definition: h2oair.hh:461
static constexpr bool isGas(int phaseIdx)
Return whether a phase is gaseous.
Definition: h2oair.hh:125
static Scalar enthalpy(const FluidState &fluidState, int phaseIdx)
Given a phase's composition, temperature and pressure, return its specific enthalpy .
Definition: h2oair.hh:646
static Scalar criticalMolarVolume(int compIdx)
Molar volume of a component at the critical point .
Definition: h2oair.hh:290
static Scalar vaporPressure(const FluidState &fluidState, int compIdx)
Vapor pressure of a component .
Definition: h2oair.hh:261
static constexpr int comp1Idx
index of the second component
Definition: h2oair.hh:94
static constexpr bool isIdealGas(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: h2oair.hh:179
static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp, Scalar pressMin, Scalar pressMax, unsigned nPress)
Initialize the fluid system's static parameters using problem specific temperature and pressure range...
Definition: h2oair.hh:343
static std::string componentName(int compIdx)
Return the human readable name of a component.
Definition: h2oair.hh:197
static constexpr int numComponents
Number of components in the fluid system.
Definition: h2oair.hh:84
static constexpr int liquidPhaseIdx
index of the liquid phase
Definition: h2oair.hh:86
static constexpr int gasPhaseIdx
index of the gas phase
Definition: h2oair.hh:87
static Scalar heatCapacity(const FluidState &fluidState, int phaseIdx)
Specific isobaric heat capacity of a fluid phase. .
Definition: h2oair.hh:739
static std::string phaseName(int phaseIdx)
Return the human readable name of a phase.
Definition: h2oair.hh:103
static Scalar thermalConductivity(const FluidState &fluidState, int phaseIdx)
Thermal conductivity of a fluid phase .
Definition: h2oair.hh:708
static Scalar componentEnthalpy(const FluidState &fluidState, int phaseIdx, int componentIdx)
Returns the specific enthalpy of a component in a specific phase.
Definition: h2oair.hh:670
static constexpr bool isIdealMixture(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: h2oair.hh:145
static Scalar criticalPressure(int compIdx)
Critical pressure of a component .
Definition: h2oair.hh:243
static constexpr int phase1Idx
index of the second phase
Definition: h2oair.hh:89
static constexpr int comp0Idx
index of the frist component
Definition: h2oair.hh:93
static constexpr int numPhases
Number of phases in the fluid system.
Definition: h2oair.hh:83
static constexpr bool isCompressible(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: h2oair.hh:163
static constexpr int H2OIdx
index of the frist component
Definition: h2oair.hh:91
H2Otype H2O
Definition: h2oair.hh:80
static Scalar fugacityCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Returns the fugacity coefficient of a component in a phase.
Definition: h2oair.hh:533
static Scalar relativeHumidity(const FluidState &fluidState)
Returns the relative humidity of the gas phase.
Definition: h2oair.hh:561
static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
The molar density of a fluid phase in .
Definition: h2oair.hh:422
static constexpr int liquidCompIdx
index of the liquid component
Definition: h2oair.hh:95
static void init()
Initialize the fluid system's static parameters generically.
Definition: h2oair.hh:322
static Scalar density(const FluidState &fluidState, const int phaseIdx)
Given a phase's composition, temperature, pressure, and the partial pressures of all components,...
Definition: h2oair.hh:373
static constexpr int AirIdx
index of the second component
Definition: h2oair.hh:92
static Scalar acentricFactor(int compIdx)
The acentric factor of a component .
Definition: h2oair.hh:301
static Scalar binaryDiffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIIdx, int compJIdx)
Given a phase's composition, temperature and pressure, return the binary diffusion coefficient for c...
Definition: h2oair.hh:588
static constexpr int phase0Idx
index of the first phase
Definition: h2oair.hh:88
static constexpr int gasCompIdx
index of the gas component
Definition: h2oair.hh:96
static constexpr bool isMiscible()
Returns whether the fluids are miscible.
Definition: h2oair.hh:117
static Scalar criticalTemperature(int compIdx)
Critical temperature of a component .
Definition: h2oair.hh:227
static Scalar diffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Definition: h2oair.hh:569
static Scalar molarMass(int compIdx)
Return the molar mass of a component .
Definition: h2oair.hh:212
Relations valid for an ideal gas.
Definition: idealgas.hh:37
static constexpr Scalar density(Scalar avgMolarMass, Scalar temperature, Scalar pressure)
The density of the gas in , depending on pressure, temperature and average molar mass of the gas.
Definition: idealgas.hh:49
static constexpr Scalar molarDensity(Scalar temperature, Scalar pressure)
The molar density of the gas , depending on pressure and temperature.
Definition: idealgas.hh:70
Fluid system base class.