3.2-git
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
41
42#include <dumux/io/name.hh>
43
44namespace Dumux {
45namespace FluidSystems {
50template<bool fastButSimplifiedRelations = false>
52{
53 static constexpr bool useH2ODensityAsLiquidMixtureDensity() { return fastButSimplifiedRelations; }
54 static constexpr bool useIdealGasDensity() { return fastButSimplifiedRelations; }
55 static constexpr bool useAirViscosityAsGasMixtureViscosity() { return fastButSimplifiedRelations; }
56};
57
69template <class Scalar,
71 class Policy = H2OAirDefaultPolicy<>,
72 bool useKelvinVaporPressure = false>
73class H2OAir
74: public Base<Scalar, H2OAir<Scalar, H2Otype, Policy> >
75{
79
80public:
81 using H2O = H2Otype;
83
84 static constexpr int numPhases = 2;
85 static constexpr int numComponents = 2;
86
87 static constexpr int liquidPhaseIdx = 0;
88 static constexpr int gasPhaseIdx = 1;
89 static constexpr int phase0Idx = liquidPhaseIdx;
90 static constexpr int phase1Idx = gasPhaseIdx;
91
92 static constexpr int H2OIdx = 0;
93 static constexpr int AirIdx = 1;
94 static constexpr int comp0Idx = H2OIdx;
95 static constexpr int comp1Idx = AirIdx;
96 static constexpr int liquidCompIdx = H2OIdx;
97 static constexpr int gasCompIdx = AirIdx;
98
104 static std::string phaseName(int phaseIdx)
105 {
106 assert(0 <= phaseIdx && phaseIdx < numPhases);
107 switch (phaseIdx)
108 {
109 case liquidPhaseIdx: return IOName::liquidPhase();
110 case gasPhaseIdx: return IOName::gaseousPhase();
111 }
112 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
113 }
114
118 static constexpr bool isMiscible()
119 { return true; }
120
126 static constexpr bool isGas(int phaseIdx)
127 {
128 assert(0 <= phaseIdx && phaseIdx < numPhases);
129 return phaseIdx == gasPhaseIdx;
130 }
131
146 static constexpr bool isIdealMixture(int phaseIdx)
147 {
148 assert(0 <= phaseIdx && phaseIdx < numPhases);
149 // we assume Henry's and Raoult's laws for the water phase and
150 // and no interaction between gas molecules of different
151 // components, so all phases are ideal mixtures!
152 return true;
153 }
154
164 static constexpr bool isCompressible(int phaseIdx)
165 {
166 assert(0 <= phaseIdx && phaseIdx < numPhases);
167 // ideal gases are always compressible
168 if (phaseIdx == gasPhaseIdx)
169 return true;
170 // the water component decides for the liquid phase...
171 return H2O::liquidIsCompressible();
172 }
173
180 static constexpr bool isIdealGas(int phaseIdx)
181 {
182 assert(0 <= phaseIdx && phaseIdx < numPhases);
183
184 // let the fluids decide
185 if (phaseIdx == gasPhaseIdx)
186 return H2O::gasIsIdeal() && Air::gasIsIdeal();
187 return false; // not a gas
188 }
189
190 /****************************************
191 * Component related static parameters
192 ****************************************/
198 static std::string componentName(int compIdx)
199 {
200 switch (compIdx)
201 {
202 case H2OIdx: return H2O::name();
203 case AirIdx: return Air::name();
204 }
205 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx);
206 }
207
213 static Scalar molarMass(int compIdx)
214 {
215 switch (compIdx)
216 {
217 case H2OIdx: return H2O::molarMass();
218 case AirIdx: return Air::molarMass();
219 }
220 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx);
221 }
222
228 static Scalar criticalTemperature(int compIdx)
229 {
230 static const Scalar TCrit[] = {
231 H2O::criticalTemperature(),
233 };
234
235 assert(0 <= compIdx && compIdx < numComponents);
236 return TCrit[compIdx];
237 }
238
244 static Scalar criticalPressure(int compIdx)
245 {
246 static const Scalar pCrit[] = {
247 H2O::criticalPressure(),
249 };
250
251 assert(0 <= compIdx && compIdx < numComponents);
252 return pCrit[compIdx];
253 }
254
261 template <class FluidState>
262 static Scalar vaporPressure(const FluidState &fluidState, int compIdx)
263 {
264 if (compIdx == H2OIdx)
265 {
266 const auto t = fluidState.temperature(H2OIdx);
267 if (!useKelvinVaporPressure)
268 return H2O::vaporPressure(t);
269 else
270 {
271 const auto pc = (fluidState.wettingPhase() == (int) H2OIdx)
272 ? fluidState.pressure(AirIdx)-fluidState.pressure(H2OIdx)
273 : fluidState.pressure(H2OIdx)-fluidState.pressure(AirIdx);
274 return H2O::vaporPressure(t)*exp( -pc * molarMass(H2OIdx)
275 / density(fluidState, H2OIdx)
277 }
278 }
279 else if (compIdx == AirIdx)
280 // return Air::vaporPressure(fluidState.temperature(AirIdx));
281 DUNE_THROW(Dune::NotImplemented, "Air::vaporPressure(t)");
282 else
283 DUNE_THROW(Dune::NotImplemented, "Invalid component index " << compIdx);
284 }
285
291 static Scalar criticalMolarVolume(int compIdx)
292 {
293 DUNE_THROW(Dune::NotImplemented,
294 "H2OAirFluidSystem::criticalMolarVolume()");
295 }
296
302 static Scalar acentricFactor(int compIdx)
303 {
304 static const Scalar accFac[] = {
305 H2O::acentricFactor(),
306 Air::acentricFactor()
307 };
308
309 assert(0 <= compIdx && compIdx < numComponents);
310 return accFac[compIdx];
311 }
312
313 /****************************************
314 * thermodynamic relations
315 ****************************************/
316
323 static void init()
324 {
325 init(/*tempMin=*/273.15,
326 /*tempMax=*/623.15,
327 /*numTemp=*/100,
328 /*pMin=*/-10.,
329 /*pMax=*/20e6,
330 /*numP=*/200);
331 }
332
344 static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp,
345 Scalar pressMin, Scalar pressMax, unsigned nPress)
346 {
347 std::cout << "The H2O-air fluid system was configured with the following policy:\n";
348 std::cout << " - use H2O density as liquid mixture density: " << std::boolalpha << Policy::useH2ODensityAsLiquidMixtureDensity() << "\n";
349 std::cout << " - use ideal gas density: " << std::boolalpha << Policy::useIdealGasDensity() << "\n";
350 std::cout << " - use air viscosity as gas mixture viscosity: " << std::boolalpha << Policy::useAirViscosityAsGasMixtureViscosity() << std::endl;
351
352 if (H2O::isTabulated)
353 {
354 H2O::init(tempMin, tempMax, nTemp,
355 pressMin, pressMax, nPress);
356 }
357 }
358
359 using Base::density;
373 template <class FluidState>
374 static Scalar density(const FluidState &fluidState,
375 const int phaseIdx)
376 {
377 assert(0 <= phaseIdx && phaseIdx < numPhases);
378
379 const Scalar T = fluidState.temperature(phaseIdx);
380 const Scalar p = fluidState.pressure(phaseIdx);
381
382 if (phaseIdx == phase0Idx)
383 {
384 if (Policy::useH2ODensityAsLiquidMixtureDensity())
385 // assume pure water
386 return H2O::liquidDensity(T, p);
387 else
388 {
389 // See: Eq. (7) in Class et al. (2002a)
390 // This assumes each gas molecule displaces exactly one
391 // molecule in the liquid.
392 return H2O::liquidMolarDensity(T, p)
393 * (H2O::molarMass()*fluidState.moleFraction(liquidPhaseIdx, H2OIdx)
394 + Air::molarMass()*fluidState.moleFraction(liquidPhaseIdx, AirIdx));
395 }
396 }
397 else if (phaseIdx == gasPhaseIdx)
398 {
399 if (Policy::useIdealGasDensity())
400 // for the gas phase assume an ideal gas
401 {
402 const Scalar averageMolarMass = fluidState.averageMolarMass(gasPhaseIdx);
403 return IdealGas::density(averageMolarMass, T, p);
404 }
405
406 return H2O::gasDensity(T, fluidState.partialPressure(gasPhaseIdx, H2OIdx))
407 + Air::gasDensity(T, fluidState.partialPressure(gasPhaseIdx, AirIdx));
408 }
409 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
410 }
411
412 using Base::molarDensity;
422 template <class FluidState>
423 static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
424 {
425 assert(0 <= phaseIdx && phaseIdx < numPhases);
426
427 const Scalar T = fluidState.temperature(phaseIdx);
428 const Scalar p = fluidState.pressure(phaseIdx);
429
430 if (phaseIdx == phase0Idx)
431 {
432 // assume pure water or that each gas molecule displaces exactly one
433 // molecule in the liquid.
434 return H2O::liquidMolarDensity(T, p);
435 }
436 else if (phaseIdx == phase1Idx)
437 {
438 if (Policy::useIdealGasDensity())
439 // for the gas phase assume an ideal gas
440 { return IdealGas::molarDensity(T, p); }
441
442 return H2O::gasMolarDensity(T, fluidState.partialPressure(phase1Idx, H2OIdx))
443 + Air::gasMolarDensity(T, fluidState.partialPressure(phase1Idx, AirIdx));
444 }
445 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
446 }
447
448 using Base::viscosity;
461 template <class FluidState>
462 static Scalar viscosity(const FluidState &fluidState,
463 int phaseIdx)
464 {
465 assert(0 <= phaseIdx && phaseIdx < numPhases);
466
467 Scalar T = fluidState.temperature(phaseIdx);
468 Scalar p = fluidState.pressure(phaseIdx);
469
470 if (phaseIdx == liquidPhaseIdx)
471 {
472 // assume pure water for the liquid phase
473 return H2O::liquidViscosity(T, p);
474 }
475 else if (phaseIdx == gasPhaseIdx)
476 {
477 if(Policy::useAirViscosityAsGasMixtureViscosity()){
478 return Air::gasViscosity(T, p);
479 }
480 else //using a complicated version of this fluid system
481 {
482 // Wilke method (Reid et al.):
483 Scalar muResult = 0;
484 const Scalar mu[numComponents] = {
485 H2O::gasViscosity(T,
486 H2O::vaporPressure(T)),
488 };
489
490 // molar masses
491 const Scalar M[numComponents] = {
492 H2O::molarMass(),
494 };
495
496 for (int i = 0; i < numComponents; ++i)
497 {
498 Scalar divisor = 0;
499 using std::sqrt;
500 for (int j = 0; j < numComponents; ++j)
501 {
502 // 1 + (mu[i]/mu[j]^1/2 * (M[i]/M[j])^1/4)
503 Scalar phiIJ = 1 + sqrt(mu[i]/mu[j] * sqrt(M[j]/M[i]));
504 phiIJ = phiIJ * phiIJ / sqrt(8*(1 + M[i]/M[j]));
505 divisor += fluidState.moleFraction(phaseIdx, j)*phiIJ;
506 }
507 muResult += fluidState.moleFraction(phaseIdx, i)*mu[i] / divisor;
508 }
509 return muResult;
510 }
511 }
512 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
513 }
514
534 template <class FluidState>
535 static Scalar fugacityCoefficient(const FluidState &fluidState,
536 int phaseIdx,
537 int compIdx)
538 {
539 assert(0 <= phaseIdx && phaseIdx < numPhases);
540 assert(0 <= compIdx && compIdx < numComponents);
541
542 Scalar T = fluidState.temperature(phaseIdx);
543 Scalar p = fluidState.pressure(phaseIdx);
544
545 if (phaseIdx == liquidPhaseIdx) {
546 if (compIdx == H2OIdx)
547 return vaporPressure(fluidState, compIdx)/p;
548 return BinaryCoeff::H2O_Air::henry(T)/p;
549 }
550
551 // for the gas phase, assume an ideal gas when it comes to
552 // fugacity (-> fugacity == partial pressure)
553 return 1.0;
554 }
555
562 template <class FluidState>
563 static Scalar relativeHumidity(const FluidState &fluidState)
564 {
565 return fluidState.partialPressure(gasPhaseIdx, comp0Idx)
566 / H2O::vaporPressure(fluidState.temperature(gasPhaseIdx));
567 }
568
570 template <class FluidState>
571 static Scalar diffusionCoefficient(const FluidState &fluidState,
572 int phaseIdx,
573 int compIdx)
574 {
575 DUNE_THROW(Dune::NotImplemented, "FluidSystems::H2OAir::diffusionCoefficient()");
576 }
577
589 template <class FluidState>
590 static Scalar binaryDiffusionCoefficient(const FluidState &fluidState,
591 int phaseIdx,
592 int compIIdx,
593 int compJIdx)
594 {
595 using std::swap;
596 if (compIIdx > compJIdx)
597 swap(compIIdx, compJIdx);
598
599 Scalar T = fluidState.temperature(phaseIdx);
600 Scalar p = fluidState.pressure(phaseIdx);
601
602 // we are in the liquid phase
603 if (phaseIdx == liquidPhaseIdx)
604 {
605 if (compIIdx == H2OIdx && compJIdx == AirIdx)
607 else
608 DUNE_THROW(Dune::InvalidStateException,
609 "Binary diffusion coefficient of components "
610 << compIIdx << " and " << compJIdx
611 << " in phase " << phaseIdx << " is undefined!\n");
612 }
613
614 // we are in the gas phase
615 else if (phaseIdx == gasPhaseIdx)
616 {
617 if (compIIdx == H2OIdx && compJIdx == AirIdx)
619 else
620 DUNE_THROW(Dune::InvalidStateException,
621 "Binary diffusion coefficient of components "
622 << compIIdx << " and " << compJIdx
623 << " in phase " << phaseIdx << " is undefined!\n");
624 }
625
626 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
627 }
628
629 using Base::enthalpy;
647 template <class FluidState>
648 static Scalar enthalpy(const FluidState &fluidState,
649 int phaseIdx)
650 {
651 Scalar T = fluidState.temperature(phaseIdx);
652 Scalar p = fluidState.pressure(phaseIdx);
655
656 if (phaseIdx == liquidPhaseIdx)
657 {
658 return H2O::liquidEnthalpy(T, p);
659 }
660
661 else if (phaseIdx == gasPhaseIdx)
662 {
663 Scalar result = 0.0;
664 result +=
665 H2O::gasEnthalpy(T, p) *
666 fluidState.massFraction(gasPhaseIdx, H2OIdx);
667
668 result +=
669 Air::gasEnthalpy(T, p) *
670 fluidState.massFraction(gasPhaseIdx, AirIdx);
671 return result;
672 }
673 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
674 }
675
683 template <class FluidState>
684 static Scalar componentEnthalpy(const FluidState &fluidState,
685 int phaseIdx,
686 int componentIdx)
687 {
688 const Scalar T = fluidState.temperature(phaseIdx);
689 const Scalar p = fluidState.pressure(phaseIdx);
690
691 if (phaseIdx == liquidPhaseIdx)
692 {
693 // the liquid enthalpy is constant
694 return H2O::liquidEnthalpy(T, p);
695 }
696 else if (phaseIdx == gasPhaseIdx)
697 {
698 if (componentIdx == H2OIdx)
699 {
700 return H2O::gasEnthalpy(T, p);
701 }
702 else if (componentIdx == AirIdx)
703 {
704 return Air::gasEnthalpy(T, p);
705 }
706 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx);
707 }
708 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
709 }
710
721 template <class FluidState>
722 static Scalar thermalConductivity(const FluidState &fluidState,
723 int phaseIdx)
724 {
725 assert(0 <= phaseIdx && phaseIdx < numPhases);
726
727 const Scalar temperature = fluidState.temperature(phaseIdx) ;
728 const Scalar pressure = fluidState.pressure(phaseIdx);
729 if (phaseIdx == liquidPhaseIdx)
730 {
731 return H2O::liquidThermalConductivity(temperature, pressure);
732 }
733 else if (phaseIdx == gasPhaseIdx)
734 {
736 }
737 else
738 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
739 }
740
741 using Base::heatCapacity;
752 template <class FluidState>
753 static Scalar heatCapacity(const FluidState &fluidState,
754 int phaseIdx)
755 {
756 const Scalar temperature = fluidState.temperature(phaseIdx);
757 const Scalar pressure = fluidState.pressure(phaseIdx);
758 if (phaseIdx == liquidPhaseIdx)
759 {
760 // influence of air is neglected
761 return H2O::liquidHeatCapacity(temperature, pressure);
762 }
763 else if (phaseIdx == gasPhaseIdx)
764 {
765 return Air::gasHeatCapacity(temperature, pressure) * fluidState.moleFraction(gasPhaseIdx, AirIdx)
766 + H2O::gasHeatCapacity(temperature, pressure) * fluidState.moleFraction(gasPhaseIdx, H2OIdx);
767 }
768 else
769 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
770 }
771};
772
773} // end namespace FluidSystems
774} // end namespace Dumux
775
776#endif
Some exceptions thrown in DuMux
Some templates to wrap the valgrind macros.
A collection of input/output field names for common physical quantities.
Binary coefficients for water and air.
A simple class for the air fluid properties.
Material properties of pure water .
Tabulates all thermodynamic properties of a given untabulated chemical species.
Relations valid for an ideal gas.
bool CheckDefined(const T &value)
Make valgrind complain if the object occupied by an object is undefined.
Definition: valgrind.hh:72
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
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:46
static Scalar gasDensity(Scalar temperature, Scalar pressure)
The density of Air at a given pressure and temperature.
Definition: air.hh:84
static constexpr Scalar molarMass()
The molar mass in of Air.
Definition: air.hh:61
static const Scalar gasHeatCapacity(Scalar temperature, Scalar pressure)
Specific isobaric heat capacity of pure air.
Definition: air.hh:305
static Scalar criticalPressure()
Returns the critical pressure of Air.
Definition: air.hh:73
static Scalar gasViscosity(Scalar temperature, Scalar pressure)
The dynamic viscosity of Air at a given pressure and temperature.
Definition: air.hh:186
static Scalar gasThermalConductivity(Scalar temperature, Scalar pressure)
Thermal conductivity of air.
Definition: air.hh:342
static constexpr bool gasIsIdeal()
Returns true, the gas phase is assumed to be ideal.
Definition: air.hh:108
static Scalar gasEnthalpy(Scalar temperature, Scalar pressure)
Specific enthalpy of Air with 273.15 as basis.
Definition: air.hh:268
static Scalar criticalTemperature()
Returns the critical temperature of Air.
Definition: air.hh:67
static std::string name()
A human readable name for Air.
Definition: air.hh:53
static Scalar gasMolarDensity(Scalar temperature, Scalar pressure)
The molar density of air in , depending on pressure and temperature.
Definition: air.hh:96
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:52
static constexpr bool useAirViscosityAsGasMixtureViscosity()
Definition: h2oair.hh:55
static constexpr bool useIdealGasDensity()
Definition: h2oair.hh:54
static constexpr bool useH2ODensityAsLiquidMixtureDensity()
Definition: h2oair.hh:53
A compositional two-phase fluid system with water and air as components in both, the liquid and the g...
Definition: h2oair.hh:75
static Scalar viscosity(const FluidState &fluidState, int phaseIdx)
Calculate the dynamic viscosity of a fluid phase .
Definition: h2oair.hh:462
static constexpr bool isGas(int phaseIdx)
Return whether a phase is gaseous.
Definition: h2oair.hh:126
static Scalar enthalpy(const FluidState &fluidState, int phaseIdx)
Given a phase's composition, temperature and pressure, return its specific enthalpy .
Definition: h2oair.hh:648
static Scalar criticalMolarVolume(int compIdx)
Molar volume of a component at the critical point .
Definition: h2oair.hh:291
static Scalar vaporPressure(const FluidState &fluidState, int compIdx)
Vapor pressure of a component .
Definition: h2oair.hh:262
static constexpr int comp1Idx
index of the second component
Definition: h2oair.hh:95
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:180
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:344
static std::string componentName(int compIdx)
Return the human readable name of a component.
Definition: h2oair.hh:198
static constexpr int numComponents
Number of components in the fluid system.
Definition: h2oair.hh:85
static constexpr int liquidPhaseIdx
index of the liquid phase
Definition: h2oair.hh:87
static constexpr int gasPhaseIdx
index of the gas phase
Definition: h2oair.hh:88
static Scalar heatCapacity(const FluidState &fluidState, int phaseIdx)
Specific isobaric heat capacity of a fluid phase. .
Definition: h2oair.hh:753
static std::string phaseName(int phaseIdx)
Return the human readable name of a phase.
Definition: h2oair.hh:104
static Scalar thermalConductivity(const FluidState &fluidState, int phaseIdx)
Thermal conductivity of a fluid phase .
Definition: h2oair.hh:722
static Scalar componentEnthalpy(const FluidState &fluidState, int phaseIdx, int componentIdx)
Returns the specific enthalpy of a component in a specific phase.
Definition: h2oair.hh:684
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:146
static Scalar criticalPressure(int compIdx)
Critical pressure of a component .
Definition: h2oair.hh:244
static constexpr int phase1Idx
index of the second phase
Definition: h2oair.hh:90
static constexpr int comp0Idx
index of the frist component
Definition: h2oair.hh:94
static constexpr int numPhases
Number of phases in the fluid system.
Definition: h2oair.hh:84
static constexpr bool isCompressible(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: h2oair.hh:164
static constexpr int H2OIdx
index of the frist component
Definition: h2oair.hh:92
H2Otype H2O
Definition: h2oair.hh:81
static Scalar fugacityCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Returns the fugacity coefficient of a component in a phase.
Definition: h2oair.hh:535
static Scalar relativeHumidity(const FluidState &fluidState)
Returns the relative humidity of the gas phase.
Definition: h2oair.hh:563
static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
The molar density of a fluid phase in .
Definition: h2oair.hh:423
static constexpr int liquidCompIdx
index of the liquid component
Definition: h2oair.hh:96
static void init()
Initialize the fluid system's static parameters generically.
Definition: h2oair.hh:323
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:374
static constexpr int AirIdx
index of the second component
Definition: h2oair.hh:93
static Scalar acentricFactor(int compIdx)
The acentric factor of a component .
Definition: h2oair.hh:302
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:590
static constexpr int phase0Idx
index of the first phase
Definition: h2oair.hh:89
static constexpr int gasCompIdx
index of the gas component
Definition: h2oair.hh:97
static constexpr bool isMiscible()
Returns whether the fluids are miscible.
Definition: h2oair.hh:118
static Scalar criticalTemperature(int compIdx)
Critical temperature of a component .
Definition: h2oair.hh:228
static Scalar diffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Definition: h2oair.hh:571
static Scalar molarMass(int compIdx)
Return the molar mass of a component .
Definition: h2oair.hh:213
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.