H2OAirFluidSystem.hpp
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  Copyright (C) 2009-2013 by Andreas Lauser
5 
6  This file is part of the Open Porous Media project (OPM).
7 
8  OPM is free software: you can redistribute it and/or modify
9  it under the terms of the GNU General Public License as published by
10  the Free Software Foundation, either version 2 of the License, or
11  (at your option) any later version.
12 
13  OPM is distributed in the hope that it will be useful,
14  but WITHOUT ANY WARRANTY; without even the implied warranty of
15  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  GNU General Public License for more details.
17 
18  You should have received a copy of the GNU General Public License
19  along with OPM. If not, see <http://www.gnu.org/licenses/>.
20 */
25 #ifndef OPM_H2O_AIR_SYSTEM_HPP
26 #define OPM_H2O_AIR_SYSTEM_HPP
27 
28 #include "BaseFluidSystem.hpp"
29 #include "NullParameterCache.hpp"
30 
37 
38 #include <opm/common/Exceptions.hpp>
39 #include <opm/common/ErrorMacros.hpp>
40 
41 #include <iostream>
42 #include <cassert>
43 
44 namespace Opm {
45 namespace FluidSystems {
46 
56 template <class Scalar,
57  //class H2Otype = Opm::SimpleH2O<Scalar>,
59  bool useComplexRelations = true>
60 class H2OAir
61  : public BaseFluidSystem<Scalar, H2OAir<Scalar, H2Otype, useComplexRelations> >
62 {
66 
67 public:
70 
72  typedef H2Otype H2O;
75 
77  static const int numPhases = 2;
78 
80  static const int liquidPhaseIdx = 0;
82  static const int gasPhaseIdx = 1;
83 
85  static const char *phaseName(unsigned phaseIdx)
86  {
87  switch (phaseIdx) {
88  case liquidPhaseIdx: return "liquid";
89  case gasPhaseIdx: return "gas";
90  };
91  OPM_THROW(std::logic_error, "Invalid phase index " << phaseIdx);
92  }
93 
95  static bool isLiquid(unsigned phaseIdx)
96  {
97  //assert(0 <= phaseIdx && phaseIdx < numPhases);
98  return phaseIdx != gasPhaseIdx;
99  }
100 
102  static bool isCompressible(unsigned phaseIdx)
103  {
104  //assert(0 <= phaseIdx && phaseIdx < numPhases);
105  return (phaseIdx == gasPhaseIdx)
106  // ideal gases are always compressible
107  ? true
108  :
109  // the water component decides for the liquid phase...
111  }
112 
114  static bool isIdealGas(unsigned phaseIdx)
115  {
116  return
117  (phaseIdx == gasPhaseIdx)
119  : false;
120  }
121 
123  static bool isIdealMixture(unsigned /*phaseIdx*/)
124  {
125  //assert(0 <= phaseIdx && phaseIdx < numPhases);
126  // we assume Henry's and Rault's laws for the water phase and
127  // and no interaction between gas molecules of different
128  // components, so all phases are ideal mixtures!
129  return true;
130  }
131 
132  /****************************************
133  * Component related static parameters
134  ****************************************/
135 
137  static const int numComponents = 2;
138 
140  static const int H2OIdx = 0;
142  static const int AirIdx = 1;
143 
145  static const char *componentName(unsigned compIdx)
146  {
147  switch (compIdx)
148  {
149  case H2OIdx: return H2O::name();
150  case AirIdx: return Air::name();
151  };
152  OPM_THROW(std::logic_error, "Invalid component index " << compIdx);
153  }
154 
156  static Scalar molarMass(unsigned compIdx)
157  {
158  return
159  (compIdx == H2OIdx)
160  ? H2O::molarMass()
161  : (compIdx == AirIdx)
162  ? Air::molarMass()
163  : 1e100;
164  }
165 
166 
172  static Scalar criticalTemperature(unsigned compIdx)
173  {
174  return
175  (compIdx == H2OIdx)
177  : (compIdx == AirIdx)
179  : 1e100;
180  }
181 
187  static Scalar criticalPressure(unsigned compIdx)
188  {
189  return
190  (compIdx == H2OIdx)
192  : (compIdx == AirIdx)
194  : 1e100;
195  }
196 
202  static Scalar acentricFactor(unsigned compIdx)
203  {
204  return
205  (compIdx == H2OIdx)
207  : (compIdx == AirIdx)
208  ? Air::acentricFactor()
209  : 1e100;
210  }
211 
212  /****************************************
213  * thermodynamic relations
214  ****************************************/
215 
222  static void init()
223  {
224  if (H2O::isTabulated)
225  init(/*tempMin=*/273.15,
226  /*tempMax=*/623.15,
227  /*numTemp=*/100,
228  /*pMin=*/-10,
229  /*pMax=*/20e6,
230  /*numP=*/200);
231  }
232 
244  static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp,
245  Scalar pressMin, Scalar pressMax, unsigned nPress)
246  {
247  if (H2O::isTabulated) {
248  H2O::init(tempMin, tempMax, nTemp,
249  pressMin, pressMax, nPress);
250  }
251  }
252 
254  template <class FluidState, class LhsEval = typename FluidState::Scalar>
255  static LhsEval density(const FluidState &fluidState,
256  const ParameterCache &/*paramCache*/,
257  unsigned phaseIdx)
258  {
260  typedef Opm::MathToolbox<LhsEval> LhsToolbox;
261 
262  assert(0 <= phaseIdx && phaseIdx < numPhases);
263 
264  const auto& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
265  LhsEval p;
266  if (isCompressible(phaseIdx))
267  p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
268  else {
269  // random value which will hopefully cause things to blow
270  // up if it is used in a calculation!
271  p = - 1e100;
273  }
274 
275 
276  LhsEval sumMoleFrac = 0;
277  for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
278  sumMoleFrac += FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(phaseIdx, compIdx));
279 
280  if (phaseIdx == liquidPhaseIdx)
281  {
282  if (!useComplexRelations)
283  // assume pure water
284  return H2O::liquidDensity(T, p);
285  else
286  {
287  // See: Ochs 2008 (2.6)
288  const LhsEval& rholH2O = H2O::liquidDensity(T, p);
289  const LhsEval& clH2O = rholH2O/H2O::molarMass();
290 
291  const auto& xlH2O = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(liquidPhaseIdx, H2OIdx));
292  const auto& xlAir = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(liquidPhaseIdx, AirIdx));
293 
294  return clH2O*(H2O::molarMass()*xlH2O + Air::molarMass()*xlAir)/sumMoleFrac;
295  }
296  }
297  else if (phaseIdx == gasPhaseIdx)
298  {
299  if (!useComplexRelations)
300  // for the gas phase assume an ideal gas
301  return
303  * FsToolbox::template toLhs<LhsEval>(fluidState.averageMolarMass(gasPhaseIdx))
304  / LhsToolbox::max(1e-5, sumMoleFrac);
305 
306  LhsEval partialPressureH2O =
307  FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(gasPhaseIdx, H2OIdx))
308  *FsToolbox::template toLhs<LhsEval>(fluidState.pressure(gasPhaseIdx));
309 
310  LhsEval partialPressureAir =
311  FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(gasPhaseIdx, AirIdx))
312  *FsToolbox::template toLhs<LhsEval>(fluidState.pressure(gasPhaseIdx));
313 
314  return H2O::gasDensity(T, partialPressureH2O) + Air::gasDensity(T, partialPressureAir);
315  }
316  OPM_THROW(std::logic_error, "Invalid phase index " << phaseIdx);
317  }
318 
320  template <class FluidState, class LhsEval = typename FluidState::Scalar>
321  static LhsEval viscosity(const FluidState &fluidState,
322  const ParameterCache &/*paramCache*/,
323  unsigned phaseIdx)
324  {
325  typedef Opm::MathToolbox<LhsEval> LhsToolbox;
327 
328  assert(0 <= phaseIdx && phaseIdx < numPhases);
329 
330  const auto& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
331  const auto& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
332 
333  if (phaseIdx == liquidPhaseIdx)
334  {
335  // assume pure water for the liquid phase
336  // TODO: viscosity of mixture
337  // couldn't find a way to solve the mixture problem
338  return H2O::liquidViscosity(T, p);
339  }
340  else if (phaseIdx == gasPhaseIdx)
341  {
342  if(!useComplexRelations){
343  return Air::gasViscosity(T, p);
344  }
345  else //using a complicated version of this fluid system
346  {
347  /* Wilke method. See:
348  *
349  * See: R. Reid, et al.: The Properties of Gases and Liquids,
350  * 4th edition, McGraw-Hill, 1987, 407-410 or
351  * 5th edition, McGraw-Hill, 2000, p. 9.21/22
352  *
353  */
354 
355  LhsEval muResult = 0;
356  const LhsEval mu[numComponents] = {
358  Air::gasViscosity(T, p)
359  };
360 
361  // molar masses
362  const Scalar M[numComponents] = {
363  H2O::molarMass(),
365  };
366 
367  for (unsigned i = 0; i < numComponents; ++i) {
368  LhsEval divisor = 0;
369  for (unsigned j = 0; j < numComponents; ++j) {
370  LhsEval phiIJ =
371  1 +
372  LhsToolbox::sqrt(mu[i]/mu[j]) * // 1 + (mu[i]/mu[j]^1/2
373  std::pow(M[j]/M[i], 1./4.0); // (M[i]/M[j])^1/4
374 
375  phiIJ *= phiIJ;
376  phiIJ /= std::sqrt(8*(1 + M[i]/M[j]));
377  divisor += FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(phaseIdx, j))*phiIJ;
378  }
379  const auto& xAlphaI = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(phaseIdx, i));
380  muResult += xAlphaI*mu[i]/divisor;
381  }
382  return muResult;
383  }
384  }
385  OPM_THROW(std::logic_error, "Invalid phase index " << phaseIdx);
386  }
387 
389  template <class FluidState, class LhsEval = typename FluidState::Scalar>
390  static LhsEval fugacityCoefficient(const FluidState &fluidState,
391  const ParameterCache &/*paramCache*/,
392  unsigned phaseIdx,
393  unsigned compIdx)
394  {
396 
397  assert(0 <= phaseIdx && phaseIdx < numPhases);
398  assert(0 <= compIdx && compIdx < numComponents);
399 
400  const auto& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
401  const auto& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
402 
403  if (phaseIdx == liquidPhaseIdx) {
404  if (compIdx == H2OIdx)
405  return H2O::vaporPressure(T)/p;
407  }
408 
409  // for the gas phase, assume an ideal gas when it comes to
410  // fugacity (-> fugacity == partial pressure)
411  return 1.0;
412  }
413 
415  template <class FluidState, class LhsEval = typename FluidState::Scalar>
416  static LhsEval binaryDiffusionCoefficient(const FluidState &fluidState,
417  const ParameterCache &/*paramCache*/,
418  unsigned phaseIdx,
419  unsigned /*compIdx*/)
420  {
422 
423  const auto& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
424  const auto& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
425 
426  if (phaseIdx == liquidPhaseIdx)
428 
429  assert(phaseIdx == gasPhaseIdx);
431  }
432 
434  template <class FluidState, class LhsEval = typename FluidState::Scalar>
435  static LhsEval enthalpy(const FluidState &fluidState,
436  const ParameterCache &/*paramCache*/,
437  unsigned phaseIdx)
438  {
440 
441  const auto& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
442  const auto& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
445 
446  if (phaseIdx == liquidPhaseIdx)
447  {
448  // TODO: correct way to deal with the solutes???
449  return H2O::liquidEnthalpy(T, p);
450  }
451 
452  else if (phaseIdx == gasPhaseIdx)
453  {
454  LhsEval result = 0.0;
455  result +=
456  H2O::gasEnthalpy(T, p) *
457  FsToolbox::template toLhs<LhsEval>(fluidState.massFraction(gasPhaseIdx, H2OIdx));
458 
459  result +=
460  Air::gasEnthalpy(T, p) *
461  FsToolbox::template toLhs<LhsEval>(fluidState.massFraction(gasPhaseIdx, AirIdx));
462  return result;
463  }
464  OPM_THROW(std::logic_error, "Invalid phase index " << phaseIdx);
465  }
466 
468  template <class FluidState, class LhsEval = typename FluidState::Scalar>
469  static LhsEval thermalConductivity(const FluidState &fluidState,
470  const ParameterCache &/*paramCache*/,
471  unsigned phaseIdx)
472  {
474 
475  assert(0 <= phaseIdx && phaseIdx < numPhases);
476 
477  const LhsEval& temperature =
478  FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
479  const LhsEval& pressure =
480  FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
481 
482  if (phaseIdx == liquidPhaseIdx)
483  return H2O::liquidThermalConductivity(temperature, pressure);
484  else { // gas phase
485  const LhsEval& lambdaDryAir = Air::gasThermalConductivity(temperature, pressure);
486 
487  if (useComplexRelations){
488  const LhsEval& xAir =
489  FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(phaseIdx, AirIdx));
490  const LhsEval& xH2O =
491  FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(phaseIdx, H2OIdx));
492  LhsEval lambdaAir = xAir*lambdaDryAir;
493 
494  // Assuming Raoult's, Daltons law and ideal gas
495  // in order to obtain the partial density of water in the air phase
496  LhsEval partialPressure = pressure*xH2O;
497 
498  LhsEval lambdaH2O =
499  xH2O*H2O::gasThermalConductivity(temperature, partialPressure);
500  return lambdaAir + lambdaH2O;
501  }
502  else
503  return lambdaDryAir; // conductivity of Nitrogen [W / (m K ) ]
504  }
505  }
506 };
507 
508 } // namespace FluidSystems
509 
510 } // namespace Opm
511 
512 #endif
static bool isCompressible(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: H2OAirFluidSystem.hpp:102
void SetUndefined(const T &value OPM_UNUSED)
Make the memory on which an object resides undefined in valgrind runs.
Definition: Valgrind.hpp:138
static bool gasIsIdeal()
Returns true iff the gas phase is assumed to be ideal.
Definition: Air.hpp:73
static bool isIdealGas(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: H2OAirFluidSystem.hpp:114
static Evaluation liquidDiffCoeff(const Evaluation &temperature, const Evaluation &)
Diffusion coefficent for molecular nitrogen in liquid water.
Definition: H2O_Air.hpp:107
bool CheckDefined(const T &value OPM_UNUSED)
Make valgrind complain if any of the memory occupied by an object is undefined.
Definition: Valgrind.hpp:74
Definition: MathToolbox.hpp:39
static Evaluation vaporPressure(Evaluation temperature)
The vapor pressure in of pure water at a given temperature.
Definition: H2O.hpp:131
Definition: Air_Mesitylene.hpp:31
static Scalar criticalPressure()
Returns the critical pressure of .
Definition: Air.hpp:93
Opm::Air< Scalar > Air
The type of the air component used for this fluid system.
Definition: H2OAirFluidSystem.hpp:74
static LhsEval enthalpy(const FluidState &fluidState, const ParameterCache &, unsigned phaseIdx)
Given a phase's composition, temperature, pressure and density, calculate its specific enthalpy [J/kg...
Definition: H2OAirFluidSystem.hpp:435
static Evaluation gasDiffCoeff(const Evaluation &temperature, const Evaluation &pressure)
Binary diffusion coefficent for molecular water and air.
Definition: H2O_Air.hpp:74
A parameter cache which does nothing.
static const int liquidPhaseIdx
The index of the liquid phase.
Definition: H2OAirFluidSystem.hpp:80
static Scalar criticalTemperature()
Returns the critical temperature of .
Definition: Air.hpp:87
static const char * name()
A human readable name for the .
Definition: Air.hpp:61
Some templates to wrap the valgrind client request macros.
Evaluation< Scalar, VarSetTag, numVars > max(const Evaluation< Scalar, VarSetTag, numVars > &x1, const Evaluation< Scalar, VarSetTag, numVars > &x2)
Definition: Math.hpp:114
Evaluation< Scalar, VarSetTag, numVars > sqrt(const Evaluation< Scalar, VarSetTag, numVars > &x)
Definition: Math.hpp:278
static void init(Scalar, Scalar, unsigned, Scalar, Scalar, unsigned)
A default routine for initialization, not needed for components and must not be called.
Definition: Component.hpp:61
static Evaluation liquidViscosity(const Evaluation &temperature, const Evaluation &pressure)
The dynamic viscosity of pure water.
Definition: H2O.hpp:815
The base class for all fluid systems.
Definition: BaseFluidSystem.hpp:43
A parameter cache which does nothing.
Definition: NullParameterCache.hpp:36
static Evaluation gasDensity(const Evaluation &temperature, const Evaluation &pressure)
The density of steam in at a given pressure and temperature.
Definition: H2O.hpp:556
H2Otype H2O
The type of the water component used for this fluid system.
Definition: H2OAirFluidSystem.hpp:72
static Evaluation gasThermalConductivity(const Evaluation &temperature, const Evaluation &pressure)
Thermal conductivity of water (IAPWS) .
Definition: H2O.hpp:862
A simple class implementing the fluid properties of air.
static const int gasPhaseIdx
The index of the gas phase.
Definition: H2OAirFluidSystem.hpp:82
static Scalar criticalTemperature(unsigned compIdx)
Critical temperature of a component [K].
Definition: H2OAirFluidSystem.hpp:172
static const int numComponents
Number of chemical species in the fluid system.
Definition: H2OAirFluidSystem.hpp:137
NullParameterCache ParameterCache
The type of the fluid system's parameter cache.
Definition: H2OAirFluidSystem.hpp:69
static const bool isTabulated
Definition: Component.hpp:47
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: H2OAirFluidSystem.hpp:244
static bool isLiquid(unsigned phaseIdx)
Return whether a phase is liquid.
Definition: H2OAirFluidSystem.hpp:95
Material properties of pure water .
static Evaluation molarDensity(const Evaluation &temperature, const Evaluation &pressure)
The molar density of the gas , depending on pressure and temperature.
Definition: IdealGas.hpp:65
static Evaluation gasViscosity(const Evaluation &temperature, const Evaluation &)
The dynamic viscosity of at a given pressure and temperature.
Definition: Air.hpp:138
static Evaluation gasDensity(const Evaluation &temperature, const Evaluation &pressure)
The density of at a given pressure and temperature [kg/m^3].
Definition: Air.hpp:103
static Evaluation liquidDensity(const Evaluation &temperature, const Evaluation &pressure)
The density of pure water in at a given pressure and temperature.
Definition: H2O.hpp:685
static bool gasIsIdeal()
Returns true iff the gas phase is assumed to be ideal.
Definition: H2O.hpp:623
static bool liquidIsCompressible()
Returns true iff the liquid phase is assumed to be compressible.
Definition: H2O.hpp:540
static LhsEval thermalConductivity(const FluidState &fluidState, const ParameterCache &, unsigned phaseIdx)
Thermal conductivity of a fluid phase [W/(m K)].
Definition: H2OAirFluidSystem.hpp:469
static Evaluation henry(const Evaluation &temperature)
Henry coefficent for air in liquid water.
Definition: H2O_Air.hpp:55
A fluid system with a liquid and a gaseous phase and water and air as components. ...
Definition: H2OAirFluidSystem.hpp:60
static Scalar molarMass()
The molar mass in of .
Definition: Air.hpp:81
A generic class which tabulates all thermodynamic properties of a given component.
Definition: TabulatedComponent.hpp:56
static Scalar molarMass(unsigned compIdx)
Return the molar mass of a component in [kg/mol].
Definition: H2OAirFluidSystem.hpp:156
static Evaluation gasThermalConductivity(const Evaluation &, const Evaluation &)
Specific heat conductivity of steam .
Definition: Air.hpp:223
static const int AirIdx
The index of the air component.
Definition: H2OAirFluidSystem.hpp:142
static const Scalar acentricFactor()
The acentric factor of water.
Definition: H2O.hpp:85
Relations valid for an ideal gas.
static Evaluation liquidEnthalpy(const Evaluation &temperature, const Evaluation &pressure)
Specific enthalpy of liquid water .
Definition: H2O.hpp:229
A generic class which tabulates all thermodynamic properties of a given component.
static const char * componentName(unsigned compIdx)
Return the human readable name of a component.
Definition: H2OAirFluidSystem.hpp:145
static const Scalar criticalPressure()
Returns the critical pressure of water.
Definition: H2O.hpp:97
static const Scalar criticalTemperature()
Returns the critical temperature of water.
Definition: H2O.hpp:91
static LhsEval fugacityCoefficient(const FluidState &fluidState, const ParameterCache &, unsigned phaseIdx, unsigned compIdx)
Calculate the fugacity coefficient [Pa] of an individual component in a fluid phase.
Definition: H2OAirFluidSystem.hpp:390
Relations valid for an ideal gas.
Definition: IdealGas.hpp:35
static Evaluation gasViscosity(const Evaluation &temperature, const Evaluation &pressure)
The dynamic viscosity of steam.
Definition: H2O.hpp:790
static LhsEval density(const FluidState &fluidState, const ParameterCache &, unsigned phaseIdx)
Calculate the density [kg/m^3] of a fluid phase.
Definition: H2OAirFluidSystem.hpp:255
static const char * name()
A human readable name for the water.
Definition: H2O.hpp:73
static Scalar acentricFactor(unsigned compIdx)
The acentric factor of a component [].
Definition: H2OAirFluidSystem.hpp:202
Evaluation< Scalar, VarSetTag, numVars > pow(const Evaluation< Scalar, VarSetTag, numVars > &base, Scalar exp)
Definition: Math.hpp:312
static LhsEval binaryDiffusionCoefficient(const FluidState &fluidState, const ParameterCache &, unsigned phaseIdx, unsigned)
Calculate the binary molecular diffusion coefficient for a component in a fluid phase [mol^2 * s / (k...
Definition: H2OAirFluidSystem.hpp:416
static const int H2OIdx
The index of the water component.
Definition: H2OAirFluidSystem.hpp:140
static bool isIdealMixture(unsigned)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: H2OAirFluidSystem.hpp:123
static Evaluation gasEnthalpy(const Evaluation &temperature, const Evaluation &)
Specific enthalpy of liquid water with 273.15 K as basis. See: W. Kays, M. Crawford, B. Weigand Convective heat and mass transfer, 4th edition (2005) p. 431ff.
Definition: Air.hpp:186
static Evaluation gasEnthalpy(const Evaluation &temperature, const Evaluation &pressure)
Specific enthalpy of water steam .
Definition: H2O.hpp:176
static Scalar criticalPressure(unsigned compIdx)
Critical pressure of a component [Pa].
Definition: H2OAirFluidSystem.hpp:187
static void init()
Initialize the fluid system's static parameters.
Definition: H2OAirFluidSystem.hpp:222
static LhsEval viscosity(const FluidState &fluidState, const ParameterCache &, unsigned phaseIdx)
Calculate the dynamic viscosity of a fluid phase [Pa*s].
Definition: H2OAirFluidSystem.hpp:321
static const char * phaseName(unsigned phaseIdx)
Return the human readable name of a fluid phase.
Definition: H2OAirFluidSystem.hpp:85
static Evaluation liquidThermalConductivity(const Evaluation &temperature, const Evaluation &pressure)
Thermal conductivity of water (IAPWS) .
Definition: H2O.hpp:842
Binary coefficients for water and nitrogen.
The base class for all fluid systems.
A simple class implementing the fluid properties of air.
Definition: Air.hpp:47
static const int numPhases
Number of fluid phases in the fluid system.
Definition: H2OAirFluidSystem.hpp:77
static const Scalar molarMass()
The molar mass in of water.
Definition: H2O.hpp:79