H2OAirXyleneFluidSystem.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) 2011-2013 by Andreas Lauser
5  Copyright (C) 2011 by Benjamin Faigle
6 
7  This file is part of the Open Porous Media project (OPM).
8 
9  OPM is free software: you can redistribute it and/or modify
10  it under the terms of the GNU General Public License as published by
11  the Free Software Foundation, either version 2 of the License, or
12  (at your option) any later version.
13 
14  OPM is distributed in the hope that it will be useful,
15  but WITHOUT ANY WARRANTY; without even the implied warranty of
16  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17  GNU General Public License for more details.
18 
19  You should have received a copy of the GNU General Public License
20  along with OPM. If not, see <http://www.gnu.org/licenses/>.
21 */
26 #ifndef OPM_H2O_AIR_XYLENE_FLUID_SYSTEM_HPP
27 #define OPM_H2O_AIR_XYLENE_FLUID_SYSTEM_HPP
28 
34 
38 
39 #include "BaseFluidSystem.hpp"
40 #include "NullParameterCache.hpp"
41 
42 namespace Opm {
43 namespace FluidSystems {
44 
50 template <class Scalar>
52  : public BaseFluidSystem<Scalar, H2OAirXylene<Scalar> >
53 {
56 
57 public:
60 
67 
69  static const int numPhases = 3;
71  static const int numComponents = 3;
72 
74  static const int waterPhaseIdx = 0;
76  static const int naplPhaseIdx = 1;
78  static const int gasPhaseIdx = 2;
79 
81  static const int H2OIdx = 0;
83  static const int NAPLIdx = 1;
85  static const int airIdx = 2;
86 
88  static void init()
89  { }
90 
92  static bool isLiquid(unsigned phaseIdx)
93  {
94  //assert(0 <= phaseIdx && phaseIdx < numPhases);
95  return phaseIdx != gasPhaseIdx;
96  }
97 
99  static bool isIdealGas(unsigned phaseIdx)
100  { return phaseIdx == gasPhaseIdx && H2O::gasIsIdeal() && Air::gasIsIdeal() && NAPL::gasIsIdeal(); }
101 
103  static bool isIdealMixture(unsigned /*phaseIdx*/)
104  {
105  //assert(0 <= phaseIdx && phaseIdx < numPhases);
106 
107  // we assume Henry's and Rault's laws for the water phase and
108  // and no interaction between gas molecules of different
109  // components, so all phases are ideal mixtures!
110  return true;
111  }
112 
114  static bool isCompressible(unsigned phaseIdx)
115  {
116  return
117  (phaseIdx == gasPhaseIdx)
118  // gases are always compressible
119  ? true
120  : (phaseIdx == waterPhaseIdx)
121  // the water component decides for the water phase...
123  // the NAPL component decides for the napl phase...
125  }
126 
128  static const char *phaseName(unsigned phaseIdx)
129  {
130  switch (phaseIdx) {
131  case waterPhaseIdx: return "water";
132  case naplPhaseIdx: return "napl";
133  case gasPhaseIdx: return "gas";
134  };
135  OPM_THROW(std::logic_error, "Invalid phase index " << phaseIdx);
136  }
137 
139  static const char *componentName(unsigned compIdx)
140  {
141  switch (compIdx) {
142  case H2OIdx: return H2O::name();
143  case airIdx: return Air::name();
144  case NAPLIdx: return NAPL::name();
145  };
146  OPM_THROW(std::logic_error, "Invalid component index " << compIdx);
147  }
148 
150  static Scalar molarMass(unsigned compIdx)
151  {
152  return
153  (compIdx == H2OIdx)
154  // gases are always compressible
155  ? H2O::molarMass()
156  : (compIdx == airIdx)
157  // the water component decides for the water comp...
158  ? Air::molarMass()
159  // the NAPL component decides for the napl comp...
160  : (compIdx == NAPLIdx)
161  ? NAPL::molarMass()
162  : 1e100;
163  }
164 
166  template <class FluidState, class LhsEval = typename FluidState::Scalar>
167  static LhsEval density(const FluidState &fluidState,
168  const ParameterCache &/*paramCache*/,
169  unsigned phaseIdx)
170  {
172 
173  if (phaseIdx == waterPhaseIdx) {
174  const auto& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
175  const auto& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
176 
177  // See: Ochs 2008
178  // \todo: proper citation
179  const LhsEval& rholH2O = H2O::liquidDensity(T, p);
180  const LhsEval& clH2O = rholH2O/H2O::molarMass();
181 
182  const auto& xwH2O = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(waterPhaseIdx, H2OIdx));
183  const auto& xwAir = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(waterPhaseIdx, airIdx));
184  const auto& xwNapl = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(waterPhaseIdx, NAPLIdx));
185  // this assumes each dissolved molecule displaces exactly one
186  // water molecule in the liquid
187  return clH2O*(H2O::molarMass()*xwH2O + Air::molarMass()*xwAir + NAPL::molarMass()*xwNapl);
188  }
189  else if (phaseIdx == naplPhaseIdx) {
190  // assume pure NAPL for the NAPL phase
191  const auto& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
192  return NAPL::liquidDensity(T, LhsEval(1e100));
193  }
194 
195  assert (phaseIdx == gasPhaseIdx);
196  const auto& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
197  const auto& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
198 
199  const LhsEval& pH2O = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(gasPhaseIdx, H2OIdx))*p;
200  const LhsEval& pAir = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(gasPhaseIdx, airIdx))*p;
201  const LhsEval& pNAPL = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(gasPhaseIdx, NAPLIdx))*p;
202  return
203  H2O::gasDensity(T, pH2O) +
204  Air::gasDensity(T, pAir) +
205  NAPL::gasDensity(T, pNAPL);
206  }
207 
209  template <class FluidState, class LhsEval = typename FluidState::Scalar>
210  static LhsEval viscosity(const FluidState &fluidState,
211  const ParameterCache &/*paramCache*/,
212  unsigned phaseIdx)
213  {
215 
216  const auto& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
217  const auto& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
218 
219  if (phaseIdx == waterPhaseIdx) {
220  // assume pure water viscosity
221  return H2O::liquidViscosity(T, p);
222  }
223  else if (phaseIdx == naplPhaseIdx) {
224  // assume pure NAPL viscosity
225  return NAPL::liquidViscosity(T, p);
226  }
227 
228  assert (phaseIdx == gasPhaseIdx);
229 
230  /* Wilke method. See:
231  *
232  * See: R. Reid, et al.: The Properties of Gases and Liquids,
233  * 4th edition, McGraw-Hill, 1987, 407-410
234  * 5th edition, McGraw-Hill, 20001, p. 9.21/22
235  *
236  * in this case, we use a simplified version in order to avoid
237  * computationally costly evaluation of sqrt and pow functions and
238  * divisions
239  * -- compare e.g. with Promo Class p. 32/33
240  */
241  const LhsEval mu[numComponents] = {
245  };
246  // molar masses
247  const Scalar M[numComponents] = {
248  H2O::molarMass(),
249  Air::molarMass(),
251  };
252 
253  const auto& xgAir = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(gasPhaseIdx, airIdx));
254  const auto& xgH2O = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(gasPhaseIdx, H2OIdx));
255  const auto& xgNapl = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(gasPhaseIdx, NAPLIdx));
256 
257  const LhsEval& xgAW = xgAir + xgH2O;
258  const LhsEval& muAW = (mu[airIdx]*xgAir + mu[H2OIdx]*xgH2O)/ xgAW;
259 
260  const LhsEval& MAW = (xgAir*Air::molarMass() + xgH2O*H2O::molarMass())/xgAW;
261 
262  Scalar phiCAW = 0.3; // simplification for this particular system
263  /* actually like this
264  * Scalar phiCAW = std::pow(1.+std::sqrt(mu[NAPLIdx]/muAW)*std::pow(MAW/M[NAPLIdx],0.25),2)
265  * / std::sqrt(8.*(1.+M[NAPLIdx]/MAW));
266  */
267  const LhsEval& phiAWC = phiCAW * muAW*M[NAPLIdx]/(mu[NAPLIdx]*MAW);
268 
269  return (xgAW*muAW)/(xgAW+xgNapl*phiAWC) + (xgNapl*mu[NAPLIdx])/(xgNapl + xgAW*phiCAW);
270  }
271 
273  template <class FluidState, class LhsEval = typename FluidState::Scalar>
274  static LhsEval diffusionCoefficient(const FluidState &fluidState,
275  const ParameterCache &/*paramCache*/,
276  unsigned phaseIdx,
277  unsigned compIdx)
278  {
280 
281  if (phaseIdx==gasPhaseIdx) {
282  const auto& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
283  const auto& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
284 
285  const LhsEval& diffAC = Opm::BinaryCoeff::Air_Xylene::gasDiffCoeff(T, p);
286  const LhsEval& diffWC = Opm::BinaryCoeff::H2O_Xylene::gasDiffCoeff(T, p);
287  const LhsEval& diffAW = Opm::BinaryCoeff::H2O_Air::gasDiffCoeff(T, p);
288 
289  const LhsEval& xga = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(gasPhaseIdx, airIdx));
290  const LhsEval& xgw = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(gasPhaseIdx, H2OIdx));
291  const LhsEval& xgc = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(gasPhaseIdx, NAPLIdx));
292 
293  if (compIdx==NAPLIdx) return (1.- xgw)/(xga/diffAW + xgc/diffWC);
294  else if (compIdx==H2OIdx) return (1.- xgc)/(xgw/diffWC + xga/diffAC);
295  else if (compIdx==airIdx) OPM_THROW(std::logic_error,
296  "Diffusivity of air in the gas phase "
297  "is constraint by sum of diffusive fluxes = 0 !\n");
298  } else if (phaseIdx==waterPhaseIdx){
299  Scalar diffACl = 1.e-9; // BinaryCoeff::Air_Xylene::liquidDiffCoeff(temperature, pressure);
300  Scalar diffWCl = 1.e-9; // BinaryCoeff::H2O_Xylene::liquidDiffCoeff(temperature, pressure);
301  Scalar diffAWl = 1.e-9; // BinaryCoeff::H2O_Air::liquidDiffCoeff(temperature, pressure);
302 
303  const LhsEval& xwa = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(waterPhaseIdx, airIdx));
304  const LhsEval& xww = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(waterPhaseIdx, H2OIdx));
305  const LhsEval& xwc = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(waterPhaseIdx, NAPLIdx));
306 
307  switch (compIdx) {
308  case NAPLIdx:
309  return (1.- xww)/(xwa/diffAWl + xwc/diffWCl);
310  case airIdx:
311  return (1.- xwc)/(xww/diffWCl + xwa/diffACl);
312  case H2OIdx:
313  OPM_THROW(std::logic_error,
314  "Diffusivity of water in the water phase "
315  "is constraint by sum of diffusive fluxes = 0 !\n");
316  };
317  } else if (phaseIdx==naplPhaseIdx) {
318 
319  OPM_THROW(std::logic_error,
320  "Diffusion coefficients of "
321  "substances in liquid phase are undefined!\n");
322  }
323  return 0;
324  }
325 
327  template <class FluidState, class LhsEval = typename FluidState::Scalar>
328  static LhsEval fugacityCoefficient(const FluidState &fluidState,
329  const ParameterCache &/*paramCache*/,
330  unsigned phaseIdx,
331  unsigned compIdx)
332  {
334 
335  assert(0 <= phaseIdx && phaseIdx < numPhases);
336  assert(0 <= compIdx && compIdx < numComponents);
337 
338  const auto& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
339  const auto& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
340 
341  if (phaseIdx == waterPhaseIdx) {
342  if (compIdx == H2OIdx)
343  return H2O::vaporPressure(T)/p;
344  else if (compIdx == airIdx)
346  else if (compIdx == NAPLIdx)
348  }
349 
350  // for the NAPL phase, we assume currently that nothing is
351  // dissolved. this means that the affinity of the NAPL
352  // component to the NAPL phase is much higher than for the
353  // other components, i.e. the fugacity cofficient is much
354  // smaller.
355  if (phaseIdx == naplPhaseIdx) {
356  const LhsEval& phiNapl = NAPL::vaporPressure(T)/p;
357  if (compIdx == NAPLIdx)
358  return phiNapl;
359  else if (compIdx == airIdx)
360  return 1e6*phiNapl;
361  else if (compIdx == H2OIdx)
362  return 1e6*phiNapl;
363  }
364 
365  // for the gas phase, assume an ideal gas when it comes to
366  // fugacity (-> fugacity == partial pressure)
367  assert(phaseIdx == gasPhaseIdx);
368  return 1.0;
369  }
370 
372  template <class FluidState, class LhsEval = typename FluidState::Scalar>
373  static LhsEval enthalpy(const FluidState &fluidState,
374  const ParameterCache &/*paramCache*/,
375  unsigned phaseIdx)
376  {
378 
379  const auto& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
380  const auto& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
381 
382  if (phaseIdx == waterPhaseIdx) {
383  return H2O::liquidEnthalpy(T, p);
384  }
385  else if (phaseIdx == naplPhaseIdx) {
386  return NAPL::liquidEnthalpy(T, p);
387  }
388  else if (phaseIdx == gasPhaseIdx) { // gas phase enthalpy depends strongly on composition
389  const LhsEval& hgc = NAPL::gasEnthalpy(T, p);
390  const LhsEval& hgw = H2O::gasEnthalpy(T, p);
391  const LhsEval& hga = Air::gasEnthalpy(T, p);
392 
393  LhsEval result = 0;
394  result += hgw * FsToolbox::template toLhs<LhsEval>(fluidState.massFraction(gasPhaseIdx, H2OIdx));
395  result += hga * FsToolbox::template toLhs<LhsEval>(fluidState.massFraction(gasPhaseIdx, airIdx));
396  result += hgc * FsToolbox::template toLhs<LhsEval>(fluidState.massFraction(gasPhaseIdx, NAPLIdx));
397 
398  return result;
399  }
400  OPM_THROW(std::logic_error, "Invalid phase index " << phaseIdx);
401  }
402 
403 private:
404  template <class LhsEval>
405  static LhsEval waterPhaseDensity_(const LhsEval& T,
406  const LhsEval& pw,
407  const LhsEval& xww,
408  const LhsEval& xwa,
409  const LhsEval& xwc)
410  {
411  const LhsEval& rholH2O = H2O::liquidDensity(T, pw);
412  const LhsEval& clH2O = rholH2O/H2O::molarMass();
413 
414  // this assumes each dissolved molecule displaces exactly one
415  // water molecule in the liquid
416  return clH2O*(xww*H2O::molarMass() + xwa*Air::molarMass() + xwc*NAPL::molarMass());
417  }
418 
419  template <class LhsEval>
420  static LhsEval gasPhaseDensity_(const LhsEval& T,
421  const LhsEval& pg,
422  const LhsEval& xgw,
423  const LhsEval& xga,
424  const LhsEval& xgc)
425  { return H2O::gasDensity(T, pg*xgw) + Air::gasDensity(T, pg*xga) + NAPL::gasDensity(T, pg*xgc); }
426 
427  template <class LhsEval>
428  static LhsEval NAPLPhaseDensity_(const LhsEval& T, const LhsEval& pn)
429  { return NAPL::liquidDensity(T, pn); }
430 };
431 } // namespace FluidSystems
432 } // namespace Opm
433 
434 #endif
static const int gasPhaseIdx
The index of the gas phase.
Definition: H2OAirXyleneFluidSystem.hpp:78
static bool gasIsIdeal()
Returns true iff the gas phase is assumed to be ideal.
Definition: Air.hpp:73
static bool liquidIsCompressible()
Returns true iff the liquid phase is assumed to be compressible.
Definition: Xylene.hpp:289
static const int NAPLIdx
The index of the NAPL component.
Definition: H2OAirXyleneFluidSystem.hpp:83
Material properties of pure water .
Definition: H2O.hpp:60
static LhsEval density(const FluidState &fluidState, const ParameterCache &, unsigned phaseIdx)
Calculate the density [kg/m^3] of a fluid phase.
Definition: H2OAirXyleneFluidSystem.hpp:167
static bool isCompressible(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: H2OAirXyleneFluidSystem.hpp:114
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 bool isLiquid(unsigned phaseIdx)
Return whether a phase is liquid.
Definition: H2OAirXyleneFluidSystem.hpp:92
static Evaluation liquidEnthalpy(const Evaluation &temperature, const Evaluation &pressure)
Specific enthalpy of the pure component in liquid.
Definition: Xylene.hpp:151
static Evaluation gasDiffCoeff(const Evaluation &temperature, const Evaluation &pressure)
Binary diffusion coefficent for molecular water and air.
Definition: H2O_Air.hpp:74
static void init()
Initialize the fluid system's static parameters.
Definition: H2OAirXyleneFluidSystem.hpp:88
static bool isIdealMixture(unsigned)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: H2OAirXyleneFluidSystem.hpp:103
A parameter cache which does nothing.
static const char * name()
A human readable name for the xylene.
Definition: Xylene.hpp:55
static const int naplPhaseIdx
The index of the NAPL phase.
Definition: H2OAirXyleneFluidSystem.hpp:76
static Evaluation gasEnthalpy(const Evaluation &temperature, const Evaluation &pressure)
Specific enthalpy of the pure component in gas.
Definition: Xylene.hpp:215
Opm::Xylene< Scalar > NAPL
The type of the xylene/napl component.
Definition: H2OAirXyleneFluidSystem.hpp:64
static const char * name()
A human readable name for the .
Definition: Air.hpp:61
static Scalar molarMass(unsigned compIdx)
Return the molar mass of a component in [kg/mol].
Definition: H2OAirXyleneFluidSystem.hpp:150
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
static Evaluation henry(const Evaluation &)
Henry coefficent for xylene in liquid water.
Definition: H2O_Xylene.hpp:49
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
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: H2OAirXyleneFluidSystem.hpp:373
A fluid system with water, gas and NAPL as phases and water, air and NAPL (contaminant) as components...
Definition: H2OAirXyleneFluidSystem.hpp:51
A simple class implementing the fluid properties of air.
Binary coefficients for water and xylene.
static LhsEval diffusionCoefficient(const FluidState &fluidState, const ParameterCache &, unsigned phaseIdx, unsigned compIdx)
Calculate the binary molecular diffusion coefficient for a component in a fluid phase [mol^2 * s / (k...
Definition: H2OAirXyleneFluidSystem.hpp:274
static const int H2OIdx
The index of the water component.
Definition: H2OAirXyleneFluidSystem.hpp:81
Material properties of pure water .
static Evaluation vaporPressure(const Evaluation &temperature)
The saturation vapor pressure in of pure xylene at a given temperature according to Antoine after Be...
Definition: Xylene.hpp:95
static const int airIdx
The index of the air pseudo-component.
Definition: H2OAirXyleneFluidSystem.hpp:85
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 Evaluation henry(const Evaluation &temperature)
Henry coefficent for air in liquid water.
Definition: H2O_Air.hpp:55
static const int numPhases
Number of fluid phases in the fluid system.
Definition: H2OAirXyleneFluidSystem.hpp:69
Component for Xylene.
static Evaluation liquidDensity(const Evaluation &temperature, const Evaluation &pressure)
The density of the liquid component at a given pressure in and temperature in . ...
Definition: Xylene.hpp:271
static Scalar molarMass()
The molar mass in of .
Definition: Air.hpp:81
static Evaluation simpleGasViscosity(const Evaluation &temperature, const Evaluation &)
Definition: Air.hpp:163
static const int numComponents
Number of chemical species in the fluid system.
Definition: H2OAirXyleneFluidSystem.hpp:71
Opm::Air< Scalar > Air
The type of the air component.
Definition: H2OAirXyleneFluidSystem.hpp:66
static Scalar molarMass()
The molar mass in of xylene.
Definition: Xylene.hpp:61
static Evaluation liquidViscosity(Evaluation temperature, const Evaluation &)
The dynamic liquid viscosity of the pure component.
Definition: Xylene.hpp:317
static Evaluation gasDensity(const Evaluation &temperature, const Evaluation &pressure)
The density in of the component at a given pressure in and temperature in .
Definition: Xylene.hpp:224
Opm::H2O< Scalar > H2O
The type of the water component.
Definition: H2OAirXyleneFluidSystem.hpp:62
static const char * componentName(unsigned compIdx)
Return the human readable name of a component.
Definition: H2OAirXyleneFluidSystem.hpp:139
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 Evaluation gasViscosity(Evaluation temperature, const Evaluation &)
The dynamic viscosity of the pure component at a given pressure in and temperature in ...
Definition: Xylene.hpp:296
static LhsEval viscosity(const FluidState &fluidState, const ParameterCache &, unsigned phaseIdx)
Calculate the dynamic viscosity of a fluid phase [Pa*s].
Definition: H2OAirXyleneFluidSystem.hpp:210
NullParameterCache ParameterCache
The type of the fluid system's parameter cache.
Definition: H2OAirXyleneFluidSystem.hpp:59
static Evaluation gasViscosity(const Evaluation &temperature, const Evaluation &pressure)
The dynamic viscosity of steam.
Definition: H2O.hpp:790
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: H2OAirXyleneFluidSystem.hpp:328
static const char * phaseName(unsigned phaseIdx)
Return the human readable name of a fluid phase.
Definition: H2OAirXyleneFluidSystem.hpp:128
static const char * name()
A human readable name for the water.
Definition: H2O.hpp:73
Component for Xylene.
Definition: Xylene.hpp:46
static bool gasIsIdeal()
Returns true iff the gas phase is assumed to be ideal.
Definition: Xylene.hpp:283
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
Binary coefficients for water and xylene.
Binary coefficients for water and nitrogen.
static const int waterPhaseIdx
The index of the water phase.
Definition: H2OAirXyleneFluidSystem.hpp:74
The base class for all fluid systems.
A simple class implementing the fluid properties of air.
Definition: Air.hpp:47
static Evaluation gasDiffCoeff(Evaluation temperature, Evaluation pressure)
Binary diffusion coefficent [m^2/s] for molecular water and xylene.
Definition: H2O_Xylene.hpp:64
static Evaluation gasDiffCoeff(Evaluation temperature, Evaluation pressure)
Binary diffusion coefficent [m^2/s] for air and xylene. method according to Wilke and Lee see Handboo...
Definition: Air_Xylene.hpp:55
static bool isIdealGas(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: H2OAirXyleneFluidSystem.hpp:99
static const Scalar molarMass()
The molar mass in of water.
Definition: H2O.hpp:79