All Classes Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
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  This file is part of the Open Porous Media project (OPM).
5 
6  OPM 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 2 of the License, or
9  (at your option) any later version.
10 
11  OPM 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 OPM. If not, see <http://www.gnu.org/licenses/>.
18 
19  Consult the COPYING file in the top-level source directory of this
20  module for the precise wording of the license and the list of
21  copyright holders.
22 */
27 #ifndef OPM_H2O_N2_FLUID_SYSTEM_HPP
28 #define OPM_H2O_N2_FLUID_SYSTEM_HPP
29 
30 #include "BaseFluidSystem.hpp"
31 #include "NullParameterCache.hpp"
32 
39 #include <opm/common/Valgrind.hpp>
40 
41 #include <opm/common/Exceptions.hpp>
42 #include <opm/common/ErrorMacros.hpp>
43 
44 #include <iostream>
45 #include <cassert>
46 
47 namespace Opm {
48 namespace FluidSystems {
49 
55 template <class Scalar, bool useComplexRelations = true>
56 class H2ON2
57  : public BaseFluidSystem<Scalar, H2ON2<Scalar, useComplexRelations> >
58 {
61 
62  // convenience typedefs
64  typedef Opm::H2O<Scalar> IapwsH2O;
66  typedef Opm::N2<Scalar> SimpleN2;
67 
68 public:
70  template <class Evaluation>
72 
73  /****************************************
74  * Fluid phase related static parameters
75  ****************************************/
76 
78  static const int numPhases = 2;
79 
81  static const int liquidPhaseIdx = 0;
83  static const int gasPhaseIdx = 1;
84 
86  static const char* phaseName(unsigned phaseIdx)
87  {
88  static const char* name[] = {
89  "liquid",
90  "gas"
91  };
92 
93  assert(0 <= phaseIdx && phaseIdx < numPhases);
94  return name[phaseIdx];
95  }
96 
98  static bool isLiquid(unsigned phaseIdx)
99  {
100  //assert(0 <= phaseIdx && phaseIdx < numPhases);
101  return phaseIdx != gasPhaseIdx;
102  }
103 
105  static bool isCompressible(unsigned phaseIdx)
106  {
107  //assert(0 <= phaseIdx && phaseIdx < numPhases);
108  // gases are always compressible
109  return
110  (phaseIdx == gasPhaseIdx)
111  ? true
112  :H2O::liquidIsCompressible();// the water component decides for the liquid phase...
113  }
114 
116  static bool isIdealGas(unsigned phaseIdx)
117  {
118  //assert(0 <= phaseIdx && phaseIdx < numPhases);
119 
120  return
121  (phaseIdx == gasPhaseIdx)
122  ? H2O::gasIsIdeal() && N2::gasIsIdeal() // let the components decide
123  : false; // not a gas
124  }
125 
127  static bool isIdealMixture(unsigned /*phaseIdx*/)
128  {
129  //assert(0 <= phaseIdx && phaseIdx < numPhases);
130  // we assume Henry's and Rault's laws for the water phase and
131  // and no interaction between gas molecules of different
132  // components, so all phases are ideal mixtures!
133  return true;
134  }
135 
136  /****************************************
137  * Component related static parameters
138  ****************************************/
139 
141  static const int numComponents = 2;
142 
144  static const int H2OIdx = 0;
146  static const int N2Idx = 1;
147 
149  typedef TabulatedH2O H2O;
150  //typedef SimpleH2O H2O;
151  //typedef IapwsH2O H2O;
152 
154  typedef SimpleN2 N2;
155 
157  static const char* componentName(unsigned compIdx)
158  {
159  static const char* name[] = {
160  H2O::name(),
161  N2::name()
162  };
163 
164  assert(0 <= compIdx && compIdx < numComponents);
165  return name[compIdx];
166  }
167 
169  static Scalar molarMass(unsigned compIdx)
170  {
171  //assert(0 <= compIdx && compIdx < numComponents);
172  return (compIdx == H2OIdx)
173  ? H2O::molarMass()
174  : (compIdx == N2Idx)
175  ? N2::molarMass()
176  : 1e30;
177  }
178 
184  static Scalar criticalTemperature(unsigned compIdx)
185  {
186  return (compIdx == H2OIdx)
188  : (compIdx == N2Idx)
190  : 1e30;
191  }
192 
198  static Scalar criticalPressure(unsigned compIdx)
199  {
200  return (compIdx == H2OIdx)
202  : (compIdx == N2Idx)
204  : 1e30;
205  }
206 
212  static Scalar acentricFactor(unsigned compIdx)
213  {
214  return (compIdx == H2OIdx)
216  : (compIdx == N2Idx)
217  ? N2::acentricFactor()
218  : 1e30;
219  }
220 
221  /****************************************
222  * thermodynamic relations
223  ****************************************/
224 
231  static void init()
232  {
233  init(/*tempMin=*/273.15,
234  /*tempMax=*/623.15,
235  /*numTemp=*/50,
236  /*pMin=*/0.0,
237  /*pMax=*/20e6,
238  /*numP=*/50);
239  }
240 
252  static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp,
253  Scalar pressMin, Scalar pressMax, unsigned nPress)
254  {
255  if (H2O::isTabulated) {
256  TabulatedH2O::init(tempMin, tempMax, nTemp,
257  pressMin, pressMax, nPress);
258  }
259  }
260 
268  template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
269  static LhsEval density(const FluidState& fluidState,
270  const ParameterCache<ParamCacheEval>& /*paramCache*/,
271  unsigned phaseIdx)
272  {
273  assert(0 <= phaseIdx && phaseIdx < numPhases);
274 
275  const auto& T = Opm::decay<LhsEval>(fluidState.temperature(phaseIdx));
276  const auto& p = Opm::decay<LhsEval>(fluidState.pressure(phaseIdx));
277 
278  LhsEval sumMoleFrac = 0;
279  for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
280  sumMoleFrac += Opm::decay<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 = Opm::decay<LhsEval>(fluidState.moleFraction(liquidPhaseIdx, H2OIdx));
294  const auto& xlN2 = Opm::decay<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  * Opm::decay<LhsEval>(fluidState.averageMolarMass(gasPhaseIdx))
310  / Opm::max(1e-5, sumMoleFrac);
311 
312  // assume ideal mixture: steam and nitrogen don't "see" each
313  // other
314  const auto& xgH2O = Opm::decay<LhsEval>(fluidState.moleFraction(gasPhaseIdx, H2OIdx));
315  const auto& xgN2 = Opm::decay<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)/Opm::max(1e-5, sumMoleFrac);
319  }
320 
322  template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
323  static LhsEval viscosity(const FluidState& fluidState,
324  const ParameterCache<ParamCacheEval>& /*paramCache*/,
325  unsigned phaseIdx)
326  {
327  assert(0 <= phaseIdx && phaseIdx < numPhases);
328 
329  const auto& T = Opm::decay<LhsEval>(fluidState.temperature(phaseIdx));
330  const auto& p = Opm::decay<LhsEval>(fluidState.pressure(phaseIdx));
331 
332  // liquid phase
333  if (phaseIdx == liquidPhaseIdx)
334  // assume pure water for the liquid phase
335  return H2O::liquidViscosity(T, p);
336 
337  // gas phase
338  assert(phaseIdx == gasPhaseIdx);
339 
340  if (!useComplexRelations)
341  // assume pure nitrogen for the gas phase
342  return N2::gasViscosity(T, p);
343  else {
344  /* Wilke method. See:
345  *
346  * See: R. Reid, et al.: The Properties of Gases and Liquids,
347  * 4th edition, McGraw-Hill, 1987, 407-410
348  * 5th edition, McGraw-Hill, 20001, p. 9.21/22
349  */
350  LhsEval muResult = 0;
351  const LhsEval mu[numComponents] = {
353  N2::gasViscosity(T, p)
354  };
355 
356  LhsEval sumx = 0.0;
357  for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
358  sumx += Opm::decay<LhsEval>(fluidState.moleFraction(phaseIdx, compIdx));
359  sumx = Opm::max(1e-10, sumx);
360 
361  for (unsigned i = 0; i < numComponents; ++i) {
362  LhsEval divisor = 0;
363  for (unsigned j = 0; j < numComponents; ++j) {
364  LhsEval phiIJ = 1 + Opm::sqrt(mu[i]/mu[j]) * std::pow(molarMass(j)/molarMass(i), 1/4.0);
365  phiIJ *= phiIJ;
366  phiIJ /= std::sqrt(8*(1 + molarMass(i)/molarMass(j)));
367  divisor +=
368  Opm::decay<LhsEval>(fluidState.moleFraction(phaseIdx, j))
369  /sumx*phiIJ;
370  }
371  muResult +=
372  Opm::decay<LhsEval>(fluidState.moleFraction(phaseIdx, i))
373  /sumx*mu[i]/divisor;
374  }
375  return muResult;
376  }
377  }
378 
380  template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
381  static LhsEval fugacityCoefficient(const FluidState& fluidState,
382  const ParameterCache<ParamCacheEval>& /*paramCache*/,
383  unsigned phaseIdx,
384  unsigned compIdx)
385  {
386  assert(0 <= phaseIdx && phaseIdx < numPhases);
387  assert(0 <= compIdx && compIdx < numComponents);
388 
389  const auto& T = Opm::decay<LhsEval>(fluidState.temperature(phaseIdx));
390  const auto& p = Opm::decay<LhsEval>(fluidState.pressure(phaseIdx));
391 
392  // liquid phase
393  if (phaseIdx == liquidPhaseIdx) {
394  if (compIdx == H2OIdx)
395  return H2O::vaporPressure(T)/p;
397  }
398 
399  assert(phaseIdx == gasPhaseIdx);
400 
401  // for the gas phase, assume an ideal gas when it comes to
402  // fugacity (-> fugacity == partial pressure)
403  return 1.0;
404  }
405 
407  template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
408  static LhsEval diffusionCoefficient(const FluidState& fluidState,
409  const ParameterCache<ParamCacheEval>& /*paramCache*/,
410  unsigned phaseIdx,
411  unsigned /*compIdx*/)
412 
413  {
414  const auto& T = Opm::decay<LhsEval>(fluidState.temperature(phaseIdx));
415  const auto& p = Opm::decay<LhsEval>(fluidState.pressure(phaseIdx));
416 
417  // liquid phase
418  if (phaseIdx == liquidPhaseIdx)
420 
421  // gas phase
422  assert(phaseIdx == gasPhaseIdx);
424  }
425 
427  template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
428  static LhsEval enthalpy(const FluidState& fluidState,
429  const ParameterCache<ParamCacheEval>& /*paramCache*/,
430  unsigned phaseIdx)
431  {
432  const auto& T = Opm::decay<LhsEval>(fluidState.temperature(phaseIdx));
433  const auto& p = Opm::decay<LhsEval>(fluidState.pressure(phaseIdx));
434  Valgrind::CheckDefined(T);
435  Valgrind::CheckDefined(p);
436 
437  // liquid phase
438  if (phaseIdx == liquidPhaseIdx) {
439  // TODO: correct way to deal with the solutes???
440  return H2O::liquidEnthalpy(T, p);
441  }
442 
443  // gas phase
444  assert(phaseIdx == gasPhaseIdx);
445 
446  // assume ideal mixture: Molecules of one component don't
447  // "see" the molecules of the other component, which means
448  // that the total specific enthalpy is the sum of the
449  // "partial specific enthalpies" of the components.
450  const auto& XgH2O = Opm::decay<LhsEval>(fluidState.massFraction(gasPhaseIdx, H2OIdx));
451  const auto& XgN2 = Opm::decay<LhsEval>(fluidState.massFraction(gasPhaseIdx, N2Idx));
452 
453  LhsEval hH2O = XgH2O*H2O::gasEnthalpy(T, p);
454  LhsEval hN2 = XgN2*N2::gasEnthalpy(T, p);
455  return hH2O + hN2;
456  }
457 
459  template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
460  static LhsEval thermalConductivity(const FluidState& fluidState,
461  const ParameterCache<ParamCacheEval>& /*paramCache*/,
462  unsigned phaseIdx)
463  {
464  assert(0 <= phaseIdx && phaseIdx < numPhases);
465 
466  const auto& T = Opm::decay<LhsEval>(fluidState.temperature(phaseIdx));
467  const auto& p = Opm::decay<LhsEval>(fluidState.pressure(phaseIdx));
468  if (phaseIdx == liquidPhaseIdx) // liquid phase
469  return H2O::liquidThermalConductivity(T, p);
470 
471  // gas phase
472  assert(phaseIdx == gasPhaseIdx);
473 
474  if (useComplexRelations){
475  // return the sum of the partial conductivity of Nitrogen and Steam
476  const auto& xH2O = Opm::decay<LhsEval>(fluidState.moleFraction(phaseIdx, H2OIdx));
477  const auto& xN2 = Opm::decay<LhsEval>(fluidState.moleFraction(phaseIdx, N2Idx));
478 
479  // Assuming Raoult's, Daltons law and ideal gas in order to obtain the
480  // partial pressures in the gas phase
481  const auto& lambdaN2 = N2::gasThermalConductivity(T, p*xN2);
482  const auto& lambdaH2O = H2O::gasThermalConductivity(T, p*xH2O);
483 
484  return lambdaN2 + lambdaH2O;
485  }
486  else
487  // return the conductivity of dry Nitrogen
488  return N2::gasThermalConductivity(T, p);
489  }
490 
492  template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
493  static LhsEval heatCapacity(const FluidState& fluidState,
494  const ParameterCache<ParamCacheEval>& /*paramCache*/,
495  unsigned phaseIdx)
496  {
497  const auto& T = Opm::decay<LhsEval>(fluidState.temperature(phaseIdx));
498  const auto& p = Opm::decay<LhsEval>(fluidState.pressure(phaseIdx));
499  const auto& xAlphaH2O = Opm::decay<LhsEval>(fluidState.moleFraction(phaseIdx, H2OIdx));
500  const auto& xAlphaN2 = Opm::decay<LhsEval>(fluidState.moleFraction(phaseIdx, N2Idx));
501  const auto& XAlphaH2O = Opm::decay<LhsEval>(fluidState.massFraction(phaseIdx, H2OIdx));
502  const auto& XAlphaN2 = Opm::decay<LhsEval>(fluidState.massFraction(phaseIdx, N2Idx));
503 
504  if (phaseIdx == liquidPhaseIdx)
505  return H2O::liquidHeatCapacity(T, p);
506 
507  assert(phaseIdx == gasPhaseIdx);
508 
509  // for the gas phase, assume ideal mixture, i.e. molecules of
510  // one component don't "see" the molecules of the other
511  // component
512  LhsEval c_pN2;
513  LhsEval c_pH2O;
514  // let the water and nitrogen components do things their own way
515  if (useComplexRelations) {
516  c_pN2 = N2::gasHeatCapacity(T, p*xAlphaN2);
517  c_pH2O = H2O::gasHeatCapacity(T, p*xAlphaH2O);
518  }
519  else {
520  // assume an ideal gas for both components. See:
521  //
522  // http://en.wikipedia.org/wiki/Heat_capacity
523  Scalar c_vN2molar = Opm::Constants<Scalar>::R*2.39;
524  Scalar c_pN2molar = Opm::Constants<Scalar>::R + c_vN2molar;
525 
526  Scalar c_vH2Omolar = Opm::Constants<Scalar>::R*3.37; // <- correct??
527  Scalar c_pH2Omolar = Opm::Constants<Scalar>::R + c_vH2Omolar;
528 
529  c_pN2 = c_pN2molar/molarMass(N2Idx);
530  c_pH2O = c_pH2Omolar/molarMass(H2OIdx);
531  }
532 
533  // mingle both components together. this assumes that there is no "cross
534  // interaction" between both flavors of molecules.
535  return XAlphaH2O*c_pH2O + XAlphaN2*c_pN2;
536  }
537 };
538 
539 } // namespace FluidSystems
540 
541 } // namespace Opm
542 
543 #endif
Material properties of pure water .
Definition: H2O.hpp:61
static const char * name()
A human readable name for nitrogen.
Definition: N2.hpp:56
static Evaluation henry(const Evaluation &temperature)
Henry coefficent for molecular nitrogen in liquid water.
Definition: H2O_N2.hpp:52
static Evaluation liquidDensity(const Evaluation &temperature, const Evaluation &pressure)
The density of liquid at a given pressure and temperature .
Definition: TabulatedComponent.hpp:435
static Evaluation liquidEnthalpy(const Evaluation &temperature, const Evaluation &pressure)
Specific enthalpy of the liquid .
Definition: TabulatedComponent.hpp:290
static bool isIdealGas(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: H2ON2FluidSystem.hpp:116
static LhsEval enthalpy(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Given a phase&#39;s composition, temperature, pressure and density, calculate its specific enthalpy [J/kg...
Definition: H2ON2FluidSystem.hpp:428
static Evaluation liquidHeatCapacity(const Evaluation &temperature, const Evaluation &pressure)
Specific isobaric heat capacity of the liquid .
Definition: TabulatedComponent.hpp:324
static Scalar molarMass()
The molar mass in of the component.
Definition: TabulatedComponent.hpp:224
static const char * phaseName(unsigned phaseIdx)
Return the human readable name of a fluid phase.
Definition: H2ON2FluidSystem.hpp:86
static const int N2Idx
The component index of molecular nitrogen.
Definition: H2ON2FluidSystem.hpp:146
A simple version of pure water.
static LhsEval diffusionCoefficient(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx, unsigned)
Calculate the binary molecular diffusion coefficient for a component in a fluid phase [mol^2 * s / (k...
Definition: H2ON2FluidSystem.hpp:408
Relations valid for an ideal gas.
SimpleN2 N2
The component for pure nitrogen.
Definition: H2ON2FluidSystem.hpp:154
Properties of pure molecular nitrogen .
Definition: N2.hpp:48
static LhsEval heatCapacity(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Specific isobaric heat capacity of a fluid phase [J/kg].
Definition: H2ON2FluidSystem.hpp:493
static const char * componentName(unsigned compIdx)
Return the human readable name of a component.
Definition: H2ON2FluidSystem.hpp:157
static LhsEval viscosity(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Calculate the dynamic viscosity of a fluid phase [Pa*s].
Definition: H2ON2FluidSystem.hpp:323
static const int gasPhaseIdx
Index of the gas phase.
Definition: H2ON2FluidSystem.hpp:83
static Scalar molarMass(unsigned compIdx)
Return the molar mass of a component in [kg/mol].
Definition: H2ON2FluidSystem.hpp:169
static const Scalar acentricFactor()
The acentric factor of water.
Definition: H2O.hpp:86
static const int H2OIdx
The component index of water.
Definition: H2ON2FluidSystem.hpp:144
static bool isLiquid(unsigned phaseIdx)
Return whether a phase is liquid.
Definition: H2ON2FluidSystem.hpp:98
static const int liquidPhaseIdx
Index of the liquid phase.
Definition: H2ON2FluidSystem.hpp:81
static Evaluation liquidDiffCoeff(const Evaluation &temperature, const Evaluation &)
Diffusion coefficent for molecular nitrogen in liquid water.
Definition: H2O_N2.hpp:102
static void init()
Initialize the fluid system&#39;s static parameters.
Definition: H2ON2FluidSystem.hpp:231
static Evaluation gasThermalConductivity(const Evaluation &, const Evaluation &)
Specific heat conductivity of steam .
Definition: N2.hpp:302
static Scalar criticalPressure()
Returns the critical pressure of molecular nitrogen.
Definition: N2.hpp:74
static Evaluation gasHeatCapacity(const Evaluation &temperature, const Evaluation &)
Specific isobaric heat capacity of pure nitrogen gas.
Definition: N2.hpp:234
static Evaluation gasDensity(const Evaluation &temperature, const Evaluation &pressure)
The density of gas at a given pressure and temperature.
Definition: N2.hpp:135
static Evaluation gasEnthalpy(const Evaluation &temperature, const Evaluation &)
Specific enthalpy of pure nitrogen gas.
Definition: N2.hpp:176
static LhsEval density(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Calculate the density [kg/m^3] of a fluid phase.
Definition: H2ON2FluidSystem.hpp:269
A parameter cache which does nothing.
Material properties of pure water .
static Evaluation gasViscosity(const Evaluation &temperature, const Evaluation &pressure)
The dynamic viscosity of gas.
Definition: TabulatedComponent.hpp:452
static Evaluation gasDensity(const Evaluation &temperature, const Evaluation &pressure)
The density of gas at a given pressure and temperature .
Definition: TabulatedComponent.hpp:417
A central place for various physical constants occuring in some equations.
Definition: Constants.hpp:40
A two-phase fluid system with water and nitrogen as components.
Definition: H2ON2FluidSystem.hpp:56
static Scalar criticalPressure(unsigned compIdx)
Critical pressure of a component [Pa].
Definition: H2ON2FluidSystem.hpp:198
Relations valid for an ideal gas.
Definition: IdealGas.hpp:37
static Evaluation liquidThermalConductivity(const Evaluation &temperature, const Evaluation &pressure)
The thermal conductivity of liquid water .
Definition: TabulatedComponent.hpp:503
static Scalar acentricFactor(unsigned compIdx)
The acentric factor of a component [].
Definition: H2ON2FluidSystem.hpp:212
static const char * name()
A human readable name for the component.
Definition: TabulatedComponent.hpp:218
static Evaluation vaporPressure(const Evaluation &temperature)
The vapor pressure in of the component at a given temperature.
Definition: TabulatedComponent.hpp:258
static Scalar molarMass()
The molar mass in of molecular nitrogen.
Definition: N2.hpp:62
static const int numComponents
Number of chemical species in the fluid system.
Definition: H2ON2FluidSystem.hpp:141
static Scalar criticalPressure()
Returns the critical pressure in of the component.
Definition: TabulatedComponent.hpp:236
static Scalar criticalTemperature(unsigned compIdx)
Critical temperature of a component [K].
Definition: H2ON2FluidSystem.hpp:184
static Evaluation gasViscosity(const Evaluation &temperature, const Evaluation &)
The dynamic viscosity of at a given pressure and temperature.
Definition: N2.hpp:266
static Evaluation gasHeatCapacity(const Evaluation &temperature, const Evaluation &pressure)
Specific isobaric heat capacity of the gas .
Definition: TabulatedComponent.hpp:307
static Scalar criticalTemperature()
Returns the critical temperature in of the component.
Definition: TabulatedComponent.hpp:230
static LhsEval thermalConductivity(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Thermal conductivity of a fluid phase [W/(m K)].
Definition: H2ON2FluidSystem.hpp:460
Properties of pure molecular nitrogen .
static LhsEval fugacityCoefficient(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx, unsigned compIdx)
Calculate the fugacity coefficient [Pa] of an individual component in a fluid phase.
Definition: H2ON2FluidSystem.hpp:381
A generic class which tabulates all thermodynamic properties of a given component.
The base class for all fluid systems.
Definition: BaseFluidSystem.hpp:43
A generic class which tabulates all thermodynamic properties of a given component.
Definition: TabulatedComponent.hpp:58
static bool isIdealMixture(unsigned)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: H2ON2FluidSystem.hpp:127
static Evaluation molarDensity(const Evaluation &temperature, const Evaluation &pressure)
The molar density of the gas , depending on pressure and temperature.
Definition: IdealGas.hpp:67
static bool isCompressible(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: H2ON2FluidSystem.hpp:105
static Evaluation gasDiffCoeff(const Evaluation &temperature, const Evaluation &pressure)
Binary diffusion coefficent for molecular water and nitrogen.
Definition: H2O_N2.hpp:70
A parameter cache which does nothing.
Definition: NullParameterCache.hpp:39
TabulatedH2O H2O
The component for pure water.
Definition: H2ON2FluidSystem.hpp:149
static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp, Scalar pressMin, Scalar pressMax, unsigned nPress)
Initialize the fluid system&#39;s static parameters using problem specific temperature and pressure range...
Definition: H2ON2FluidSystem.hpp:252
Scalar Scalar
The type used for scalar quantities.
Definition: BaseFluidSystem.hpp:49
static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp, Scalar pressMin, Scalar pressMax, unsigned nPress)
Initialize the tables.
Definition: TabulatedComponent.hpp:75
static bool gasIsIdeal()
Returns true iff the gas phase is assumed to be ideal.
Definition: N2.hpp:150
static Evaluation gasEnthalpy(const Evaluation &temperature, const Evaluation &pressure)
Specific enthalpy of the gas .
Definition: TabulatedComponent.hpp:273
static Scalar criticalTemperature()
Returns the critical temperature of molecular nitrogen.
Definition: N2.hpp:68
static const int numPhases
Number of fluid phases in the fluid system.
Definition: H2ON2FluidSystem.hpp:78
static bool gasIsIdeal()
Returns true iff the gas phase is assumed to be ideal.
Definition: TabulatedComponent.hpp:405
static Evaluation liquidViscosity(const Evaluation &temperature, const Evaluation &pressure)
The dynamic viscosity of liquid.
Definition: TabulatedComponent.hpp:469
static Evaluation gasThermalConductivity(const Evaluation &temperature, const Evaluation &pressure)
The thermal conductivity of gaseous water .
Definition: TabulatedComponent.hpp:486
The base class for all fluid systems.
static bool liquidIsCompressible()
Returns true iff the liquid phase is assumed to be compressible.
Definition: TabulatedComponent.hpp:399
Binary coefficients for water and nitrogen.