3.6-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
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 viscosityIsConstant(int phaseIdx)
180 {
181 // water decides for the liquid phase
182 if (phaseIdx == liquidPhaseIdx)
183 return H2O::liquidViscosityIsConstant();
184 // air decides if policy is enabled
185 else if (phaseIdx == gasPhaseIdx && Policy::useAirViscosityAsGasMixtureViscosity())
187 // in general it depends on the mixture
188 else
189 return false;
190 }
191
198 static constexpr bool isIdealGas(int phaseIdx)
199 {
200 assert(0 <= phaseIdx && phaseIdx < numPhases);
201
202 // let the fluids decide
203 if (phaseIdx == gasPhaseIdx)
204 return H2O::gasIsIdeal() && Air::gasIsIdeal();
205 return false; // not a gas
206 }
207
208 /****************************************
209 * Component related static parameters
210 ****************************************/
216 static std::string componentName(int compIdx)
217 {
218 switch (compIdx)
219 {
220 case H2OIdx: return H2O::name();
221 case AirIdx: return Air::name();
222 }
223 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx);
224 }
225
231 static Scalar molarMass(int compIdx)
232 {
233 switch (compIdx)
234 {
235 case H2OIdx: return H2O::molarMass();
236 case AirIdx: return Air::molarMass();
237 }
238 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx);
239 }
240
246 static Scalar criticalTemperature(int compIdx)
247 {
248 static const Scalar TCrit[] = {
249 H2O::criticalTemperature(),
251 };
252
253 assert(0 <= compIdx && compIdx < numComponents);
254 return TCrit[compIdx];
255 }
256
262 static Scalar criticalPressure(int compIdx)
263 {
264 static const Scalar pCrit[] = {
265 H2O::criticalPressure(),
267 };
268
269 assert(0 <= compIdx && compIdx < numComponents);
270 return pCrit[compIdx];
271 }
272
279 template <class FluidState>
280 static Scalar vaporPressure(const FluidState &fluidState, int compIdx)
281 {
282 if (compIdx == H2OIdx)
283 {
284 const auto t = fluidState.temperature(H2OIdx);
285 if (!useKelvinVaporPressure)
286 return H2O::vaporPressure(t);
287 else
288 {
289 const auto pc = (fluidState.wettingPhase() == (int) H2OIdx)
290 ? fluidState.pressure(AirIdx)-fluidState.pressure(H2OIdx)
291 : fluidState.pressure(H2OIdx)-fluidState.pressure(AirIdx);
292 return H2O::vaporPressure(t)*exp( -pc * molarMass(H2OIdx)
293 / density(fluidState, H2OIdx)
295 }
296 }
297 else if (compIdx == AirIdx)
298 // return Air::vaporPressure(fluidState.temperature(AirIdx));
299 DUNE_THROW(Dune::NotImplemented, "Air::vaporPressure(t)");
300 else
301 DUNE_THROW(Dune::NotImplemented, "Invalid component index " << compIdx);
302 }
303
309 static Scalar criticalMolarVolume(int compIdx)
310 {
311 DUNE_THROW(Dune::NotImplemented,
312 "H2OAirFluidSystem::criticalMolarVolume()");
313 }
314
320 static Scalar acentricFactor(int compIdx)
321 {
322 static const Scalar accFac[] = {
323 H2O::acentricFactor(),
324 Air::acentricFactor()
325 };
326
327 assert(0 <= compIdx && compIdx < numComponents);
328 return accFac[compIdx];
329 }
330
331 /****************************************
332 * thermodynamic relations
333 ****************************************/
334
341 static void init()
342 {
343 init(/*tempMin=*/273.15,
344 /*tempMax=*/623.15,
345 /*numTemp=*/100,
346 /*pMin=*/-10.,
347 /*pMax=*/20e6,
348 /*numP=*/200);
349 }
350
362 static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp,
363 Scalar pressMin, Scalar pressMax, unsigned nPress)
364 {
365 std::cout << "The H2O-air fluid system was configured with the following policy:\n";
366 std::cout << " - use H2O density as liquid mixture density: " << std::boolalpha << Policy::useH2ODensityAsLiquidMixtureDensity() << "\n";
367 std::cout << " - use ideal gas density: " << std::boolalpha << Policy::useIdealGasDensity() << "\n";
368 std::cout << " - use air viscosity as gas mixture viscosity: " << std::boolalpha << Policy::useAirViscosityAsGasMixtureViscosity() << std::endl;
369
370 if (H2O::isTabulated)
371 {
372 H2O::init(tempMin, tempMax, nTemp,
373 pressMin, pressMax, nPress);
374 }
375 }
376
377 using Base::density;
391 template <class FluidState>
392 static Scalar density(const FluidState &fluidState,
393 const int phaseIdx)
394 {
395 assert(0 <= phaseIdx && phaseIdx < numPhases);
396
397 const Scalar T = fluidState.temperature(phaseIdx);
398 const Scalar p = fluidState.pressure(phaseIdx);
399
400 if (phaseIdx == phase0Idx)
401 {
402 if (Policy::useH2ODensityAsLiquidMixtureDensity())
403 // assume pure water
404 return H2O::liquidDensity(T, p);
405 else
406 {
407 // See: Eq. (7) in Class et al. (2002a)
408 // This assumes each gas molecule displaces exactly one
409 // molecule in the liquid.
410 return H2O::liquidMolarDensity(T, p)
411 * (H2O::molarMass()*fluidState.moleFraction(liquidPhaseIdx, H2OIdx)
412 + Air::molarMass()*fluidState.moleFraction(liquidPhaseIdx, AirIdx));
413 }
414 }
415 else if (phaseIdx == gasPhaseIdx)
416 {
417 if (Policy::useIdealGasDensity())
418 // for the gas phase assume an ideal gas
419 {
420 const Scalar averageMolarMass = fluidState.averageMolarMass(gasPhaseIdx);
421 return IdealGas::density(averageMolarMass, T, p);
422 }
423
424 return H2O::gasDensity(T, fluidState.partialPressure(gasPhaseIdx, H2OIdx))
425 + Air::gasDensity(T, fluidState.partialPressure(gasPhaseIdx, AirIdx));
426 }
427 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
428 }
429
430 using Base::molarDensity;
440 template <class FluidState>
441 static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
442 {
443 assert(0 <= phaseIdx && phaseIdx < numPhases);
444
445 const Scalar T = fluidState.temperature(phaseIdx);
446 const Scalar p = fluidState.pressure(phaseIdx);
447
448 if (phaseIdx == phase0Idx)
449 {
450 // assume pure water or that each gas molecule displaces exactly one
451 // molecule in the liquid.
452 return H2O::liquidMolarDensity(T, p);
453 }
454 else if (phaseIdx == phase1Idx)
455 {
456 if (Policy::useIdealGasDensity())
457 // for the gas phase assume an ideal gas
458 { return IdealGas::molarDensity(T, p); }
459
460 return H2O::gasMolarDensity(T, fluidState.partialPressure(phase1Idx, H2OIdx))
461 + Air::gasMolarDensity(T, fluidState.partialPressure(phase1Idx, AirIdx));
462 }
463 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
464 }
465
466 using Base::viscosity;
479 template <class FluidState>
480 static Scalar viscosity(const FluidState &fluidState,
481 int phaseIdx)
482 {
483 assert(0 <= phaseIdx && phaseIdx < numPhases);
484
485 Scalar T = fluidState.temperature(phaseIdx);
486 Scalar p = fluidState.pressure(phaseIdx);
487
488 if (phaseIdx == liquidPhaseIdx)
489 {
490 // assume pure water for the liquid phase
491 return H2O::liquidViscosity(T, p);
492 }
493 else if (phaseIdx == gasPhaseIdx)
494 {
495 if(Policy::useAirViscosityAsGasMixtureViscosity()){
496 return Air::gasViscosity(T, p);
497 }
498 else //using a complicated version of this fluid system
499 {
500 // Wilke method (Reid et al.):
501 Scalar muResult = 0;
502 const Scalar mu[numComponents] = {
505 };
506
507 // molar masses
508 const Scalar M[numComponents] = {
509 H2O::molarMass(),
511 };
512
513 for (int i = 0; i < numComponents; ++i)
514 {
515 Scalar divisor = 0;
516 using std::sqrt;
517 for (int j = 0; j < numComponents; ++j)
518 {
519 // 1 + (mu[i]/mu[j]^1/2 * (M[i]/M[j])^1/4)
520 Scalar phiIJ = 1 + sqrt(mu[i]/mu[j] * sqrt(M[j]/M[i]));
521 phiIJ = phiIJ * phiIJ / sqrt(8*(1 + M[i]/M[j]));
522 divisor += fluidState.moleFraction(phaseIdx, j)*phiIJ;
523 }
524 muResult += fluidState.moleFraction(phaseIdx, i)*mu[i] / divisor;
525 }
526 return muResult;
527 }
528 }
529 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
530 }
531
551 template <class FluidState>
552 static Scalar fugacityCoefficient(const FluidState &fluidState,
553 int phaseIdx,
554 int compIdx)
555 {
556 assert(0 <= phaseIdx && phaseIdx < numPhases);
557 assert(0 <= compIdx && compIdx < numComponents);
558
559 Scalar T = fluidState.temperature(phaseIdx);
560 Scalar p = fluidState.pressure(phaseIdx);
561
562 if (phaseIdx == liquidPhaseIdx) {
563 if (compIdx == H2OIdx)
564 return vaporPressure(fluidState, compIdx)/p;
565 return BinaryCoeff::H2O_Air::henry(T)/p;
566 }
567
568 // for the gas phase, assume an ideal gas when it comes to
569 // fugacity (-> fugacity == partial pressure)
570 return 1.0;
571 }
572
579 template <class FluidState>
580 static Scalar relativeHumidity(const FluidState &fluidState)
581 {
582 return fluidState.partialPressure(gasPhaseIdx, comp0Idx)
583 / H2O::vaporPressure(fluidState.temperature(gasPhaseIdx));
584 }
585
587 template <class FluidState>
588 static Scalar diffusionCoefficient(const FluidState &fluidState,
589 int phaseIdx,
590 int compIdx)
591 {
592 DUNE_THROW(Dune::NotImplemented, "FluidSystems::H2OAir::diffusionCoefficient()");
593 }
594
606 template <class FluidState>
607 static Scalar binaryDiffusionCoefficient(const FluidState &fluidState,
608 int phaseIdx,
609 int compIIdx,
610 int compJIdx)
611 {
612 using std::swap;
613 if (compIIdx > compJIdx)
614 swap(compIIdx, compJIdx);
615
616 Scalar T = fluidState.temperature(phaseIdx);
617 Scalar p = fluidState.pressure(phaseIdx);
618
619 // we are in the liquid phase
620 if (phaseIdx == liquidPhaseIdx)
621 {
622 if (compIIdx == H2OIdx && compJIdx == AirIdx)
624 else
625 DUNE_THROW(Dune::InvalidStateException,
626 "Binary diffusion coefficient of components "
627 << compIIdx << " and " << compJIdx
628 << " in phase " << phaseIdx << " is undefined!\n");
629 }
630
631 // we are in the gas phase
632 else if (phaseIdx == gasPhaseIdx)
633 {
634 if (compIIdx == H2OIdx && compJIdx == AirIdx)
636 else
637 DUNE_THROW(Dune::InvalidStateException,
638 "Binary diffusion coefficient of components "
639 << compIIdx << " and " << compJIdx
640 << " in phase " << phaseIdx << " is undefined!\n");
641 }
642
643 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
644 }
645
646 using Base::enthalpy;
664 template <class FluidState>
665 static Scalar enthalpy(const FluidState &fluidState,
666 int phaseIdx)
667 {
668 const Scalar T = fluidState.temperature(phaseIdx);
669 const Scalar p = fluidState.pressure(phaseIdx);
670
671 if (phaseIdx == liquidPhaseIdx)
672 return H2O::liquidEnthalpy(T, p);
673
674 else if (phaseIdx == gasPhaseIdx)
675 return H2O::gasEnthalpy(T, p)*fluidState.massFraction(gasPhaseIdx, H2OIdx)
676 + Air::gasEnthalpy(T, p)*fluidState.massFraction(gasPhaseIdx, AirIdx);
677
678 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
679 }
680
688 template <class FluidState>
689 static Scalar componentEnthalpy(const FluidState &fluidState,
690 int phaseIdx,
691 int componentIdx)
692 {
693 const Scalar T = fluidState.temperature(phaseIdx);
694 const Scalar p = fluidState.pressure(phaseIdx);
695
696 if (phaseIdx == liquidPhaseIdx)
697 {
698 // the liquid enthalpy is constant
699 return H2O::liquidEnthalpy(T, p);
700 }
701 else if (phaseIdx == gasPhaseIdx)
702 {
703 if (componentIdx == H2OIdx)
704 {
705 return H2O::gasEnthalpy(T, p);
706 }
707 else if (componentIdx == AirIdx)
708 {
709 return Air::gasEnthalpy(T, p);
710 }
711 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx);
712 }
713 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
714 }
715
726 template <class FluidState>
727 static Scalar thermalConductivity(const FluidState &fluidState,
728 int phaseIdx)
729 {
730 assert(0 <= phaseIdx && phaseIdx < numPhases);
731
732 const Scalar temperature = fluidState.temperature(phaseIdx) ;
733 const Scalar pressure = fluidState.pressure(phaseIdx);
734 if (phaseIdx == liquidPhaseIdx)
735 {
736 return H2O::liquidThermalConductivity(temperature, pressure);
737 }
738 else if (phaseIdx == gasPhaseIdx)
739 {
741 }
742 else
743 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
744 }
745
746 using Base::heatCapacity;
757 template <class FluidState>
758 static Scalar heatCapacity(const FluidState &fluidState,
759 int phaseIdx)
760 {
761 const Scalar temperature = fluidState.temperature(phaseIdx);
762 const Scalar pressure = fluidState.pressure(phaseIdx);
763 if (phaseIdx == liquidPhaseIdx)
764 {
765 // influence of air is neglected
766 return H2O::liquidHeatCapacity(temperature, pressure);
767 }
768 else if (phaseIdx == gasPhaseIdx)
769 {
770 return Air::gasHeatCapacity(temperature, pressure) * fluidState.moleFraction(gasPhaseIdx, AirIdx)
771 + H2O::gasHeatCapacity(temperature, pressure) * fluidState.moleFraction(gasPhaseIdx, H2OIdx);
772 }
773 else
774 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
775 }
776};
777
778} // end namespace FluidSystems
779} // end namespace Dumux
780
781#endif
Some exceptions thrown in DuMux
A collection of input/output field names for common physical quantities.
Relations valid for an ideal gas.
Material properties of pure water .
Tabulates all thermodynamic properties of a given untabulated chemical species.
A simple class for the air fluid properties.
Binary coefficients for water and air.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
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:973
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:313
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:193
static Scalar gasThermalConductivity(Scalar temperature, Scalar pressure)
Thermal conductivity of air.
Definition: air.hh:351
static constexpr bool gasViscosityIsConstant()
Returns true if the gas phase viscosity is constant.
Definition: air.hh:115
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:276
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 component.
Definition: tabulatedcomponent.hh:684
A central place for various physical constants occurring 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:480
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:665
static Scalar criticalMolarVolume(int compIdx)
Molar volume of a component at the critical point .
Definition: h2oair.hh:309
static Scalar vaporPressure(const FluidState &fluidState, int compIdx)
Vapor pressure of a component .
Definition: h2oair.hh:280
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:198
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:362
static std::string componentName(int compIdx)
Return the human readable name of a component.
Definition: h2oair.hh:216
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:758
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:727
static Scalar componentEnthalpy(const FluidState &fluidState, int phaseIdx, int componentIdx)
Returns the specific enthalpy of a component in a specific phase.
Definition: h2oair.hh:689
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:262
static constexpr int phase1Idx
index of the second phase
Definition: h2oair.hh:89
static constexpr int comp0Idx
index of the first 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 first component
Definition: h2oair.hh:91
static constexpr bool viscosityIsConstant(int phaseIdx)
Returns true if and only if a fluid phase is assumed to have a constant viscosity.
Definition: h2oair.hh:179
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:552
static Scalar relativeHumidity(const FluidState &fluidState)
Returns the relative humidity of the gas phase.
Definition: h2oair.hh:580
static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
The molar density of a fluid phase in .
Definition: h2oair.hh:441
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:341
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:392
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:320
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:607
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:246
static Scalar diffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Definition: h2oair.hh:588
static Scalar molarMass(int compIdx)
Return the molar mass of a component .
Definition: h2oair.hh:231
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.