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