H2ON2FluidSystem.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_N2_FLUID_SYSTEM_HPP
26 #define OPM_H2O_N2_FLUID_SYSTEM_HPP
27 
28 #include "BaseFluidSystem.hpp"
29 #include "NullParameterCache.hpp"
30 
38 
39 #include <opm/common/Exceptions.hpp>
40 #include <opm/common/ErrorMacros.hpp>
41 
42 #include <iostream>
43 #include <cassert>
44 
45 namespace Opm {
46 namespace FluidSystems {
47 
53 template <class Scalar, bool useComplexRelations = true>
54 class H2ON2
55  : public BaseFluidSystem<Scalar, H2ON2<Scalar, useComplexRelations> >
56 {
59 
60  // convenience typedefs
62  typedef Opm::H2O<Scalar> IapwsH2O;
64  typedef Opm::N2<Scalar> SimpleN2;
65 
66 public:
69 
70  /****************************************
71  * Fluid phase related static parameters
72  ****************************************/
73 
75  static const int numPhases = 2;
76 
78  static const int liquidPhaseIdx = 0;
80  static const int gasPhaseIdx = 1;
81 
83  static const char *phaseName(unsigned phaseIdx)
84  {
85  static const char *name[] = {
86  "liquid",
87  "gas"
88  };
89 
90  assert(0 <= phaseIdx && phaseIdx < numPhases);
91  return name[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  // gases are always compressible
106  return
107  (phaseIdx == gasPhaseIdx)
108  ? true
109  :H2O::liquidIsCompressible();// the water component decides for the liquid phase...
110  }
111 
113  static bool isIdealGas(unsigned phaseIdx)
114  {
115  //assert(0 <= phaseIdx && phaseIdx < numPhases);
116 
117  return
118  (phaseIdx == gasPhaseIdx)
119  ? H2O::gasIsIdeal() && N2::gasIsIdeal() // let the components decide
120  : false; // not a gas
121  }
122 
124  static bool isIdealMixture(unsigned /*phaseIdx*/)
125  {
126  //assert(0 <= phaseIdx && phaseIdx < numPhases);
127  // we assume Henry's and Rault's laws for the water phase and
128  // and no interaction between gas molecules of different
129  // components, so all phases are ideal mixtures!
130  return true;
131  }
132 
133  /****************************************
134  * Component related static parameters
135  ****************************************/
136 
138  static const int numComponents = 2;
139 
141  static const int H2OIdx = 0;
143  static const int N2Idx = 1;
144 
146  typedef TabulatedH2O H2O;
147  //typedef SimpleH2O H2O;
148  //typedef IapwsH2O H2O;
149 
151  typedef SimpleN2 N2;
152 
154  static const char *componentName(unsigned compIdx)
155  {
156  static const char *name[] = {
157  H2O::name(),
158  N2::name()
159  };
160 
161  assert(0 <= compIdx && compIdx < numComponents);
162  return name[compIdx];
163  }
164 
166  static Scalar molarMass(unsigned compIdx)
167  {
168  //assert(0 <= compIdx && compIdx < numComponents);
169  return (compIdx == H2OIdx)
170  ? H2O::molarMass()
171  : (compIdx == N2Idx)
172  ? N2::molarMass()
173  : 1e100;
174  }
175 
181  static Scalar criticalTemperature(unsigned compIdx)
182  {
183  return (compIdx == H2OIdx)
185  : (compIdx == N2Idx)
187  : 1e100;
188  }
189 
195  static Scalar criticalPressure(unsigned compIdx)
196  {
197  return (compIdx == H2OIdx)
199  : (compIdx == N2Idx)
201  : 1e100;
202  }
203 
209  static Scalar acentricFactor(unsigned compIdx)
210  {
211  return (compIdx == H2OIdx)
213  : (compIdx == N2Idx)
214  ? N2::acentricFactor()
215  : 1e100;
216  }
217 
218  /****************************************
219  * thermodynamic relations
220  ****************************************/
221 
228  static void init()
229  {
230  init(/*tempMin=*/273.15,
231  /*tempMax=*/623.15,
232  /*numTemp=*/100,
233  /*pMin=*/0.0,
234  /*pMax=*/20e6,
235  /*numP=*/200);
236  }
237 
249  static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp,
250  Scalar pressMin, Scalar pressMax, unsigned nPress)
251  {
252  if (H2O::isTabulated) {
253  TabulatedH2O::init(tempMin, tempMax, nTemp,
254  pressMin, pressMax, nPress);
255  }
256  }
257 
265  template <class FluidState, class LhsEval = typename FluidState::Scalar>
266  static LhsEval density(const FluidState &fluidState,
267  const ParameterCache &/*paramCache*/,
268  unsigned phaseIdx)
269  {
270  assert(0 <= phaseIdx && phaseIdx < numPhases);
271 
272  typedef Opm::MathToolbox<LhsEval> LhsToolbox;
274 
275  const auto& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
276  const auto& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
277 
278  LhsEval sumMoleFrac = 0;
279  for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
280  sumMoleFrac += FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(phaseIdx, compIdx));
281 
282  // liquid phase
283  if (phaseIdx == liquidPhaseIdx) {
284  if (!useComplexRelations)
285  // assume pure water
286  return H2O::liquidDensity(T, p);
287  else
288  {
289  // See: Ochs 2008
290  const auto& rholH2O = H2O::liquidDensity(T, p);
291  const auto& clH2O = rholH2O/H2O::molarMass();
292 
293  const auto& xlH2O = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(liquidPhaseIdx, H2OIdx));
294  const auto& xlN2 = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(liquidPhaseIdx, N2Idx));
295 
296  // this assumes each nitrogen molecule displaces exactly one
297  // water molecule in the liquid
298  return clH2O*(H2O::molarMass()*xlH2O + N2::molarMass()*xlN2)/sumMoleFrac;
299  }
300  }
301 
302  // gas phase
303  assert(phaseIdx == gasPhaseIdx);
304 
305  if (!useComplexRelations)
306  // for the gas phase assume an ideal gas
307  return
309  * FsToolbox::template toLhs<LhsEval>(fluidState.averageMolarMass(gasPhaseIdx))
310  / LhsToolbox::max(1e-5, sumMoleFrac);
311 
312  // assume ideal mixture: steam and nitrogen don't "see" each
313  // other
314  const auto& xgH2O = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(gasPhaseIdx, H2OIdx));
315  const auto& xgN2 = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(gasPhaseIdx, N2Idx));
316  const auto& rho_gH2O = H2O::gasDensity(T, p*xgH2O);
317  const auto& rho_gN2 = N2::gasDensity(T, p*xgN2);
318  return (rho_gH2O + rho_gN2)/LhsToolbox::max(1e-5, sumMoleFrac);
319  }
320 
322  template <class FluidState, class LhsEval = typename FluidState::Scalar>
323  static LhsEval viscosity(const FluidState &fluidState,
324  const ParameterCache &/*paramCache*/,
325  unsigned phaseIdx)
326  {
327  assert(0 <= phaseIdx && phaseIdx < numPhases);
328 
329  typedef Opm::MathToolbox<LhsEval> LhsToolbox;
331 
332  const auto& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
333  const auto& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
334 
335  // liquid phase
336  if (phaseIdx == liquidPhaseIdx)
337  // assume pure water for the liquid phase
338  return H2O::liquidViscosity(T, p);
339 
340  // gas phase
341  assert(phaseIdx == gasPhaseIdx);
342 
343  if (!useComplexRelations)
344  // assume pure nitrogen for the gas phase
345  return N2::gasViscosity(T, p);
346  else {
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
351  * 5th edition, McGraw-Hill, 20001, p. 9.21/22
352  */
353  LhsEval muResult = 0;
354  const LhsEval mu[numComponents] = {
356  N2::gasViscosity(T, p)
357  };
358 
359  LhsEval sumx = 0.0;
360  for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
361  sumx += FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(phaseIdx, compIdx));
362  sumx = LhsToolbox::max(1e-10, sumx);
363 
364  for (unsigned i = 0; i < numComponents; ++i) {
365  LhsEval divisor = 0;
366  for (unsigned j = 0; j < numComponents; ++j) {
367  LhsEval phiIJ = 1 + LhsToolbox::sqrt(mu[i]/mu[j]) * std::pow(molarMass(j)/molarMass(i), 1/4.0);
368  phiIJ *= phiIJ;
369  phiIJ /= std::sqrt(8*(1 + molarMass(i)/molarMass(j)));
370  divisor +=
371  FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(phaseIdx, j))
372  /sumx*phiIJ;
373  }
374  muResult +=
375  FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(phaseIdx, i))
376  /sumx*mu[i]/divisor;
377  }
378  return muResult;
379  }
380  }
381 
383  template <class FluidState, class LhsEval = typename FluidState::Scalar>
384  static LhsEval fugacityCoefficient(const FluidState &fluidState,
385  const ParameterCache &/*paramCache*/,
386  unsigned phaseIdx,
387  unsigned compIdx)
388  {
389  assert(0 <= phaseIdx && phaseIdx < numPhases);
390  assert(0 <= compIdx && compIdx < numComponents);
391 
393 
394  const auto& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
395  const auto& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
396 
397  // liquid phase
398  if (phaseIdx == liquidPhaseIdx) {
399  if (compIdx == H2OIdx)
400  return H2O::vaporPressure(T)/p;
402  }
403 
404  assert(phaseIdx == gasPhaseIdx);
405 
406  // for the gas phase, assume an ideal gas when it comes to
407  // fugacity (-> fugacity == partial pressure)
408  return 1.0;
409  }
410 
412  template <class FluidState, class LhsEval = typename FluidState::Scalar>
413  static LhsEval diffusionCoefficient(const FluidState &fluidState,
414  const ParameterCache &/*paramCache*/,
415  unsigned phaseIdx,
416  unsigned /*compIdx*/)
417 
418  {
420 
421  const auto& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
422  const auto& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
423 
424  // liquid phase
425  if (phaseIdx == liquidPhaseIdx)
427 
428  // gas phase
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  // liquid phase
447  if (phaseIdx == liquidPhaseIdx) {
448  // TODO: correct way to deal with the solutes???
449  return H2O::liquidEnthalpy(T, p);
450  }
451 
452  // gas phase
453  assert(phaseIdx == gasPhaseIdx);
454 
455  // assume ideal mixture: Molecules of one component don't
456  // "see" the molecules of the other component, which means
457  // that the total specific enthalpy is the sum of the
458  // "partial specific enthalpies" of the components.
459  const auto& XgH2O = FsToolbox::template toLhs<LhsEval>(fluidState.massFraction(gasPhaseIdx, H2OIdx));
460  const auto& XgN2 = FsToolbox::template toLhs<LhsEval>(fluidState.massFraction(gasPhaseIdx, N2Idx));
461 
462  LhsEval hH2O = XgH2O*H2O::gasEnthalpy(T, p);
463  LhsEval hN2 = XgN2*N2::gasEnthalpy(T, p);
464  return hH2O + hN2;
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  {
473  assert(0 <= phaseIdx && phaseIdx < numPhases);
474 
476 
477  const auto& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
478  const auto& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
479  if (phaseIdx == liquidPhaseIdx) // liquid phase
480  return H2O::liquidThermalConductivity(T, p);
481 
482  // gas phase
483  assert(phaseIdx == gasPhaseIdx);
484 
485  if (useComplexRelations){
486  // return the sum of the partial conductivity of Nitrogen and Steam
487  const auto& xH2O = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(phaseIdx, H2OIdx));
488  const auto& xN2 = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(phaseIdx, N2Idx));
489 
490  // Assuming Raoult's, Daltons law and ideal gas in order to obtain the
491  // partial pressures in the gas phase
492  const auto& lambdaN2 = N2::gasThermalConductivity(T, p*xN2);
493  const auto& lambdaH2O = H2O::gasThermalConductivity(T, p*xH2O);
494 
495  return lambdaN2 + lambdaH2O;
496  }
497  else
498  // return the conductivity of dry Nitrogen
499  return N2::gasThermalConductivity(T, p);
500  }
501 
503  template <class FluidState, class LhsEval = typename FluidState::Scalar>
504  static LhsEval heatCapacity(const FluidState &fluidState,
505  const ParameterCache &/*paramCache*/,
506  unsigned phaseIdx)
507  {
509 
510  const auto& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
511  const auto& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
512  const auto& xAlphaH2O = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(phaseIdx, H2OIdx));
513  const auto& xAlphaN2 = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(phaseIdx, N2Idx));
514  const auto& XAlphaH2O = FsToolbox::template toLhs<LhsEval>(fluidState.massFraction(phaseIdx, H2OIdx));
515  const auto& XAlphaN2 = FsToolbox::template toLhs<LhsEval>(fluidState.massFraction(phaseIdx, N2Idx));
516 
517  if (phaseIdx == liquidPhaseIdx)
518  return H2O::liquidHeatCapacity(T, p);
519 
520  assert(phaseIdx == gasPhaseIdx);
521 
522  // for the gas phase, assume ideal mixture, i.e. molecules of
523  // one component don't "see" the molecules of the other
524  // component
525  LhsEval c_pN2;
526  LhsEval c_pH2O;
527  // let the water and nitrogen components do things their own way
528  if (useComplexRelations) {
529  c_pN2 = N2::gasHeatCapacity(T, p*xAlphaN2);
530  c_pH2O = H2O::gasHeatCapacity(T, p*xAlphaH2O);
531  }
532  else {
533  // assume an ideal gas for both components. See:
534  //
535  // http://en.wikipedia.org/wiki/Heat_capacity
536  Scalar c_vN2molar = Opm::Constants<Scalar>::R*2.39;
537  Scalar c_pN2molar = Opm::Constants<Scalar>::R + c_vN2molar;
538 
539  Scalar c_vH2Omolar = Opm::Constants<Scalar>::R*3.37; // <- correct??
540  Scalar c_pH2Omolar = Opm::Constants<Scalar>::R + c_vH2Omolar;
541 
542  c_pN2 = c_pN2molar/molarMass(N2Idx);
543  c_pH2O = c_pH2Omolar/molarMass(H2OIdx);
544  }
545 
546  // mingle both components together. this assumes that there is no "cross
547  // interaction" between both flavors of molecules.
548  return XAlphaH2O*c_pH2O + XAlphaN2*c_pN2;
549  }
550 };
551 
552 } // namespace FluidSystems
553 
554 } // namespace Opm
555 
556 #endif
static bool isIdealMixture(unsigned)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: H2ON2FluidSystem.hpp:124
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: H2ON2FluidSystem.hpp:249
static LhsEval diffusionCoefficient(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: H2ON2FluidSystem.hpp:413
static const char * name()
A human readable name for the component.
Definition: TabulatedComponent.hpp:216
static Scalar molarMass()
The molar mass in of molecular nitrogen.
Definition: N2.hpp:63
static Scalar criticalTemperature()
Returns the critical temperature of molecular nitrogen.
Definition: N2.hpp:69
static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp, Scalar pressMin, Scalar pressMax, unsigned nPress)
Initialize the tables.
Definition: TabulatedComponent.hpp:73
Material properties of pure water .
Definition: H2O.hpp:60
static const char * componentName(unsigned compIdx)
Return the human readable name of a component.
Definition: H2ON2FluidSystem.hpp:154
Properties of pure molecular nitrogen .
Definition: N2.hpp:49
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
static Evaluation vaporPressure(const Evaluation &temperature)
The vapor pressure in of the component at a given temperature.
Definition: TabulatedComponent.hpp:256
static Evaluation gasHeatCapacity(const Evaluation &temperature, const Evaluation &pressure)
Specific isobaric heat capacity of the gas .
Definition: TabulatedComponent.hpp:311
static Scalar criticalPressure()
Returns the critical pressure of molecular nitrogen.
Definition: N2.hpp:75
static Evaluation gasViscosity(const Evaluation &temperature, const Evaluation &)
The dynamic viscosity of at a given pressure and temperature.
Definition: N2.hpp:269
static bool isLiquid(unsigned phaseIdx)
Return whether a phase is liquid.
Definition: H2ON2FluidSystem.hpp:95
static const int liquidPhaseIdx
Index of the liquid phase.
Definition: H2ON2FluidSystem.hpp:78
Definition: MathToolbox.hpp:39
Definition: Air_Mesitylene.hpp:31
static Evaluation liquidDensity(const Evaluation &temperature, const Evaluation &pressure)
The density of liquid at a given pressure and temperature .
Definition: TabulatedComponent.hpp:449
NullParameterCache ParameterCache
The type of the fluid system's parameter cache.
Definition: H2ON2FluidSystem.hpp:68
static Scalar acentricFactor(unsigned compIdx)
The acentric factor of a component [].
Definition: H2ON2FluidSystem.hpp:209
A parameter cache which does nothing.
static Evaluation gasDensity(const Evaluation &temperature, const Evaluation &pressure)
The density of gas at a given pressure and temperature.
Definition: N2.hpp:138
Binary coefficients for water and nitrogen.
static LhsEval viscosity(const FluidState &fluidState, const ParameterCache &, unsigned phaseIdx)
Calculate the dynamic viscosity of a fluid phase [Pa*s].
Definition: H2ON2FluidSystem.hpp:323
static Evaluation gasThermalConductivity(const Evaluation &temperature, const Evaluation &pressure)
The thermal conductivity of gaseous water .
Definition: TabulatedComponent.hpp:506
static Evaluation henry(const Evaluation &temperature)
Henry coefficent for molecular nitrogen in liquid water.
Definition: H2O_N2.hpp:52
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 Evaluation liquidHeatCapacity(const Evaluation &temperature, const Evaluation &pressure)
Specific isobaric heat capacity of the liquid .
Definition: TabulatedComponent.hpp:330
TabulatedH2O H2O
The component for pure water.
Definition: H2ON2FluidSystem.hpp:146
The base class for all fluid systems.
Definition: BaseFluidSystem.hpp:43
A parameter cache which does nothing.
Definition: NullParameterCache.hpp:36
Properties of pure molecular nitrogen .
static Scalar molarMass()
The molar mass in of the component.
Definition: TabulatedComponent.hpp:222
static Evaluation gasEnthalpy(const Evaluation &temperature, const Evaluation &)
Specific enthalpy of pure nitrogen gas.
Definition: N2.hpp:179
static bool isIdealGas(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: H2ON2FluidSystem.hpp:113
static Scalar criticalPressure(unsigned compIdx)
Critical pressure of a component [Pa].
Definition: H2ON2FluidSystem.hpp:195
static Scalar molarMass(unsigned compIdx)
Return the molar mass of a component in [kg/mol].
Definition: H2ON2FluidSystem.hpp:166
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: H2ON2FluidSystem.hpp:384
static bool gasIsIdeal()
Returns true iff the gas phase is assumed to be ideal.
Definition: N2.hpp:153
static const char * phaseName(unsigned phaseIdx)
Return the human readable name of a fluid phase.
Definition: H2ON2FluidSystem.hpp:83
Material properties of pure water .
static LhsEval heatCapacity(const FluidState &fluidState, const ParameterCache &, unsigned phaseIdx)
Specific isobaric heat capacity of a fluid phase [J/kg].
Definition: H2ON2FluidSystem.hpp:504
static bool isCompressible(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: H2ON2FluidSystem.hpp:102
static Evaluation liquidDiffCoeff(const Evaluation &temperature, const Evaluation &)
Diffusion coefficent for molecular nitrogen in liquid water.
Definition: H2O_N2.hpp:102
A two-phase fluid system with water and nitrogen as components.
Definition: H2ON2FluidSystem.hpp:54
static Scalar criticalTemperature(unsigned compIdx)
Critical temperature of a component [K].
Definition: H2ON2FluidSystem.hpp:181
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 gasThermalConductivity(const Evaluation &, const Evaluation &)
Specific heat conductivity of steam .
Definition: N2.hpp:307
static Scalar criticalPressure()
Returns the critical pressure in of the component.
Definition: TabulatedComponent.hpp:234
A simple version of pure water.
static const int H2OIdx
The component index of water.
Definition: H2ON2FluidSystem.hpp:141
static const int N2Idx
The component index of molecular nitrogen.
Definition: H2ON2FluidSystem.hpp:143
static Scalar criticalTemperature()
Returns the critical temperature in of the component.
Definition: TabulatedComponent.hpp:228
static bool gasIsIdeal()
Returns true iff the gas phase is assumed to be ideal.
Definition: TabulatedComponent.hpp:417
A generic class which tabulates all thermodynamic properties of a given component.
Definition: TabulatedComponent.hpp:56
static LhsEval density(const FluidState &fluidState, const ParameterCache &, unsigned phaseIdx)
Calculate the density [kg/m^3] of a fluid phase.
Definition: H2ON2FluidSystem.hpp:266
static Evaluation liquidEnthalpy(const Evaluation &temperature, const Evaluation &pressure)
Specific enthalpy of the liquid .
Definition: TabulatedComponent.hpp:292
static const int numComponents
Number of chemical species in the fluid system.
Definition: H2ON2FluidSystem.hpp:138
static const bool isTabulated
Definition: TabulatedComponent.hpp:61
static const Scalar acentricFactor()
The acentric factor of water.
Definition: H2O.hpp:85
Relations valid for an ideal gas.
A generic class which tabulates all thermodynamic properties of a given component.
static const int gasPhaseIdx
Index of the gas phase.
Definition: H2ON2FluidSystem.hpp:80
static LhsEval thermalConductivity(const FluidState &fluidState, const ParameterCache &, unsigned phaseIdx)
Thermal conductivity of a fluid phase [W/(m K)].
Definition: H2ON2FluidSystem.hpp:469
static bool liquidIsCompressible()
Returns true iff the liquid phase is assumed to be compressible.
Definition: TabulatedComponent.hpp:411
Relations valid for an ideal gas.
Definition: IdealGas.hpp:35
static Evaluation gasEnthalpy(const Evaluation &temperature, const Evaluation &pressure)
Specific enthalpy of the gas .
Definition: TabulatedComponent.hpp:273
static const int numPhases
Number of fluid phases in the fluid system.
Definition: H2ON2FluidSystem.hpp:75
Evaluation< Scalar, VarSetTag, numVars > pow(const Evaluation< Scalar, VarSetTag, numVars > &base, Scalar exp)
Definition: Math.hpp:312
static Evaluation gasHeatCapacity(const Evaluation &temperature, const Evaluation &)
Specific isobaric heat capacity of pure nitrogen gas.
Definition: N2.hpp:237
static Evaluation gasDensity(const Evaluation &temperature, const Evaluation &pressure)
The density of gas at a given pressure and temperature .
Definition: TabulatedComponent.hpp:429
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: H2ON2FluidSystem.hpp:435
static void init()
Initialize the fluid system's static parameters.
Definition: H2ON2FluidSystem.hpp:228
The base class for all fluid systems.
static Evaluation gasViscosity(const Evaluation &temperature, const Evaluation &pressure)
The dynamic viscosity of gas.
Definition: TabulatedComponent.hpp:468
SimpleN2 N2
The component for pure nitrogen.
Definition: H2ON2FluidSystem.hpp:151
static Evaluation liquidViscosity(const Evaluation &temperature, const Evaluation &pressure)
The dynamic viscosity of liquid.
Definition: TabulatedComponent.hpp:487
static Evaluation liquidThermalConductivity(const Evaluation &temperature, const Evaluation &pressure)
The thermal conductivity of liquid water .
Definition: TabulatedComponent.hpp:525
A central place for various physical constants occuring in some equations.
Definition: Constants.hpp:39
static Evaluation gasDiffCoeff(const Evaluation &temperature, const Evaluation &pressure)
Binary diffusion coefficent for molecular water and nitrogen.
Definition: H2O_N2.hpp:70
static const char * name()
A human readable name for nitrogen.
Definition: N2.hpp:57