H2ON2LiquidPhaseFluidSystem.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) 2012 by Bernd Flemisch
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_N2_LIQUIDPHASE_FLUID_SYSTEM_HPP
27 #define OPM_H2O_N2_LIQUIDPHASE_FLUID_SYSTEM_HPP
28 
29 #include "BaseFluidSystem.hpp"
30 #include "NullParameterCache.hpp"
31 
39 
40 #include <opm/common/Exceptions.hpp>
41 #include <opm/common/ErrorMacros.hpp>
42 
43 #include <iostream>
44 #include <cassert>
45 
46 namespace Opm {
47 namespace FluidSystems {
48 
55 template <class Scalar, bool useComplexRelations = true>
57  : public BaseFluidSystem<Scalar, H2ON2LiquidPhase<Scalar, useComplexRelations> >
58 {
61 
62  // convenience typedefs
63  typedef Opm::H2O<Scalar> IapwsH2O;
65  typedef Opm::N2<Scalar> SimpleN2;
66 
67 public:
70 
71  /****************************************
72  * Fluid phase related static parameters
73  ****************************************/
74 
76  static const int numPhases = 1;
77 
79  static const int liquidPhaseIdx = 0;
80 
82  static const char *phaseName(OPM_OPTIM_UNUSED unsigned phaseIdx)
83  {
84  assert(phaseIdx == liquidPhaseIdx);
85 
86  return "liquid";
87  }
88 
90  static bool isLiquid(unsigned /*phaseIdx*/)
91  {
92  //assert(phaseIdx == liquidPhaseIdx);
93  return true; //only water phase present
94  }
95 
97  static bool isCompressible(unsigned /*phaseIdx*/)
98  {
99  //assert(0 <= phaseIdx && phaseIdx < numPhases);
100  // the water component decides for the liquid phase...
101  return H2O::liquidIsCompressible();
102  }
103 
105  static bool isIdealGas(unsigned /*phaseIdx*/)
106  {
107  //assert(0 <= phaseIdx && phaseIdx < numPhases);
108  return false; // not a gas (only liquid phase present)
109  }
110 
112  static bool isIdealMixture(unsigned /*phaseIdx*/)
113  {
114  //assert(0 <= phaseIdx && phaseIdx < numPhases);
115  // we assume Henry's and Rault's laws for the water phase and
116  // and no interaction between gas molecules of different
117  // components, so all phases are ideal mixtures!
118  return true;
119  }
120 
121  /****************************************
122  * Component related static parameters
123  ****************************************/
124 
126  static const int numComponents = 2;
127 
129  static const int H2OIdx = 0;
131  static const int N2Idx = 1;
132 
134  typedef TabulatedH2O H2O;
135  //typedef SimpleH2O H2O;
136  //typedef IapwsH2O H2O;
137 
139  typedef SimpleN2 N2;
140 
142  static const char *componentName(unsigned compIdx)
143  {
144  static const char *name[] = {
145  H2O::name(),
146  N2::name()
147  };
148 
149  assert(0 <= compIdx && compIdx < numComponents);
150  return name[compIdx];
151  }
152 
154  static Scalar molarMass(unsigned compIdx)
155  {
156  //assert(0 <= compIdx && compIdx < numComponents);
157  return (compIdx == H2OIdx)
158  ? H2O::molarMass()
159  : (compIdx == N2Idx)
160  ? N2::molarMass()
161  : 1e100;
162  }
163 
169  static Scalar criticalTemperature(unsigned compIdx)
170  {
171  //assert(0 <= compIdx && compIdx < numComponents);
172  return (compIdx == H2OIdx)
174  : (compIdx == N2Idx)
176  : 1e100;
177  }
178 
184  static Scalar criticalPressure(unsigned compIdx)
185  {
186  //assert(0 <= compIdx && compIdx < numComponents);
187  return (compIdx == H2OIdx)
189  : (compIdx == N2Idx)
191  : 1e100;
192  }
193 
199  static Scalar acentricFactor(unsigned compIdx)
200  {
201  //assert(0 <= compIdx && compIdx < numComponents);
202  return (compIdx == H2OIdx)
204  : (compIdx == N2Idx)
205  ? N2::acentricFactor()
206  : 1e100;
207  }
208 
209  /****************************************
210  * thermodynamic relations
211  ****************************************/
212 
219  static void init()
220  {
221  init(/*tempMin=*/273.15,
222  /*tempMax=*/623.15,
223  /*numTemp=*/100,
224  /*pMin=*/0.0,
225  /*pMax=*/20e6,
226  /*numP=*/200);
227  }
228 
240  static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp,
241  Scalar pressMin, Scalar pressMax, unsigned nPress)
242  {
243  if (H2O::isTabulated) {
244  TabulatedH2O::init(tempMin, tempMax, nTemp,
245  pressMin, pressMax, nPress);
246  }
247  }
248 
250  template <class FluidState, class LhsEval = typename FluidState::Scalar>
251  static LhsEval density(const FluidState &fluidState,
252  const ParameterCache &/*paramCache*/,
253  unsigned phaseIdx)
254  {
256 
257  assert(0 <= phaseIdx && phaseIdx < numPhases);
258 
259  const auto& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
260  const auto& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
261 
262  LhsEval sumMoleFrac = 0;
263  for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
264  sumMoleFrac += FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(phaseIdx, compIdx));
265 
266  assert(phaseIdx == liquidPhaseIdx);
267 
268  if (!useComplexRelations)
269  // assume pure water
270  return H2O::liquidDensity(T, p);
271  else
272  {
273  // See: Ochs 2008
274  const LhsEval& rholH2O = H2O::liquidDensity(T, p);
275  const LhsEval& clH2O = rholH2O/H2O::molarMass();
276 
277  const auto& xlH2O = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(liquidPhaseIdx, H2OIdx));
278  const auto& xlN2 = FsToolbox::template toLhs<LhsEval>(fluidState.moleFraction(liquidPhaseIdx, N2Idx));
279 
280  // this assumes each nitrogen molecule displaces exactly one
281  // water molecule in the liquid
282  return clH2O*(H2O::molarMass()*xlH2O + N2::molarMass()*xlN2)/sumMoleFrac;
283  }
284  }
285 
287  template <class FluidState, class LhsEval = typename FluidState::Scalar>
288  static LhsEval viscosity(const FluidState &fluidState,
289  const ParameterCache &/*paramCache*/,
290  unsigned phaseIdx)
291  {
293 
294  assert(phaseIdx == liquidPhaseIdx);
295 
296  const auto& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
297  const auto& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
298 
299  // assume pure water for the liquid phase
300  return H2O::liquidViscosity(T, p);
301  }
302 
304  template <class FluidState, class LhsEval = typename FluidState::Scalar>
305  static LhsEval fugacityCoefficient(const FluidState &fluidState,
306  const ParameterCache &/*paramCache*/,
307  unsigned phaseIdx,
308  unsigned compIdx)
309  {
311 
312  assert(phaseIdx == liquidPhaseIdx);
313  assert(0 <= compIdx && compIdx < numComponents);
314 
315  const auto& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
316  const auto& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
317 
318  if (compIdx == H2OIdx)
319  return H2O::vaporPressure(T)/p;
321  }
322 
324  template <class FluidState, class LhsEval = typename FluidState::Scalar>
325  static LhsEval diffusionCoefficient(const FluidState &fluidState,
326  const ParameterCache &/*paramCache*/,
327  unsigned phaseIdx,
328  unsigned /*compIdx*/)
329 
330  {
332 
333  assert(phaseIdx == liquidPhaseIdx);
334 
335  const auto& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
336  const auto& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
337 
339  }
340 
342  template <class FluidState, class LhsEval = typename FluidState::Scalar>
343  static LhsEval enthalpy(const FluidState &fluidState,
344  const ParameterCache &/*paramCache*/,
345  unsigned phaseIdx)
346  {
348 
349  assert (phaseIdx == liquidPhaseIdx);
350 
351  const auto& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
352  const auto& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
355 
356  // TODO: way to deal with the solutes???
357  return H2O::liquidEnthalpy(T, p);
358  }
359 
361  template <class FluidState, class LhsEval = typename FluidState::Scalar>
362  static LhsEval thermalConductivity(const FluidState &fluidState,
363  const ParameterCache &/*paramCache*/,
364  const unsigned phaseIdx)
365  {
367 
368  assert(phaseIdx == liquidPhaseIdx);
369 
370  if(useComplexRelations){
371  const auto& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
372  const auto& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
373  return H2O::liquidThermalConductivity(T, p);
374  }
375  else
376  return 0.578078; // conductivity of water[W / (m K ) ] IAPWS evaluated at p=.1 MPa, T=8C
377  }
378 
380  template <class FluidState, class LhsEval = typename FluidState::Scalar>
381  static LhsEval heatCapacity(const FluidState &fluidState,
382  const ParameterCache &/*paramCache*/,
383  unsigned phaseIdx)
384  {
386 
387  assert (phaseIdx == liquidPhaseIdx);
388 
389  const auto& T = FsToolbox::template toLhs<LhsEval>(fluidState.temperature(phaseIdx));
390  const auto& p = FsToolbox::template toLhs<LhsEval>(fluidState.pressure(phaseIdx));
391 
392  return H2O::liquidHeatCapacity(T, p);
393  }
394 };
395 
396 } // namespace FluidSystems
397 
398 } // namespace Opm
399 
400 #endif
static const char * name()
A human readable name for the component.
Definition: TabulatedComponent.hpp:216
TabulatedH2O H2O
The type of the component for pure water.
Definition: H2ON2LiquidPhaseFluidSystem.hpp:134
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
static bool isLiquid(unsigned)
Return whether a phase is liquid.
Definition: H2ON2LiquidPhaseFluidSystem.hpp:90
static bool isIdealMixture(unsigned)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: H2ON2LiquidPhaseFluidSystem.hpp:112
Material properties of pure water .
Definition: H2O.hpp:60
static bool isIdealGas(unsigned)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: H2ON2LiquidPhaseFluidSystem.hpp:105
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 Scalar criticalPressure()
Returns the critical pressure of molecular nitrogen.
Definition: N2.hpp:75
Definition: MathToolbox.hpp:39
static LhsEval viscosity(const FluidState &fluidState, const ParameterCache &, unsigned phaseIdx)
Calculate the dynamic viscosity of a fluid phase [Pa*s].
Definition: H2ON2LiquidPhaseFluidSystem.hpp:288
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 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: H2ON2LiquidPhaseFluidSystem.hpp:240
A parameter cache which does nothing.
Binary coefficients for water and nitrogen.
static const int numComponents
Number of chemical species in the fluid system.
Definition: H2ON2LiquidPhaseFluidSystem.hpp:126
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.
static Evaluation liquidHeatCapacity(const Evaluation &temperature, const Evaluation &pressure)
Specific isobaric heat capacity of the liquid .
Definition: TabulatedComponent.hpp:330
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 const char * phaseName(OPM_OPTIM_UNUSED unsigned phaseIdx)
Return the human readable name of a fluid phase.
Definition: H2ON2LiquidPhaseFluidSystem.hpp:82
static Scalar molarMass()
The molar mass in of the component.
Definition: TabulatedComponent.hpp:222
A liquid-phase-only fluid system with water and nitrogen as components.
Definition: H2ON2LiquidPhaseFluidSystem.hpp:56
static const int N2Idx
The index of the component for molecular nitrogen.
Definition: H2ON2LiquidPhaseFluidSystem.hpp:131
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: H2ON2LiquidPhaseFluidSystem.hpp:325
static LhsEval thermalConductivity(const FluidState &fluidState, const ParameterCache &, const unsigned phaseIdx)
Thermal conductivity of a fluid phase [W/(m K)].
Definition: H2ON2LiquidPhaseFluidSystem.hpp:362
Material properties of pure water .
static Evaluation liquidDiffCoeff(const Evaluation &temperature, const Evaluation &)
Diffusion coefficent for molecular nitrogen in liquid water.
Definition: H2O_N2.hpp:102
static Scalar criticalPressure()
Returns the critical pressure in of the component.
Definition: TabulatedComponent.hpp:234
A simple version of pure water.
static Scalar criticalPressure(unsigned compIdx)
Critical pressure of a component [Pa].
Definition: H2ON2LiquidPhaseFluidSystem.hpp:184
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: H2ON2LiquidPhaseFluidSystem.hpp:343
static Scalar criticalTemperature()
Returns the critical temperature in of the component.
Definition: TabulatedComponent.hpp:228
A generic class which tabulates all thermodynamic properties of a given component.
Definition: TabulatedComponent.hpp:56
static const int liquidPhaseIdx
Index of the liquid phase.
Definition: H2ON2LiquidPhaseFluidSystem.hpp:79
static Evaluation liquidEnthalpy(const Evaluation &temperature, const Evaluation &pressure)
Specific enthalpy of the liquid .
Definition: TabulatedComponent.hpp:292
#define OPM_OPTIM_UNUSED
Definition: Unused.hpp:42
static const int H2OIdx
The index of the water component.
Definition: H2ON2LiquidPhaseFluidSystem.hpp:129
static const bool isTabulated
Definition: TabulatedComponent.hpp:61
static const Scalar acentricFactor()
The acentric factor of water.
Definition: H2O.hpp:85
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: H2ON2LiquidPhaseFluidSystem.hpp:305
Relations valid for an ideal gas.
A generic class which tabulates all thermodynamic properties of a given component.
SimpleN2 N2
The type of the component for pure molecular nitrogen.
Definition: H2ON2LiquidPhaseFluidSystem.hpp:139
static const int numPhases
Number of fluid phases in the fluid system.
Definition: H2ON2LiquidPhaseFluidSystem.hpp:76
static bool liquidIsCompressible()
Returns true iff the liquid phase is assumed to be compressible.
Definition: TabulatedComponent.hpp:411
static const char * componentName(unsigned compIdx)
Return the human readable name of a component.
Definition: H2ON2LiquidPhaseFluidSystem.hpp:142
static Scalar acentricFactor(unsigned compIdx)
The acentric factor of a component [].
Definition: H2ON2LiquidPhaseFluidSystem.hpp:199
static void init()
Initialize the fluid system's static parameters.
Definition: H2ON2LiquidPhaseFluidSystem.hpp:219
static Scalar criticalTemperature(unsigned compIdx)
Critical temperature of a component [K].
Definition: H2ON2LiquidPhaseFluidSystem.hpp:169
The base class for all fluid systems.
static LhsEval density(const FluidState &fluidState, const ParameterCache &, unsigned phaseIdx)
Calculate the density [kg/m^3] of a fluid phase.
Definition: H2ON2LiquidPhaseFluidSystem.hpp:251
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 bool isCompressible(unsigned)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: H2ON2LiquidPhaseFluidSystem.hpp:97
static Scalar molarMass(unsigned compIdx)
Return the molar mass of a component in [kg/mol].
Definition: H2ON2LiquidPhaseFluidSystem.hpp:154
static LhsEval heatCapacity(const FluidState &fluidState, const ParameterCache &, unsigned phaseIdx)
Specific isobaric heat capacity of a fluid phase [J/kg].
Definition: H2ON2LiquidPhaseFluidSystem.hpp:381
static const char * name()
A human readable name for nitrogen.
Definition: N2.hpp:57
NullParameterCache ParameterCache
The type of the fluid system's parameter cache.
Definition: H2ON2LiquidPhaseFluidSystem.hpp:69