H2OAirMesityleneFluidSystem.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  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_MESITYLENE_FLUID_SYSTEM_HPP
27 #define OPM_H2O_AIR_MESITYLENE_FLUID_SYSTEM_HPP
28 
29 #include "BaseFluidSystem.hpp"
30 #include "NullParameterCache.hpp"
31 
43 
44 #include <iostream>
45 
46 namespace Opm {
47 namespace FluidSystems {
48 
54 template <class Scalar>
56  : public BaseFluidSystem<Scalar, H2OAirMesitylene<Scalar> >
57 {
60 
61  typedef Opm::H2O<Scalar> IapwsH2O;
62  typedef Opm::TabulatedComponent<Scalar, IapwsH2O, /*alongVaporPressure=*/false> TabulatedH2O;
63 
64 public:
67 
70 
73 
75  //typedef SimpleH2O H2O;
76  typedef TabulatedH2O H2O;
77  //typedef IapwsH2O H2O;
78 
80  static const int numPhases = 3;
82  static const int numComponents = 3;
83 
85  static const int waterPhaseIdx = 0;
87  static const int naplPhaseIdx = 1;
89  static const int gasPhaseIdx = 2;
90 
92  static const int H2OIdx = 0;
94  static const int NAPLIdx = 1;
96  static const int airIdx = 2;
97 
99  static void init()
100  {
101  init(/*tempMin=*/273.15,
102  /*tempMax=*/623.15,
103  /*numTemp=*/100,
104  /*pMin=*/0.0,
105  /*pMax=*/20e6,
106  /*numP=*/200);
107  }
108 
120  static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp,
121  Scalar pressMin, Scalar pressMax, unsigned nPress)
122  {
123  if (H2O::isTabulated) {
124  TabulatedH2O::init(tempMin, tempMax, nTemp,
125  pressMin, pressMax, nPress);
126  }
127  }
128 
130  static bool isLiquid(unsigned phaseIdx)
131  {
132  //assert(0 <= phaseIdx && phaseIdx < numPhases);
133  return phaseIdx != gasPhaseIdx;
134  }
135 
137  static bool isIdealGas(unsigned phaseIdx)
138  { return phaseIdx == gasPhaseIdx && H2O::gasIsIdeal() && Air::gasIsIdeal() && NAPL::gasIsIdeal(); }
139 
141  static bool isCompressible(unsigned phaseIdx)
142  {
143  //assert(0 <= phaseIdx && phaseIdx < numPhases);
144  // gases are always compressible
145  return (phaseIdx == gasPhaseIdx)
146  ? true
147  : (phaseIdx == waterPhaseIdx)
150  }
151 
153  static bool isIdealMixture(unsigned /*phaseIdx*/)
154  {
155  //assert(0 <= phaseIdx && phaseIdx < numPhases);
156  // we assume Henry's and Rault's laws for the water phase and
157  // and no interaction between gas molecules of different
158  // components, so all phases are ideal mixtures!
159  return true;
160  }
161 
163  static const char *phaseName(unsigned phaseIdx)
164  {
165  switch (phaseIdx) {
166  case waterPhaseIdx: return "water";
167  case naplPhaseIdx: return "napl";
168  case gasPhaseIdx: return "gas";
169  };
170  OPM_THROW(std::logic_error, "Invalid phase index " << phaseIdx);
171  }
172 
174  static const char *componentName(unsigned compIdx)
175  {
176  switch (compIdx) {
177  case H2OIdx: return H2O::name();
178  case airIdx: return Air::name();
179  case NAPLIdx: return NAPL::name();
180  };
181  OPM_THROW(std::logic_error, "Invalid component index " << compIdx);
182  }
183 
185  static Scalar molarMass(unsigned compIdx)
186  {
187  return
188  (compIdx == H2OIdx)
189  ? H2O::molarMass()
190  : (compIdx == airIdx)
191  ? Air::molarMass()
192  : (compIdx == NAPLIdx)
193  ? NAPL::molarMass()
194  : -1e10;
195  //OPM_THROW(std::logic_error, "Invalid component index " << compIdx);
196  }
197 
199  template <class FluidState, class LhsEval = typename FluidState::Scalar>
200  static LhsEval density(const FluidState &fluidState,
201  const ParameterCache &/*paramCache*/,
202  unsigned phaseIdx)
203  {
205 
206  const LhsEval& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
207 
208  if (phaseIdx == waterPhaseIdx) {
209  // See: Ochs 2008
210  const LhsEval& p =
212  ? FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx))
213  : 1e100;
214 
215  const LhsEval& rholH2O = H2O::liquidDensity(T, p);
216  const LhsEval& clH2O = rholH2O/H2O::molarMass();
217 
218  // this assumes each dissolved molecule displaces exactly one
219  // water molecule in the liquid
220  return
221  clH2O*(H2O::molarMass()*FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(waterPhaseIdx, H2OIdx)) +
222  Air::molarMass()*FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(waterPhaseIdx, airIdx)) +
223  NAPL::molarMass()*FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(waterPhaseIdx, NAPLIdx)));
224  }
225  else if (phaseIdx == naplPhaseIdx) {
226  // assume pure NAPL for the NAPL phase
227  const LhsEval& p =
229  ? FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx))
230  : 1e100;
231  return NAPL::liquidDensity(T, p);
232  }
233 
234  assert (phaseIdx == gasPhaseIdx);
235  const LhsEval& pg = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(gasPhaseIdx));
236  const LhsEval& pH2O = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(gasPhaseIdx, H2OIdx))*pg;
237  const LhsEval& pAir = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(gasPhaseIdx, airIdx))*pg;
238  const LhsEval& pNAPL = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(gasPhaseIdx, NAPLIdx))*pg;
239  return
240  H2O::gasDensity(T, pH2O) +
241  Air::gasDensity(T, pAir) +
242  NAPL::gasDensity(T, pNAPL);
243  }
244 
246  template <class FluidState, class LhsEval = typename FluidState::Scalar>
247  static LhsEval viscosity(const FluidState &fluidState,
248  const ParameterCache &/*paramCache*/,
249  unsigned phaseIdx)
250  {
252 
253  const LhsEval& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
254  const LhsEval& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
255 
256  if (phaseIdx == waterPhaseIdx) {
257  // assume pure water viscosity
258 
259  return H2O::liquidViscosity(T,
260  p);
261  }
262  else if (phaseIdx == naplPhaseIdx) {
263  // assume pure NAPL viscosity
264  return NAPL::liquidViscosity(T, p);
265  }
266 
267  assert (phaseIdx == gasPhaseIdx);
268 
269  /* Wilke method. See:
270  *
271  * See: R. Reid, et al.: The Properties of Gases and Liquids,
272  * 4th edition, McGraw-Hill, 1987, 407-410
273  * 5th edition, McGraw-Hill, 20001, p. 9.21/22
274  *
275  * in this case, we use a simplified version in order to avoid
276  * computationally costly evaluation of sqrt and pow functions and
277  * divisions
278  * -- compare e.g. with Promo Class p. 32/33
279  */
280  const LhsEval mu[numComponents] = {
282  Air::gasViscosity(T, p),
284  };
285  // molar masses
286  const Scalar M[numComponents] = {
287  H2O::molarMass(),
288  Air::molarMass(),
290  };
291 
292  const LhsEval& xgAir = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(gasPhaseIdx, airIdx));
293  const LhsEval& xgH2O = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(gasPhaseIdx, H2OIdx));
294  const LhsEval& xgNapl = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(gasPhaseIdx, NAPLIdx));
295  const LhsEval& xgAW = xgAir + xgH2O;
296  const LhsEval& muAW = (mu[airIdx]*xgAir + mu[H2OIdx]*xgH2O)/xgAW;
297  const LhsEval& MAW = (xgAir*Air::molarMass() + xgH2O*H2O::molarMass())/xgAW;
298 
299  Scalar phiCAW = 0.3; // simplification for this particular system
300  /* actually like this
301  * Scalar phiCAW = std::pow(1.+std::sqrt(mu[NAPLIdx]/muAW)*std::pow(MAW/M[NAPLIdx],0.25),2)
302  * / std::sqrt(8.*(1.+M[NAPLIdx]/MAW));
303  */
304  const LhsEval& phiAWC = phiCAW * muAW*M[NAPLIdx]/(mu[NAPLIdx]*MAW);
305 
306  return (xgAW*muAW)/(xgAW + xgNapl*phiAWC) + (xgNapl*mu[NAPLIdx])/(xgNapl + xgAW*phiCAW);
307  }
308 
310  template <class FluidState, class LhsEval = typename FluidState::Scalar>
311  static LhsEval diffusionCoefficient(const FluidState &/*fluidState*/,
312  const ParameterCache &/*paramCache*/,
313  unsigned /*phaseIdx*/,
314  unsigned /*compIdx*/)
315  {
316  return 0;
317 #if 0
319 
320  const LhsEval& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
321  const LhsEval& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
322  LhsEval diffCont;
323 
324  if (phaseIdx==gasPhaseIdx) {
325  const LhsEval& diffAC = Opm::BinaryCoeff::Air_Mesitylene::gasDiffCoeff(T, p);
326  const LhsEval& diffWC = Opm::BinaryCoeff::H2O_Mesitylene::gasDiffCoeff(T, p);
327  const LhsEval& diffAW = Opm::BinaryCoeff::H2O_Air::gasDiffCoeff(T, p);
328 
329  const LhsEval& xga = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(gasPhaseIdx, airIdx));
330  const LhsEval& xgw = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(gasPhaseIdx, H2OIdx));
331  const LhsEval& xgc = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(gasPhaseIdx, NAPLIdx));
332 
333  if (compIdx==NAPLIdx) return (1 - xgw)/(xga/diffAW + xgc/diffWC);
334  else if (compIdx==H2OIdx) return (1 - xgc)/(xgw/diffWC + xga/diffAC);
335  else if (compIdx==airIdx) OPM_THROW(std::logic_error,
336  "Diffusivity of air in the gas phase "
337  "is constraint by sum of diffusive fluxes = 0 !\n");
338  }
339  else if (phaseIdx==waterPhaseIdx){
340  const LhsEval& diffACl = 1.e-9; // BinaryCoeff::Air_Mesitylene::liquidDiffCoeff(temperature, pressure);
341  const LhsEval& diffWCl = 1.e-9; // BinaryCoeff::H2O_Mesitylene::liquidDiffCoeff(temperature, pressure);
342  const LhsEval& diffAWl = 1.e-9; // BinaryCoeff::H2O_Air::liquidDiffCoeff(temperature, pressure);
343 
344  const LhsEval& xwa = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(waterPhaseIdx, airIdx));
345  const LhsEval& xww = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(waterPhaseIdx, H2OIdx));
346  const LhsEval& xwc = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(waterPhaseIdx, NAPLIdx));
347 
348  switch (compIdx) {
349  case NAPLIdx:
350  diffCont = (1.- xww)/(xwa/diffAWl + xwc/diffWCl);
351  return diffCont;
352  case airIdx:
353  diffCont = (1.- xwc)/(xww/diffWCl + xwa/diffACl);
354  return diffCont;
355  case H2OIdx:
356  OPM_THROW(std::logic_error,
357  "Diffusivity of water in the water phase "
358  "is constraint by sum of diffusive fluxes = 0 !\n");
359  };
360  }
361  else if (phaseIdx==naplPhaseIdx) {
362  OPM_THROW(std::logic_error,
363  "Diffusion coefficients of "
364  "substances in liquid phase are undefined!\n");
365  }
366  return 0;
367 #endif
368  }
369 
371  template <class FluidState, class LhsEval = typename FluidState::Scalar>
372  static LhsEval fugacityCoefficient(const FluidState &fluidState,
373  const ParameterCache &/*paramCache*/,
374  unsigned phaseIdx,
375  unsigned compIdx)
376  {
378 
379  assert(0 <= phaseIdx && phaseIdx < numPhases);
380  assert(0 <= compIdx && compIdx < numComponents);
381 
382  const LhsEval& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
383  const LhsEval& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
386 
387  if (phaseIdx == waterPhaseIdx) {
388  if (compIdx == H2OIdx)
389  return H2O::vaporPressure(T)/p;
390  else if (compIdx == airIdx)
392  else if (compIdx == NAPLIdx)
394  assert(false);
395  }
396  // for the NAPL phase, we assume currently that nothing is
397  // dissolved. this means that the affinity of the NAPL
398  // component to the NAPL phase is much higher than for the
399  // other components, i.e. the fugacity cofficient is much
400  // smaller.
401  else if (phaseIdx == naplPhaseIdx) {
402  const LhsEval& phiNapl = NAPL::vaporPressure(T)/p;
403  if (compIdx == NAPLIdx)
404  return phiNapl;
405  else if (compIdx == airIdx)
406  return 1e6*phiNapl;
407  else if (compIdx == H2OIdx)
408  return 1e6*phiNapl;
409  assert(false);
410  }
411 
412  // for the gas phase, assume an ideal gas when it comes to
413  // fugacity (-> fugacity == partial pressure)
414  assert(phaseIdx == gasPhaseIdx);
415  return 1.0;
416  }
417 
418 
420  template <class FluidState, class LhsEval = typename FluidState::Scalar>
421  static LhsEval enthalpy(const FluidState &fluidState,
422  const ParameterCache &/*paramCache*/,
423  unsigned phaseIdx)
424  {
426 
427  const LhsEval& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
428  const LhsEval& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
429 
430  if (phaseIdx == waterPhaseIdx) {
431  return H2O::liquidEnthalpy(T, p);
432  }
433  else if (phaseIdx == naplPhaseIdx) {
434  return NAPL::liquidEnthalpy(T, p);
435  }
436  else if (phaseIdx == gasPhaseIdx) {
437  // gas phase enthalpy depends strongly on composition
438  LhsEval result = 0;
439  result += H2O::gasEnthalpy(T, p) * FsToolbox::template toLhs<LhsEval>(fluidState.massFraction(gasPhaseIdx, H2OIdx));
440  result += NAPL::gasEnthalpy(T, p) * FsToolbox::template toLhs<LhsEval>(fluidState.massFraction(gasPhaseIdx, airIdx));
441  result += Air::gasEnthalpy(T, p) * FsToolbox::template toLhs<LhsEval>(fluidState.massFraction(gasPhaseIdx, NAPLIdx));
442 
443  return result;
444  }
445  OPM_THROW(std::logic_error, "Invalid phase index " << phaseIdx);
446  }
447 
449  template <class FluidState, class LhsEval = typename FluidState::Scalar>
450  static LhsEval thermalConductivity(const FluidState &fluidState,
451  const ParameterCache &/*paramCache*/,
452  unsigned phaseIdx)
453  {
455 
456  assert(0 <= phaseIdx && phaseIdx < numPhases);
457 
458  if (phaseIdx == waterPhaseIdx){ // water phase
459  const LhsEval& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
460  const LhsEval& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
461 
462  return H2O::liquidThermalConductivity(T, p);
463  }
464  else if (phaseIdx == gasPhaseIdx) { // gas phase
465  const LhsEval& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
466  const LhsEval& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
467 
468  return Air::gasThermalConductivity(T, p);
469  }
470 
471  assert(phaseIdx == naplPhaseIdx);
472 
473  // Taken from:
474  //
475  // D. K. H. Briggs: "Thermal Conductivity of Liquids",
476  // Ind. Eng. Chem., 1957, 49 (3), pp 418–421
477  //
478  // Convertion to SI units:
479  // 344e-6 cal/(s cm K) = 0.0143964 J/(s m K)
480  return 0.0143964;
481  }
482 };
483 } // namespace FluidSystems
484 } // namespace Opm
485 
486 #endif
static Scalar molarMass(unsigned compIdx)
Return the molar mass of a component in [kg/mol].
Definition: H2OAirMesityleneFluidSystem.hpp:185
static const int waterPhaseIdx
The index of the water phase.
Definition: H2OAirMesityleneFluidSystem.hpp:85
Component for Mesitylene.
Definition: Mesitylene.hpp:43
static const int naplPhaseIdx
The index of the NAPL phase.
Definition: H2OAirMesityleneFluidSystem.hpp:87
static const char * name()
A human readable name for the component.
Definition: TabulatedComponent.hpp:216
static bool gasIsIdeal()
Returns true iff the gas phase is assumed to be ideal.
Definition: Air.hpp:73
static const int gasPhaseIdx
The index of the gas phase.
Definition: H2OAirMesityleneFluidSystem.hpp:89
static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp, Scalar pressMin, Scalar pressMax, unsigned nPress)
Initialize the tables.
Definition: TabulatedComponent.hpp:73
static LhsEval thermalConductivity(const FluidState &fluidState, const ParameterCache &, unsigned phaseIdx)
Thermal conductivity of a fluid phase [W/(m K)].
Definition: H2OAirMesityleneFluidSystem.hpp:450
Material properties of pure water .
Definition: H2O.hpp:60
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 bool isLiquid(unsigned phaseIdx)
Return whether a phase is liquid.
Definition: H2OAirMesityleneFluidSystem.hpp:130
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: H2OAirMesityleneFluidSystem.hpp:372
TabulatedH2O H2O
The type of the water component.
Definition: H2OAirMesityleneFluidSystem.hpp:76
Definition: MathToolbox.hpp:39
static Evaluation gasDiffCoeff(Evaluation temperature, Evaluation pressure)
Binary diffusion coefficent [m^2/s] for molecular water and mesitylene.
Definition: H2O_Mesitylene.hpp:65
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
static Evaluation gasDiffCoeff(const Evaluation &temperature, const Evaluation &pressure)
Binary diffusion coefficent for molecular water and air.
Definition: H2O_Air.hpp:74
static bool isCompressible(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: H2OAirMesityleneFluidSystem.hpp:141
A parameter cache which does nothing.
static LhsEval viscosity(const FluidState &fluidState, const ParameterCache &, unsigned phaseIdx)
Calculate the dynamic viscosity of a fluid phase [Pa*s].
Definition: H2OAirMesityleneFluidSystem.hpp:247
Binary coefficients for water and nitrogen.
static const char * name()
A human readable name for the .
Definition: Air.hpp:61
static Evaluation henry(const Evaluation &temperature)
Henry coefficent for molecular nitrogen in liquid water.
Definition: H2O_N2.hpp:52
static const int NAPLIdx
The index of the NAPL component.
Definition: H2OAirMesityleneFluidSystem.hpp:94
static void init()
Initialize the fluid system's static parameters.
Definition: H2OAirMesityleneFluidSystem.hpp:99
static Evaluation gasEnthalpy(const Evaluation &temperature, const Evaluation &pressure)
Specific enthalpy of mesitylene vapor .
Definition: Mesitylene.hpp:180
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 .
NullParameterCache ParameterCache
The type of the fluid system's parameter cache.
Definition: H2OAirMesityleneFluidSystem.hpp:66
static Evaluation gasDiffCoeff(Evaluation temperature, Evaluation pressure)
Binary diffusion coefficent [m^2/s] for air and mesitylene. I used the method according to Wilke and ...
Definition: Air_Mesitylene.hpp:57
static Scalar molarMass()
The molar mass in of the component.
Definition: TabulatedComponent.hpp:222
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: H2OAirMesityleneFluidSystem.hpp:120
static bool liquidIsCompressible()
Returns true iff the liquid phase is assumed to be compressible.
Definition: Mesitylene.hpp:220
A simple class implementing the fluid properties of air.
static const char * name()
A human readable name for the mesitylene.
Definition: Mesitylene.hpp:51
Material properties of pure water .
static Evaluation gasViscosity(const Evaluation &temperature, const Evaluation &)
The dynamic viscosity of at a given pressure and temperature.
Definition: Air.hpp:138
static const int numPhases
Number of fluid phases in the fluid system.
Definition: H2OAirMesityleneFluidSystem.hpp:80
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 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: H2OAirMesityleneFluidSystem.hpp:421
static Evaluation gasDensity(const Evaluation &temperature, const Evaluation &pressure)
The density of pure mesitylene vapor at a given pressure and temperature .
Definition: Mesitylene.hpp:192
Binary coefficients for water and mesitylene.
static bool isIdealGas(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: H2OAirMesityleneFluidSystem.hpp:137
A simple version of pure water.
Opm::Air< Scalar > Air
The type of the air component.
Definition: H2OAirMesityleneFluidSystem.hpp:72
static const int airIdx
The index of the air pseudo-component.
Definition: H2OAirMesityleneFluidSystem.hpp:96
static Evaluation gasViscosity(Evaluation temperature, const Evaluation &, bool=true)
The dynamic viscosity of mesitylene vapor.
Definition: Mesitylene.hpp:231
static Scalar molarMass()
The molar mass in of .
Definition: Air.hpp:81
static bool gasIsIdeal()
Returns true iff the gas phase is assumed to be ideal.
Definition: TabulatedComponent.hpp:417
Binary coefficients for water and mesitylene.
Component for Mesitylene.
A fluid system with water, gas and NAPL as phases and water, air and mesitylene (DNAPL) as components...
Definition: H2OAirMesityleneFluidSystem.hpp:55
A generic class which tabulates all thermodynamic properties of a given component.
Definition: TabulatedComponent.hpp:56
static Evaluation liquidEnthalpy(const Evaluation &temperature, const Evaluation &pressure)
Specific enthalpy of the liquid .
Definition: TabulatedComponent.hpp:292
static Evaluation gasThermalConductivity(const Evaluation &, const Evaluation &)
Specific heat conductivity of steam .
Definition: Air.hpp:223
static const bool isTabulated
Definition: TabulatedComponent.hpp:61
static const char * componentName(unsigned compIdx)
Return the human readable name of a component.
Definition: H2OAirMesityleneFluidSystem.hpp:174
Relations valid for an ideal gas.
A generic class which tabulates all thermodynamic properties of a given component.
static Evaluation henry(const Evaluation &)
Henry coefficent for mesitylene in liquid water.
Definition: H2O_Mesitylene.hpp:50
static bool gasIsIdeal()
Returns true iff the gas phase is assumed to be ideal.
Definition: Mesitylene.hpp:214
static bool liquidIsCompressible()
Returns true iff the liquid phase is assumed to be compressible.
Definition: TabulatedComponent.hpp:411
static const int numComponents
Number of chemical species in the fluid system.
Definition: H2OAirMesityleneFluidSystem.hpp:82
static Evaluation liquidEnthalpy(const Evaluation &temperature, const Evaluation &pressure)
Specific enthalpy of liquid mesitylene .
Definition: Mesitylene.hpp:119
static Evaluation gasEnthalpy(const Evaluation &temperature, const Evaluation &pressure)
Specific enthalpy of the gas .
Definition: TabulatedComponent.hpp:273
static Evaluation liquidDensity(const Evaluation &temperature, const Evaluation &)
The density of pure mesitylene at a given pressure and temperature .
Definition: Mesitylene.hpp:202
static LhsEval diffusionCoefficient(const FluidState &, const ParameterCache &, unsigned, unsigned)
Calculate the binary molecular diffusion coefficient for a component in a fluid phase [mol^2 * s / (k...
Definition: H2OAirMesityleneFluidSystem.hpp:311
static LhsEval density(const FluidState &fluidState, const ParameterCache &, unsigned phaseIdx)
Calculate the density [kg/m^3] of a fluid phase.
Definition: H2OAirMesityleneFluidSystem.hpp:200
static Evaluation gasDensity(const Evaluation &temperature, const Evaluation &pressure)
The density of gas at a given pressure and temperature .
Definition: TabulatedComponent.hpp:429
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
Binary coefficients for water and nitrogen.
static const int H2OIdx
The index of the water component.
Definition: H2OAirMesityleneFluidSystem.hpp:92
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
static Evaluation liquidViscosity(Evaluation temperature, const Evaluation &)
The dynamic viscosity of pure mesitylene.
Definition: Mesitylene.hpp:259
A simple class implementing the fluid properties of air.
Definition: Air.hpp:47
static Scalar molarMass()
The molar mass in of mesitylene.
Definition: Mesitylene.hpp:57
static const char * phaseName(unsigned phaseIdx)
Return the human readable name of a fluid phase.
Definition: H2OAirMesityleneFluidSystem.hpp:163
Opm::Mesitylene< Scalar > NAPL
The type of the mesithylene/napl component.
Definition: H2OAirMesityleneFluidSystem.hpp:69
static bool isIdealMixture(unsigned)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: H2OAirMesityleneFluidSystem.hpp:153
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
static Evaluation vaporPressure(const Evaluation &temperature)
The saturation vapor pressure in of pure mesitylene at a given temperature according to Antoine afte...
Definition: Mesitylene.hpp:98