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 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_AIR_XYLENE_FLUID_SYSTEM_HPP
28#define OPM_H2O_AIR_XYLENE_FLUID_SYSTEM_HPP
29
35
39
40#include "BaseFluidSystem.hpp"
42
43namespace Opm {
44
50template <class Scalar>
52 : public BaseFluidSystem<Scalar, H2OAirXyleneFluidSystem<Scalar> >
53{
56
57public:
58 template <class Evaluation>
59 struct ParameterCache : public NullParameterCache<Evaluation>
60 {};
61
63 typedef ::Opm::H2O<Scalar> H2O;
67 typedef ::Opm::Air<Scalar> Air;
68
70 static const int numPhases = 3;
72 static const int numComponents = 3;
73
75 static const int waterPhaseIdx = 0;
77 static const int naplPhaseIdx = 1;
79 static const int gasPhaseIdx = 2;
80
82 static const int H2OIdx = 0;
84 static const int NAPLIdx = 1;
86 static const int airIdx = 2;
87
89 static void init()
90 { }
91
93 static bool isLiquid(unsigned phaseIdx)
94 {
95 //assert(0 <= phaseIdx && phaseIdx < numPhases);
96 return phaseIdx != gasPhaseIdx;
97 }
98
100 static bool isIdealGas(unsigned phaseIdx)
101 { return phaseIdx == gasPhaseIdx && H2O::gasIsIdeal() && Air::gasIsIdeal() && NAPL::gasIsIdeal(); }
102
104 static bool isIdealMixture(unsigned /*phaseIdx*/)
105 {
106 //assert(0 <= phaseIdx && phaseIdx < numPhases);
107
108 // we assume Henry's and Rault's laws for the water phase and
109 // and no interaction between gas molecules of different
110 // components, so all phases are ideal mixtures!
111 return true;
112 }
113
115 static bool isCompressible(unsigned phaseIdx)
116 {
117 return
118 (phaseIdx == gasPhaseIdx)
119 // gases are always compressible
120 ? true
121 : (phaseIdx == waterPhaseIdx)
122 // the water component decides for the water phase...
124 // the NAPL component decides for the napl phase...
126 }
127
129 static const char* phaseName(unsigned phaseIdx)
130 {
131 switch (phaseIdx) {
132 case waterPhaseIdx: return "water";
133 case naplPhaseIdx: return "napl";
134 case gasPhaseIdx: return "gas";
135 };
136 throw std::logic_error("Invalid phase index "+std::to_string(phaseIdx));
137 }
138
140 static const char* componentName(unsigned compIdx)
141 {
142 switch (compIdx) {
143 case H2OIdx: return H2O::name();
144 case airIdx: return Air::name();
145 case NAPLIdx: return NAPL::name();
146 };
147 throw std::logic_error("Invalid component index "+std::to_string(compIdx));
148 }
149
151 static Scalar molarMass(unsigned compIdx)
152 {
153 return
154 (compIdx == H2OIdx)
155 // gases are always compressible
157 : (compIdx == airIdx)
158 // the water component decides for the water comp...
160 // the NAPL component decides for the napl comp...
161 : (compIdx == NAPLIdx)
163 : 1e30;
164 }
165
167 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
168 static LhsEval density(const FluidState& fluidState,
169 const ParameterCache<ParamCacheEval>& /*paramCache*/,
170 unsigned phaseIdx)
171 {
172 if (phaseIdx == waterPhaseIdx) {
173 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
174 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
175
176 // See: Ochs 2008
177 // \todo: proper citation
178 const LhsEval& rholH2O = H2O::liquidDensity(T, p);
179 const LhsEval& clH2O = rholH2O/H2O::molarMass();
180
181 const auto& xwH2O = decay<LhsEval>(fluidState.moleFraction(waterPhaseIdx, H2OIdx));
182 const auto& xwAir = decay<LhsEval>(fluidState.moleFraction(waterPhaseIdx, airIdx));
183 const auto& xwNapl = decay<LhsEval>(fluidState.moleFraction(waterPhaseIdx, NAPLIdx));
184 // this assumes each dissolved molecule displaces exactly one
185 // water molecule in the liquid
186 return clH2O*(H2O::molarMass()*xwH2O + Air::molarMass()*xwAir + NAPL::molarMass()*xwNapl);
187 }
188 else if (phaseIdx == naplPhaseIdx) {
189 // assume pure NAPL for the NAPL phase
190 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
191 return NAPL::liquidDensity(T, LhsEval(1e30));
192 }
193
194 assert (phaseIdx == gasPhaseIdx);
195 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
196 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
197
198 const LhsEval& pH2O = decay<LhsEval>(fluidState.moleFraction(gasPhaseIdx, H2OIdx))*p;
199 const LhsEval& pAir = decay<LhsEval>(fluidState.moleFraction(gasPhaseIdx, airIdx))*p;
200 const LhsEval& pNAPL = decay<LhsEval>(fluidState.moleFraction(gasPhaseIdx, NAPLIdx))*p;
201 return
202 H2O::gasDensity(T, pH2O) +
203 Air::gasDensity(T, pAir) +
204 NAPL::gasDensity(T, pNAPL);
205 }
206
208 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
209 static LhsEval viscosity(const FluidState& fluidState,
210 const ParameterCache<ParamCacheEval>& /*paramCache*/,
211 unsigned phaseIdx)
212 {
213 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
214 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
215
216 if (phaseIdx == waterPhaseIdx) {
217 // assume pure water viscosity
218 return H2O::liquidViscosity(T, p);
219 }
220 else if (phaseIdx == naplPhaseIdx) {
221 // assume pure NAPL viscosity
222 return NAPL::liquidViscosity(T, p);
223 }
224
225 assert (phaseIdx == gasPhaseIdx);
226
227 /* Wilke method. See:
228 *
229 * See: R. Reid, et al.: The Properties of Gases and Liquids,
230 * 4th edition, McGraw-Hill, 1987, 407-410
231 * 5th edition, McGraw-Hill, 20001, p. 9.21/22
232 *
233 * in this case, we use a simplified version in order to avoid
234 * computationally costly evaluation of sqrt and pow functions and
235 * divisions
236 * -- compare e.g. with Promo Class p. 32/33
237 */
238 const LhsEval mu[numComponents] = {
242 };
243 // molar masses
244 const Scalar M[numComponents] = {
248 };
249
250 const auto& xgAir = decay<LhsEval>(fluidState.moleFraction(gasPhaseIdx, airIdx));
251 const auto& xgH2O = decay<LhsEval>(fluidState.moleFraction(gasPhaseIdx, H2OIdx));
252 const auto& xgNapl = decay<LhsEval>(fluidState.moleFraction(gasPhaseIdx, NAPLIdx));
253
254 const LhsEval& xgAW = xgAir + xgH2O;
255 const LhsEval& muAW = (mu[airIdx]*xgAir + mu[H2OIdx]*xgH2O)/ xgAW;
256
257 const LhsEval& MAW = (xgAir*Air::molarMass() + xgH2O*H2O::molarMass())/xgAW;
258
259 Scalar phiCAW = 0.3; // simplification for this particular system
260 /* actually like this
261 * Scalar phiCAW = std::pow(1.+std::sqrt(mu[NAPLIdx]/muAW)*std::pow(MAW/M[NAPLIdx],0.25),2)
262 * / std::sqrt(8.*(1.+M[NAPLIdx]/MAW));
263 */
264 const LhsEval& phiAWC = phiCAW * muAW*M[NAPLIdx]/(mu[NAPLIdx]*MAW);
265
266 return (xgAW*muAW)/(xgAW+xgNapl*phiAWC) + (xgNapl*mu[NAPLIdx])/(xgNapl + xgAW*phiCAW);
267 }
268
270 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
271 static LhsEval diffusionCoefficient(const FluidState& fluidState,
272 const ParameterCache<ParamCacheEval>& /*paramCache*/,
273 unsigned phaseIdx,
274 unsigned compIdx)
275 {
276 if (phaseIdx==gasPhaseIdx) {
277 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
278 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
279
280 const LhsEval& diffAC = BinaryCoeff::Air_Xylene::gasDiffCoeff(T, p);
281 const LhsEval& diffWC = BinaryCoeff::H2O_Xylene::gasDiffCoeff(T, p);
282 const LhsEval& diffAW = BinaryCoeff::H2O_Air::gasDiffCoeff(T, p);
283
284 const LhsEval& xga = decay<LhsEval>(fluidState.moleFraction(gasPhaseIdx, airIdx));
285 const LhsEval& xgw = decay<LhsEval>(fluidState.moleFraction(gasPhaseIdx, H2OIdx));
286 const LhsEval& xgc = decay<LhsEval>(fluidState.moleFraction(gasPhaseIdx, NAPLIdx));
287
288 if (compIdx==NAPLIdx) return (1.- xgw)/(xga/diffAW + xgc/diffWC);
289 else if (compIdx==H2OIdx) return (1.- xgc)/(xgw/diffWC + xga/diffAC);
290 else if (compIdx==airIdx) throw std::logic_error("Diffusivity of air in the gas phase "
291 "is constraint by sum of diffusive fluxes = 0 !\n");
292 } else if (phaseIdx==waterPhaseIdx){
293 Scalar diffACl = 1.e-9; // BinaryCoeff::Air_Xylene::liquidDiffCoeff(temperature, pressure);
294 Scalar diffWCl = 1.e-9; // BinaryCoeff::H2O_Xylene::liquidDiffCoeff(temperature, pressure);
295 Scalar diffAWl = 1.e-9; // BinaryCoeff::H2O_Air::liquidDiffCoeff(temperature, pressure);
296
297 const LhsEval& xwa = decay<LhsEval>(fluidState.moleFraction(waterPhaseIdx, airIdx));
298 const LhsEval& xww = decay<LhsEval>(fluidState.moleFraction(waterPhaseIdx, H2OIdx));
299 const LhsEval& xwc = decay<LhsEval>(fluidState.moleFraction(waterPhaseIdx, NAPLIdx));
300
301 switch (compIdx) {
302 case NAPLIdx:
303 return (1.- xww)/(xwa/diffAWl + xwc/diffWCl);
304 case airIdx:
305 return (1.- xwc)/(xww/diffWCl + xwa/diffACl);
306 case H2OIdx:
307 throw std::logic_error("Diffusivity of water in the water phase "
308 "is constraint by sum of diffusive fluxes = 0 !\n");
309 };
310 } else if (phaseIdx==naplPhaseIdx) {
311
312 throw std::logic_error("Diffusion coefficients of "
313 "substances in liquid phase are undefined!\n");
314 }
315 return 0;
316 }
317
319 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
320 static LhsEval fugacityCoefficient(const FluidState& fluidState,
321 const ParameterCache<ParamCacheEval>& /*paramCache*/,
322 unsigned phaseIdx,
323 unsigned compIdx)
324 {
325 assert(phaseIdx < numPhases);
326 assert(compIdx < numComponents);
327
328 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
329 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
330
331 if (phaseIdx == waterPhaseIdx) {
332 if (compIdx == H2OIdx)
333 return H2O::vaporPressure(T)/p;
334 else if (compIdx == airIdx)
335 return BinaryCoeff::H2O_Air::henry(T)/p;
336 else if (compIdx == NAPLIdx)
338 }
339
340 // for the NAPL phase, we assume currently that nothing is
341 // dissolved. this means that the affinity of the NAPL
342 // component to the NAPL phase is much higher than for the
343 // other components, i.e. the fugacity cofficient is much
344 // smaller.
345 if (phaseIdx == naplPhaseIdx) {
346 const LhsEval& phiNapl = NAPL::vaporPressure(T)/p;
347 if (compIdx == NAPLIdx)
348 return phiNapl;
349 else if (compIdx == airIdx)
350 return 1e6*phiNapl;
351 else if (compIdx == H2OIdx)
352 return 1e6*phiNapl;
353 }
354
355 // for the gas phase, assume an ideal gas when it comes to
356 // fugacity (-> fugacity == partial pressure)
357 assert(phaseIdx == gasPhaseIdx);
358 return 1.0;
359 }
360
362 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
363 static LhsEval enthalpy(const FluidState& fluidState,
364 const ParameterCache<ParamCacheEval>& /*paramCache*/,
365 unsigned phaseIdx)
366 {
367 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
368 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
369
370 if (phaseIdx == waterPhaseIdx) {
371 return H2O::liquidEnthalpy(T, p);
372 }
373 else if (phaseIdx == naplPhaseIdx) {
374 return NAPL::liquidEnthalpy(T, p);
375 }
376 else if (phaseIdx == gasPhaseIdx) { // gas phase enthalpy depends strongly on composition
377 const LhsEval& hgc = NAPL::gasEnthalpy(T, p);
378 const LhsEval& hgw = H2O::gasEnthalpy(T, p);
379 const LhsEval& hga = Air::gasEnthalpy(T, p);
380
381 LhsEval result = 0;
382 result += hgw * decay<LhsEval>(fluidState.massFraction(gasPhaseIdx, H2OIdx));
383 result += hga * decay<LhsEval>(fluidState.massFraction(gasPhaseIdx, airIdx));
384 result += hgc * decay<LhsEval>(fluidState.massFraction(gasPhaseIdx, NAPLIdx));
385
386 return result;
387 }
388 throw std::logic_error("Invalid phase index "+std::to_string(phaseIdx));
389 }
390
391private:
392 template <class LhsEval>
393 static LhsEval waterPhaseDensity_(const LhsEval& T,
394 const LhsEval& pw,
395 const LhsEval& xww,
396 const LhsEval& xwa,
397 const LhsEval& xwc)
398 {
399 const LhsEval& rholH2O = H2O::liquidDensity(T, pw);
400 const LhsEval& clH2O = rholH2O/H2O::molarMass();
401
402 // this assumes each dissolved molecule displaces exactly one
403 // water molecule in the liquid
404 return clH2O*(xww*H2O::molarMass() + xwa*Air::molarMass() + xwc*NAPL::molarMass());
405 }
406
407 template <class LhsEval>
408 static LhsEval gasPhaseDensity_(const LhsEval& T,
409 const LhsEval& pg,
410 const LhsEval& xgw,
411 const LhsEval& xga,
412 const LhsEval& xgc)
413 { return H2O::gasDensity(T, pg*xgw) + Air::gasDensity(T, pg*xga) + NAPL::gasDensity(T, pg*xgc); }
414
415 template <class LhsEval>
416 static LhsEval NAPLPhaseDensity_(const LhsEval& T, const LhsEval& pn)
417 { return NAPL::liquidDensity(T, pn); }
418};
419
420} // namespace Opm
421
422#endif
static Evaluation simpleGasViscosity(const Evaluation &temperature, const Evaluation &)
Definition: Air.hpp:160
static Evaluation gasDensity(const Evaluation &temperature, const Evaluation &pressure)
The density of at a given pressure and temperature [kg/m^3].
Definition: Air.hpp:102
static Scalar molarMass()
The molar mass in of .
Definition: Air.hpp:80
static const char * name()
A human readable name for the .
Definition: Air.hpp:60
static bool gasIsIdeal()
Returns true iff the gas phase is assumed to be ideal.
Definition: Air.hpp:72
static Evaluation gasEnthalpy(const Evaluation &temperature, const Evaluation &)
Specific enthalpy of liquid water with 273.15 K as basis. See: W. Kays, M. Crawford,...
Definition: Air.hpp:180
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 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:57
static Evaluation henry(const Evaluation &temperature)
Henry coefficent for air in liquid water.
Definition: H2O_Air.hpp:55
static Evaluation gasDiffCoeff(const Evaluation &temperature, const Evaluation &pressure)
Binary diffusion coefficent for molecular water and air.
Definition: H2O_Air.hpp:70
static Evaluation gasDiffCoeff(Evaluation temperature, Evaluation pressure)
Binary diffusion coefficent [m^2/s] for molecular water and xylene.
Definition: H2O_Xylene.hpp:66
static Evaluation henry(const Evaluation &)
Henry coefficent for xylene in liquid water.
Definition: H2O_Xylene.hpp:51
A fluid system with water, gas and NAPL as phases and water, air and NAPL (contaminant) as components...
Definition: H2OAirXyleneFluidSystem.hpp:53
static LhsEval viscosity(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Calculate the dynamic viscosity of a fluid phase [Pa*s].
Definition: H2OAirXyleneFluidSystem.hpp:209
static const int numPhases
Number of fluid phases in the fluid system.
Definition: H2OAirXyleneFluidSystem.hpp:70
static const int naplPhaseIdx
The index of the NAPL phase.
Definition: H2OAirXyleneFluidSystem.hpp:77
static Scalar molarMass(unsigned compIdx)
Return the molar mass of a component in [kg/mol].
Definition: H2OAirXyleneFluidSystem.hpp:151
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: H2OAirXyleneFluidSystem.hpp:363
static void init()
Initialize the fluid system's static parameters.
Definition: H2OAirXyleneFluidSystem.hpp:89
static LhsEval density(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Calculate the density [kg/m^3] of a fluid phase.
Definition: H2OAirXyleneFluidSystem.hpp:168
static LhsEval diffusionCoefficient(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx, unsigned compIdx)
Calculate the binary molecular diffusion coefficient for a component in a fluid phase [mol^2 * s / (k...
Definition: H2OAirXyleneFluidSystem.hpp:271
::Opm::H2O< Scalar > H2O
The type of the water component.
Definition: H2OAirXyleneFluidSystem.hpp:63
Xylene< Scalar > NAPL
The type of the xylene/napl component.
Definition: H2OAirXyleneFluidSystem.hpp:65
static bool isIdealMixture(unsigned)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: H2OAirXyleneFluidSystem.hpp:104
static const int waterPhaseIdx
The index of the water phase.
Definition: H2OAirXyleneFluidSystem.hpp:75
static const int numComponents
Number of chemical species in the fluid system.
Definition: H2OAirXyleneFluidSystem.hpp:72
static const char * phaseName(unsigned phaseIdx)
Return the human readable name of a fluid phase.
Definition: H2OAirXyleneFluidSystem.hpp:129
static const int H2OIdx
The index of the water component.
Definition: H2OAirXyleneFluidSystem.hpp:82
static bool isCompressible(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: H2OAirXyleneFluidSystem.hpp:115
static const int NAPLIdx
The index of the NAPL component.
Definition: H2OAirXyleneFluidSystem.hpp:84
::Opm::Air< Scalar > Air
The type of the air component.
Definition: H2OAirXyleneFluidSystem.hpp:67
static const char * componentName(unsigned compIdx)
Return the human readable name of a component.
Definition: H2OAirXyleneFluidSystem.hpp:140
static bool isLiquid(unsigned phaseIdx)
Return whether a phase is liquid.
Definition: H2OAirXyleneFluidSystem.hpp:93
static const int airIdx
The index of the air pseudo-component.
Definition: H2OAirXyleneFluidSystem.hpp:86
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: H2OAirXyleneFluidSystem.hpp:320
static const int gasPhaseIdx
The index of the gas phase.
Definition: H2OAirXyleneFluidSystem.hpp:79
static bool isIdealGas(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: H2OAirXyleneFluidSystem.hpp:100
static Evaluation liquidDensity(const Evaluation &temperature, const Evaluation &pressure, bool extrapolate=false)
The density of pure water in at a given pressure and temperature.
Definition: H2O.hpp:698
static Evaluation gasDensity(const Evaluation &temperature, const Evaluation &pressure)
The density of steam in at a given pressure and temperature.
Definition: H2O.hpp:572
static Evaluation vaporPressure(Evaluation temperature)
The vapor pressure in of pure water at a given temperature.
Definition: H2O.hpp:138
static Evaluation gasViscosity(const Evaluation &temperature, const Evaluation &pressure)
The dynamic viscosity of steam.
Definition: H2O.hpp:802
static Evaluation liquidEnthalpy(const Evaluation &temperature, const Evaluation &pressure)
Specific enthalpy of liquid water .
Definition: H2O.hpp:236
static bool gasIsIdeal()
Returns true iff the gas phase is assumed to be ideal.
Definition: H2O.hpp:638
static const Scalar molarMass()
The molar mass in of water.
Definition: H2O.hpp:80
static bool liquidIsCompressible()
Returns true iff the liquid phase is assumed to be compressible.
Definition: H2O.hpp:556
static Evaluation liquidViscosity(const Evaluation &temperature, const Evaluation &pressure, bool extrapolate=false)
The dynamic viscosity of pure water.
Definition: H2O.hpp:828
static Evaluation gasEnthalpy(const Evaluation &temperature, const Evaluation &pressure)
Specific enthalpy of water steam .
Definition: H2O.hpp:183
static const char * name()
A human readable name for the water.
Definition: H2O.hpp:74
A parameter cache which does nothing.
Definition: NullParameterCache.hpp:40
Component for Xylene.
Definition: Xylene.hpp:47
static Evaluation liquidViscosity(Evaluation temperature, const Evaluation &)
The dynamic liquid viscosity of the pure component.
Definition: Xylene.hpp:309
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:220
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 Evaluation gasEnthalpy(const Evaluation &temperature, const Evaluation &pressure)
Specific enthalpy of the pure component in gas.
Definition: Xylene.hpp:211
static bool liquidIsCompressible()
Returns true iff the liquid phase is assumed to be compressible.
Definition: Xylene.hpp:283
static const char * name()
A human readable name for the xylene.
Definition: Xylene.hpp:55
static Evaluation liquidEnthalpy(const Evaluation &temperature, const Evaluation &pressure)
Specific enthalpy of the pure component in liquid.
Definition: Xylene.hpp:149
static Scalar molarMass()
The molar mass in of xylene.
Definition: Xylene.hpp:61
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:290
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:265
static bool gasIsIdeal()
Returns true iff the gas phase is assumed to be ideal.
Definition: Xylene.hpp:277
Definition: Air_Mesitylene.hpp:34
Definition: H2OAirXyleneFluidSystem.hpp:60