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"
32
40
41#include <iostream>
42#include <cassert>
43
44namespace Opm {
45
51template <class Scalar>
53 : public BaseFluidSystem<Scalar, H2ON2FluidSystem<Scalar> >
54{
57
58 // convenience typedefs
59 typedef ::Opm::IdealGas<Scalar> IdealGas;
60 typedef ::Opm::H2O<Scalar> IapwsH2O;
62 typedef ::Opm::N2<Scalar> SimpleN2;
63
64public:
66 template <class Evaluation>
68
69 /****************************************
70 * Fluid phase related static parameters
71 ****************************************/
72
74 static const int numPhases = 2;
75
77 static const int liquidPhaseIdx = 0;
79 static const int gasPhaseIdx = 1;
80
82 static const char* phaseName(unsigned phaseIdx)
83 {
84 static const char* name[] = {
85 "liquid",
86 "gas"
87 };
88
89 assert(phaseIdx < numPhases);
90 return name[phaseIdx];
91 }
92
94 static bool isLiquid(unsigned phaseIdx)
95 {
96 //assert(0 <= phaseIdx && phaseIdx < numPhases);
97 return phaseIdx != gasPhaseIdx;
98 }
99
101 static bool isCompressible(unsigned phaseIdx)
102 {
103 //assert(0 <= phaseIdx && phaseIdx < numPhases);
104 // gases are always compressible
105 return
106 (phaseIdx == gasPhaseIdx)
107 ? true
108 :H2O::liquidIsCompressible();// the water component decides for the liquid phase...
109 }
110
112 static bool isIdealGas(unsigned phaseIdx)
113 {
114 //assert(0 <= phaseIdx && phaseIdx < numPhases);
115
116 return
117 (phaseIdx == gasPhaseIdx)
118 ? H2O::gasIsIdeal() && N2::gasIsIdeal() // let the components decide
119 : false; // not a gas
120 }
121
123 static bool isIdealMixture(unsigned /*phaseIdx*/)
124 {
125 //assert(0 <= phaseIdx && phaseIdx < numPhases);
126 // we assume Henry's and Rault's laws for the water phase and
127 // and no interaction between gas molecules of different
128 // components, so all phases are ideal mixtures!
129 return true;
130 }
131
132 /****************************************
133 * Component related static parameters
134 ****************************************/
135
137 static const int numComponents = 2;
138
140 static const int H2OIdx = 0;
142 static const int N2Idx = 1;
143
146 //typedef SimpleH2O H2O;
147 //typedef IapwsH2O H2O;
148
150 typedef SimpleN2 N2;
151
153 static const char* componentName(unsigned compIdx)
154 {
155 static const char* name[] = {
156 H2O::name(),
157 N2::name()
158 };
159
160 assert(compIdx < numComponents);
161 return name[compIdx];
162 }
163
165 static Scalar molarMass(unsigned compIdx)
166 {
167 //assert(0 <= compIdx && compIdx < numComponents);
168 return (compIdx == H2OIdx)
170 : (compIdx == N2Idx)
171 ? N2::molarMass()
172 : 1e30;
173 }
174
180 static Scalar criticalTemperature(unsigned compIdx)
181 {
182 return (compIdx == H2OIdx)
184 : (compIdx == N2Idx)
186 : 1e30;
187 }
188
194 static Scalar criticalPressure(unsigned compIdx)
195 {
196 return (compIdx == H2OIdx)
198 : (compIdx == N2Idx)
200 : 1e30;
201 }
202
208 static Scalar acentricFactor(unsigned compIdx)
209 {
210 return (compIdx == H2OIdx)
212 : (compIdx == N2Idx)
214 : 1e30;
215 }
216
217 /****************************************
218 * thermodynamic relations
219 ****************************************/
220
227 static void init()
228 {
229 init(/*tempMin=*/273.15,
230 /*tempMax=*/623.15,
231 /*numTemp=*/50,
232 /*pMin=*/0.0,
233 /*pMax=*/20e6,
234 /*numP=*/50);
235 }
236
248 static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp,
249 Scalar pressMin, Scalar pressMax, unsigned nPress)
250 {
251 if (H2O::isTabulated) {
252 TabulatedH2O::init(tempMin, tempMax, nTemp,
253 pressMin, pressMax, nPress);
254 }
255 }
256
260 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
261 static LhsEval density(const FluidState& fluidState,
262 const ParameterCache<ParamCacheEval>& /*paramCache*/,
263 unsigned phaseIdx)
264 {
265 assert(phaseIdx < numPhases);
266
267 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
268 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
269
270 LhsEval sumMoleFrac = 0;
271 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
272 sumMoleFrac += decay<LhsEval>(fluidState.moleFraction(phaseIdx, compIdx));
273
274 // liquid phase
275 if (phaseIdx == liquidPhaseIdx) {
276 // assume ideal mixture where each molecule occupies the same volume regardless
277 // of whether it is water or nitrogen.
278 const LhsEval& clH2O = H2O::liquidDensity(T, p)/H2O::molarMass();
279
280 const auto& xlH2O = decay<LhsEval>(fluidState.moleFraction(liquidPhaseIdx, H2OIdx));
281 const auto& xlN2 = decay<LhsEval>(fluidState.moleFraction(liquidPhaseIdx, N2Idx));
282
283 return clH2O*(H2O::molarMass()*xlH2O + N2::molarMass()*xlN2)/sumMoleFrac;
284 }
285
286 // gas phase
287 assert(phaseIdx == gasPhaseIdx);
288
289 // assume ideal mixture: steam and nitrogen don't "distinguish" each other
290 const auto& xgH2O = decay<LhsEval>(fluidState.moleFraction(gasPhaseIdx, H2OIdx));
291 const auto& xgN2 = decay<LhsEval>(fluidState.moleFraction(gasPhaseIdx, N2Idx));
292 const auto& rho_gH2O = H2O::gasDensity(T, p*xgH2O);
293 const auto& rho_gN2 = N2::gasDensity(T, p*xgN2);
294 return (rho_gH2O + rho_gN2)/max(1e-5, sumMoleFrac);
295 }
296
298 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
299 static LhsEval viscosity(const FluidState& fluidState,
300 const ParameterCache<ParamCacheEval>& /*paramCache*/,
301 unsigned phaseIdx)
302 {
303 assert(phaseIdx < numPhases);
304
305 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
306 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
307
308 // liquid phase
309 if (phaseIdx == liquidPhaseIdx)
310 // assume pure water for the liquid phase
311 return H2O::liquidViscosity(T, p);
312
313 // gas phase
314 assert(phaseIdx == gasPhaseIdx);
315
316 /* Wilke method. See:
317 *
318 * See: R. Reid, et al.: The Properties of Gases and Liquids,
319 * 4th edition, McGraw-Hill, 1987, 407-410
320 * 5th edition, McGraw-Hill, 20001, p. 9.21/22
321 */
322 LhsEval muResult = 0;
323 const LhsEval mu[numComponents] = {
325 N2::gasViscosity(T, p)
326 };
327
328 LhsEval sumx = 0.0;
329 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
330 sumx += decay<LhsEval>(fluidState.moleFraction(phaseIdx, compIdx));
331 sumx = max(1e-10, sumx);
332
333 for (unsigned i = 0; i < numComponents; ++i) {
334 LhsEval divisor = 0;
335 for (unsigned j = 0; j < numComponents; ++j) {
336 LhsEval phiIJ = 1 + sqrt(mu[i]/mu[j]) * std::pow(molarMass(j)/molarMass(i), 1/4.0);
337 phiIJ *= phiIJ;
338 phiIJ /= std::sqrt(8*(1 + molarMass(i)/molarMass(j)));
339 divisor +=
340 decay<LhsEval>(fluidState.moleFraction(phaseIdx, j))
341 /sumx*phiIJ;
342 }
343 muResult +=
344 decay<LhsEval>(fluidState.moleFraction(phaseIdx, i))
345 /sumx*mu[i]/divisor;
346 }
347 return muResult;
348 }
349
351 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
352 static LhsEval fugacityCoefficient(const FluidState& fluidState,
353 const ParameterCache<ParamCacheEval>& /*paramCache*/,
354 unsigned phaseIdx,
355 unsigned compIdx)
356 {
357 assert(phaseIdx < numPhases);
358 assert(compIdx < numComponents);
359
360 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
361 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
362
363 // liquid phase
364 if (phaseIdx == liquidPhaseIdx) {
365 if (compIdx == H2OIdx)
366 return H2O::vaporPressure(T)/p;
367 return BinaryCoeff::H2O_N2::henry(T)/p;
368 }
369
370 assert(phaseIdx == gasPhaseIdx);
371
372 // for the gas phase, assume an ideal gas when it comes to
373 // fugacity (-> fugacity == partial pressure)
374 return 1.0;
375 }
376
378 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
379 static LhsEval diffusionCoefficient(const FluidState& fluidState,
380 const ParameterCache<ParamCacheEval>& /*paramCache*/,
381 unsigned phaseIdx,
382 unsigned /*compIdx*/)
383
384 {
385 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
386 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
387
388 // liquid phase
389 if (phaseIdx == liquidPhaseIdx)
391
392 // gas phase
393 assert(phaseIdx == gasPhaseIdx);
395 }
396
398 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
399 static LhsEval enthalpy(const FluidState& fluidState,
400 const ParameterCache<ParamCacheEval>& /*paramCache*/,
401 unsigned phaseIdx)
402 {
403 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
404 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
407
408 // liquid phase
409 if (phaseIdx == liquidPhaseIdx) {
410 // TODO: correct way to deal with the solutes???
411 return H2O::liquidEnthalpy(T, p);
412 }
413
414 // gas phase
415 assert(phaseIdx == gasPhaseIdx);
416
417 // assume ideal mixture: Molecules of one component don't discriminate between
418 // their own kind and molecules of the other component.
419 const auto& XgH2O = decay<LhsEval>(fluidState.massFraction(gasPhaseIdx, H2OIdx));
420 const auto& XgN2 = decay<LhsEval>(fluidState.massFraction(gasPhaseIdx, N2Idx));
421
422 LhsEval hH2O = XgH2O*H2O::gasEnthalpy(T, p);
423 LhsEval hN2 = XgN2*N2::gasEnthalpy(T, p);
424 return hH2O + hN2;
425 }
426
428 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
429 static LhsEval thermalConductivity(const FluidState& fluidState,
430 const ParameterCache<ParamCacheEval>& /*paramCache*/,
431 unsigned phaseIdx)
432 {
433 assert(phaseIdx < numPhases);
434
435 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
436 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
437 if (phaseIdx == liquidPhaseIdx) // liquid phase
439
440 // gas phase
441 assert(phaseIdx == gasPhaseIdx);
442
443 // return the sum of the partial conductivity of Nitrogen and Steam
444 const auto& xH2O = decay<LhsEval>(fluidState.moleFraction(phaseIdx, H2OIdx));
445 const auto& xN2 = decay<LhsEval>(fluidState.moleFraction(phaseIdx, N2Idx));
446
447 // Assuming Raoult's, Daltons law and ideal gas in order to obtain the
448 // partial pressures in the gas phase
449 const auto& lambdaN2 = N2::gasThermalConductivity(T, p*xN2);
450 const auto& lambdaH2O = H2O::gasThermalConductivity(T, p*xH2O);
451
452 return lambdaN2 + lambdaH2O;
453 }
454
456 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
457 static LhsEval heatCapacity(const FluidState& fluidState,
458 const ParameterCache<ParamCacheEval>& /*paramCache*/,
459 unsigned phaseIdx)
460 {
461 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
462 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
463 const auto& xAlphaH2O = decay<LhsEval>(fluidState.moleFraction(phaseIdx, H2OIdx));
464 const auto& xAlphaN2 = decay<LhsEval>(fluidState.moleFraction(phaseIdx, N2Idx));
465 const auto& XAlphaH2O = decay<LhsEval>(fluidState.massFraction(phaseIdx, H2OIdx));
466 const auto& XAlphaN2 = decay<LhsEval>(fluidState.massFraction(phaseIdx, N2Idx));
467
468 if (phaseIdx == liquidPhaseIdx)
469 return H2O::liquidHeatCapacity(T, p);
470
471 assert(phaseIdx == gasPhaseIdx);
472
473 // for the gas phase, assume ideal mixture
474 LhsEval c_pN2;
475 LhsEval c_pH2O;
476
477 c_pN2 = N2::gasHeatCapacity(T, p*xAlphaN2);
478 c_pH2O = H2O::gasHeatCapacity(T, p*xAlphaH2O);
479
480 // mingle both components together. this assumes that there is no "cross
481 // interaction" between both flavors of molecules.
482 return XAlphaH2O*c_pH2O + XAlphaN2*c_pN2;
483 }
484};
485
486} // namespace Opm
487
488#endif
Some templates to wrap the valgrind client request macros.
The base class for all fluid systems.
Definition: BaseFluidSystem.hpp:44
Scalar Scalar
The type used for scalar quantities.
Definition: BaseFluidSystem.hpp:49
static Evaluation liquidDiffCoeff(const Evaluation &temperature, const Evaluation &)
Diffusion coefficent for molecular nitrogen in liquid water.
Definition: H2O_N2.hpp:102
static Evaluation gasDiffCoeff(const Evaluation &temperature, const Evaluation &pressure)
Binary diffusion coefficent for molecular water and nitrogen.
Definition: H2O_N2.hpp:70
static Evaluation henry(const Evaluation &temperature)
Henry coefficent for molecular nitrogen in liquid water.
Definition: H2O_N2.hpp:52
A two-phase fluid system with water and nitrogen as components.
Definition: H2ON2FluidSystem.hpp:54
static const int H2OIdx
The component index of water.
Definition: H2ON2FluidSystem.hpp:140
static Scalar acentricFactor(unsigned compIdx)
The acentric factor of a component [].
Definition: H2ON2FluidSystem.hpp:208
static LhsEval thermalConductivity(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Thermal conductivity of a fluid phase [W/(m K)].
Definition: H2ON2FluidSystem.hpp:429
static Scalar criticalTemperature(unsigned compIdx)
Critical temperature of a component [K].
Definition: H2ON2FluidSystem.hpp:180
TabulatedH2O H2O
The component for pure water.
Definition: H2ON2FluidSystem.hpp:145
static void init()
Initialize the fluid system's static parameters.
Definition: H2ON2FluidSystem.hpp:227
static const char * phaseName(unsigned phaseIdx)
Return the human readable name of a fluid phase.
Definition: H2ON2FluidSystem.hpp:82
static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp, Scalar pressMin, Scalar pressMax, unsigned nPress)
Initialize the fluid system's static parameters using problem specific temperature and pressure range...
Definition: H2ON2FluidSystem.hpp:248
static LhsEval density(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Calculate the density [kg/m^3] of a fluid phase.
Definition: H2ON2FluidSystem.hpp:261
static bool isIdealMixture(unsigned)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: H2ON2FluidSystem.hpp:123
static bool isLiquid(unsigned phaseIdx)
Return whether a phase is liquid.
Definition: H2ON2FluidSystem.hpp:94
static const int liquidPhaseIdx
Index of the liquid phase.
Definition: H2ON2FluidSystem.hpp:77
static const int numPhases
Number of fluid phases in the fluid system.
Definition: H2ON2FluidSystem.hpp:74
static const int gasPhaseIdx
Index of the gas phase.
Definition: H2ON2FluidSystem.hpp:79
static LhsEval heatCapacity(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Specific isobaric heat capacity of a fluid phase [J/kg].
Definition: H2ON2FluidSystem.hpp:457
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:379
static Scalar criticalPressure(unsigned compIdx)
Critical pressure of a component [Pa].
Definition: H2ON2FluidSystem.hpp:194
static LhsEval enthalpy(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Given a phase's composition, temperature, pressure and density, calculate its specific enthalpy [J/kg...
Definition: H2ON2FluidSystem.hpp:399
static const char * componentName(unsigned compIdx)
Return the human readable name of a component.
Definition: H2ON2FluidSystem.hpp:153
static const int numComponents
Number of chemical species in the fluid system.
Definition: H2ON2FluidSystem.hpp:137
static Scalar molarMass(unsigned compIdx)
Return the molar mass of a component in [kg/mol].
Definition: H2ON2FluidSystem.hpp:165
static bool isIdealGas(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: H2ON2FluidSystem.hpp:112
static LhsEval viscosity(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Calculate the dynamic viscosity of a fluid phase [Pa*s].
Definition: H2ON2FluidSystem.hpp:299
SimpleN2 N2
The component for pure nitrogen.
Definition: H2ON2FluidSystem.hpp:150
static const int N2Idx
The component index of molecular nitrogen.
Definition: H2ON2FluidSystem.hpp:142
static bool isCompressible(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: H2ON2FluidSystem.hpp:101
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:352
Material properties of pure water .
Definition: H2O.hpp:62
Relations valid for an ideal gas.
Definition: IdealGas.hpp:38
Properties of pure molecular nitrogen .
Definition: N2.hpp:49
static Evaluation gasThermalConductivity(const Evaluation &, const Evaluation &)
Specific heat conductivity of steam .
Definition: N2.hpp:303
static Evaluation gasDensity(const Evaluation &temperature, const Evaluation &pressure)
The density of gas at a given pressure and temperature.
Definition: N2.hpp:145
static bool gasIsIdeal()
Returns true iff the gas phase is assumed to be ideal.
Definition: N2.hpp:160
static Scalar criticalPressure()
Returns the critical pressure of molecular nitrogen.
Definition: N2.hpp:74
static Scalar criticalTemperature()
Returns the critical temperature of molecular nitrogen.
Definition: N2.hpp:68
static Evaluation gasViscosity(const Evaluation &temperature, const Evaluation &)
The dynamic viscosity of at a given pressure and temperature.
Definition: N2.hpp:267
static const char * name()
A human readable name for nitrogen.
Definition: N2.hpp:56
static Scalar acentricFactor()
Acentric factor of .
Definition: N2.hpp:85
static Scalar molarMass()
The molar mass in of molecular nitrogen.
Definition: N2.hpp:62
static Evaluation gasHeatCapacity(const Evaluation &temperature, const Evaluation &)
Specific isobaric heat capacity of pure nitrogen gas.
Definition: N2.hpp:236
static Evaluation gasEnthalpy(const Evaluation &temperature, const Evaluation &)
Specific enthalpy of pure nitrogen gas.
Definition: N2.hpp:186
A parameter cache which does nothing.
Definition: NullParameterCache.hpp:40
A generic class which tabulates all thermodynamic properties of a given component.
Definition: TabulatedComponent.hpp:56
static Evaluation gasThermalConductivity(const Evaluation &temperature, const Evaluation &pressure)
The thermal conductivity of gaseous water .
Definition: TabulatedComponent.hpp:495
static Evaluation gasEnthalpy(const Evaluation &temperature, const Evaluation &pressure)
Specific enthalpy of the gas .
Definition: TabulatedComponent.hpp:282
static Evaluation liquidHeatCapacity(const Evaluation &temperature, const Evaluation &pressure)
Specific isobaric heat capacity of the liquid .
Definition: TabulatedComponent.hpp:333
static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp, Scalar pressMin, Scalar pressMax, unsigned nPress)
Initialize the tables.
Definition: TabulatedComponent.hpp:72
static Scalar criticalTemperature()
Returns the critical temperature in of the component.
Definition: TabulatedComponent.hpp:227
static Scalar criticalPressure()
Returns the critical pressure in of the component.
Definition: TabulatedComponent.hpp:233
static Scalar molarMass()
The molar mass in of the component.
Definition: TabulatedComponent.hpp:221
static const bool isTabulated
Definition: TabulatedComponent.hpp:60
static bool liquidIsCompressible()
Returns true iff the liquid phase is assumed to be compressible.
Definition: TabulatedComponent.hpp:408
static Evaluation liquidViscosity(const Evaluation &temperature, const Evaluation &pressure)
The dynamic viscosity of liquid.
Definition: TabulatedComponent.hpp:478
static Evaluation gasDensity(const Evaluation &temperature, const Evaluation &pressure)
The density of gas at a given pressure and temperature .
Definition: TabulatedComponent.hpp:426
static Evaluation gasHeatCapacity(const Evaluation &temperature, const Evaluation &pressure)
Specific isobaric heat capacity of the gas .
Definition: TabulatedComponent.hpp:316
static bool gasIsIdeal()
Returns true iff the gas phase is assumed to be ideal.
Definition: TabulatedComponent.hpp:414
static Scalar acentricFactor()
Returns the acentric factor of the component.
Definition: TabulatedComponent.hpp:239
static Evaluation liquidDensity(const Evaluation &temperature, const Evaluation &pressure)
The density of liquid at a given pressure and temperature .
Definition: TabulatedComponent.hpp:444
static Evaluation vaporPressure(const Evaluation &temperature)
The vapor pressure in of the component at a given temperature.
Definition: TabulatedComponent.hpp:267
static Evaluation liquidThermalConductivity(const Evaluation &temperature, const Evaluation &pressure)
The thermal conductivity of liquid water .
Definition: TabulatedComponent.hpp:512
static Evaluation gasViscosity(const Evaluation &temperature, const Evaluation &pressure)
The dynamic viscosity of gas.
Definition: TabulatedComponent.hpp:461
static Evaluation liquidEnthalpy(const Evaluation &temperature, const Evaluation &pressure)
Specific enthalpy of the liquid .
Definition: TabulatedComponent.hpp:299
static const char * name()
A human readable name for the component.
Definition: TabulatedComponent.hpp:215
bool CheckDefined(const T &value)
Make valgrind complain if any of the memory occupied by an object is undefined.
Definition: Valgrind.hpp:74
Definition: Air_Mesitylene.hpp:34
Evaluation sqrt(const Evaluation &value)
Definition: MathToolbox.hpp:399
ReturnEval_< Evaluation1, Evaluation2 >::type max(const Evaluation1 &arg1, const Evaluation2 &arg2)
Definition: MathToolbox.hpp:341
ReturnEval_< Evaluation1, Evaluation2 >::type pow(const Evaluation1 &base, const Evaluation2 &exp)
Definition: MathToolbox.hpp:416