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