H2OAirMesityleneFluidSystem.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_MESITYLENE_FLUID_SYSTEM_HPP
28#define OPM_H2O_AIR_MESITYLENE_FLUID_SYSTEM_HPP
29
30#include "BaseFluidSystem.hpp"
32
44
45#include <iostream>
46
47namespace Opm {
48
54template <class Scalar>
56 : public BaseFluidSystem<Scalar, H2OAirMesityleneFluidSystem<Scalar> >
57{
60
61 typedef ::Opm::H2O<Scalar> IapwsH2O;
62 typedef TabulatedComponent<Scalar, IapwsH2O, /*alongVaporPressure=*/false> TabulatedH2O;
63
64public:
65 template <class Evaluation>
66 struct ParameterCache : public NullParameterCache<Evaluation>
67 {};
68
71
73 typedef ::Opm::Air<Scalar> Air;
74
76 //typedef SimpleH2O H2O;
78 //typedef IapwsH2O H2O;
79
81 static const int numPhases = 3;
83 static const int numComponents = 3;
84
86 static const int waterPhaseIdx = 0;
88 static const int naplPhaseIdx = 1;
90 static const int gasPhaseIdx = 2;
91
93 static const int H2OIdx = 0;
95 static const int NAPLIdx = 1;
97 static const int airIdx = 2;
98
100 static void init()
101 {
102 init(/*tempMin=*/273.15,
103 /*tempMax=*/623.15,
104 /*numTemp=*/50,
105 /*pMin=*/0.0,
106 /*pMax=*/20e6,
107 /*numP=*/50);
108 }
109
121 static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp,
122 Scalar pressMin, Scalar pressMax, unsigned nPress)
123 {
124 if (H2O::isTabulated) {
125 TabulatedH2O::init(tempMin, tempMax, nTemp,
126 pressMin, pressMax, nPress);
127 }
128 }
129
131 static bool isLiquid(unsigned phaseIdx)
132 {
133 //assert(0 <= phaseIdx && phaseIdx < numPhases);
134 return phaseIdx != gasPhaseIdx;
135 }
136
138 static bool isIdealGas(unsigned phaseIdx)
139 { return phaseIdx == gasPhaseIdx && H2O::gasIsIdeal() && Air::gasIsIdeal() && NAPL::gasIsIdeal(); }
140
142 static bool isCompressible(unsigned phaseIdx)
143 {
144 //assert(0 <= phaseIdx && phaseIdx < numPhases);
145 // gases are always compressible
146 return (phaseIdx == gasPhaseIdx)
147 ? true
148 : (phaseIdx == waterPhaseIdx)
151 }
152
154 static bool isIdealMixture(unsigned /*phaseIdx*/)
155 {
156 //assert(0 <= phaseIdx && phaseIdx < numPhases);
157 // we assume Henry's and Rault's laws for the water phase and
158 // and no interaction between gas molecules of different
159 // components, so all phases are ideal mixtures!
160 return true;
161 }
162
164 static const char* phaseName(unsigned phaseIdx)
165 {
166 switch (phaseIdx) {
167 case waterPhaseIdx: return "water";
168 case naplPhaseIdx: return "napl";
169 case gasPhaseIdx: return "gas";
170 };
171 throw std::logic_error("Invalid phase index "+std::to_string(phaseIdx));
172 }
173
175 static const char* componentName(unsigned compIdx)
176 {
177 switch (compIdx) {
178 case H2OIdx: return H2O::name();
179 case airIdx: return Air::name();
180 case NAPLIdx: return NAPL::name();
181 };
182 throw std::logic_error("Invalid component index "+std::to_string(compIdx));
183 }
184
186 static Scalar molarMass(unsigned compIdx)
187 {
188 return
189 (compIdx == H2OIdx)
191 : (compIdx == airIdx)
193 : (compIdx == NAPLIdx)
195 : -1e10;
196 //throw std::logic_error("Invalid component index "+std::to_string(compIdx));
197 }
198
200 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
201 static LhsEval density(const FluidState& fluidState,
202 const ParameterCache<ParamCacheEval>& /*paramCache*/,
203 unsigned phaseIdx)
204 {
205 const LhsEval& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
206
207 if (phaseIdx == waterPhaseIdx) {
208 // See: Ochs 2008
209 const LhsEval& p =
211 ? decay<LhsEval>(fluidState.pressure(phaseIdx))
212 : 1e30;
213
214 const LhsEval& rholH2O = H2O::liquidDensity(T, p);
215 const LhsEval& clH2O = rholH2O/H2O::molarMass();
216
217 // this assumes each dissolved molecule displaces exactly one
218 // water molecule in the liquid
219 return
220 clH2O*(H2O::molarMass()*decay<LhsEval>(fluidState.moleFraction(waterPhaseIdx, H2OIdx)) +
221 Air::molarMass()*decay<LhsEval>(fluidState.moleFraction(waterPhaseIdx, airIdx)) +
222 NAPL::molarMass()*decay<LhsEval>(fluidState.moleFraction(waterPhaseIdx, NAPLIdx)));
223 }
224 else if (phaseIdx == naplPhaseIdx) {
225 // assume pure NAPL for the NAPL phase
226 const LhsEval& p =
228 ? decay<LhsEval>(fluidState.pressure(phaseIdx))
229 : 1e30;
230 return NAPL::liquidDensity(T, p);
231 }
232
233 assert (phaseIdx == gasPhaseIdx);
234 const LhsEval& pg = decay<LhsEval>(fluidState.pressure(gasPhaseIdx));
235 const LhsEval& pH2O = decay<LhsEval>(fluidState.moleFraction(gasPhaseIdx, H2OIdx))*pg;
236 const LhsEval& pAir = decay<LhsEval>(fluidState.moleFraction(gasPhaseIdx, airIdx))*pg;
237 const LhsEval& pNAPL = decay<LhsEval>(fluidState.moleFraction(gasPhaseIdx, NAPLIdx))*pg;
238 return
239 H2O::gasDensity(T, pH2O) +
240 Air::gasDensity(T, pAir) +
241 NAPL::gasDensity(T, pNAPL);
242 }
243
245 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
246 static LhsEval viscosity(const FluidState& fluidState,
247 const ParameterCache<ParamCacheEval>& /*paramCache*/,
248 unsigned phaseIdx)
249 {
250 const LhsEval& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
251 const LhsEval& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
252
253 if (phaseIdx == waterPhaseIdx) {
254 // assume pure water viscosity
255
256 return H2O::liquidViscosity(T,
257 p);
258 }
259 else if (phaseIdx == naplPhaseIdx) {
260 // assume pure NAPL viscosity
261 return NAPL::liquidViscosity(T, p);
262 }
263
264 assert (phaseIdx == gasPhaseIdx);
265
266 /* Wilke method. See:
267 *
268 * See: R. Reid, et al.: The Properties of Gases and Liquids,
269 * 4th edition, McGraw-Hill, 1987, 407-410
270 * 5th edition, McGraw-Hill, 20001, p. 9.21/22
271 *
272 * in this case, we use a simplified version in order to avoid
273 * computationally costly evaluation of sqrt and pow functions and
274 * divisions
275 * -- compare e.g. with Promo Class p. 32/33
276 */
277 const LhsEval mu[numComponents] = {
279 Air::gasViscosity(T, p),
281 };
282 // molar masses
283 const Scalar M[numComponents] = {
287 };
288
289 const LhsEval& xgAir = decay<LhsEval>(fluidState.moleFraction(gasPhaseIdx, airIdx));
290 const LhsEval& xgH2O = decay<LhsEval>(fluidState.moleFraction(gasPhaseIdx, H2OIdx));
291 const LhsEval& xgNapl = decay<LhsEval>(fluidState.moleFraction(gasPhaseIdx, NAPLIdx));
292 const LhsEval& xgAW = xgAir + xgH2O;
293 const LhsEval& muAW = (mu[airIdx]*xgAir + mu[H2OIdx]*xgH2O)/xgAW;
294 const LhsEval& MAW = (xgAir*Air::molarMass() + xgH2O*H2O::molarMass())/xgAW;
295
296 Scalar phiCAW = 0.3; // simplification for this particular system
297 /* actually like this
298 * Scalar phiCAW = std::pow(1.+std::sqrt(mu[NAPLIdx]/muAW)*std::pow(MAW/M[NAPLIdx],0.25),2)
299 * / std::sqrt(8.*(1.+M[NAPLIdx]/MAW));
300 */
301 const LhsEval& phiAWC = phiCAW * muAW*M[NAPLIdx]/(mu[NAPLIdx]*MAW);
302
303 return (xgAW*muAW)/(xgAW + xgNapl*phiAWC) + (xgNapl*mu[NAPLIdx])/(xgNapl + xgAW*phiCAW);
304 }
305
307 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
308 static LhsEval diffusionCoefficient(const FluidState& /*fluidState*/,
309 const ParameterCache<ParamCacheEval>& /*paramCache*/,
310 unsigned /*phaseIdx*/,
311 unsigned /*compIdx*/)
312 {
313 return 0;
314#if 0
316
317 const LhsEval& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
318 const LhsEval& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
319 LhsEval diffCont;
320
321 if (phaseIdx==gasPhaseIdx) {
322 const LhsEval& diffAC = BinaryCoeff::Air_Mesitylene::gasDiffCoeff(T, p);
323 const LhsEval& diffWC = BinaryCoeff::H2O_Mesitylene::gasDiffCoeff(T, p);
324 const LhsEval& diffAW = BinaryCoeff::H2O_Air::gasDiffCoeff(T, p);
325
326 const LhsEval& xga = decay<LhsEval>(fluidState.moleFraction(gasPhaseIdx, airIdx));
327 const LhsEval& xgw = decay<LhsEval>(fluidState.moleFraction(gasPhaseIdx, H2OIdx));
328 const LhsEval& xgc = decay<LhsEval>(fluidState.moleFraction(gasPhaseIdx, NAPLIdx));
329
330 if (compIdx==NAPLIdx) return (1 - xgw)/(xga/diffAW + xgc/diffWC);
331 else if (compIdx==H2OIdx) return (1 - xgc)/(xgw/diffWC + xga/diffAC);
332 else if (compIdx==airIdx) throw std::logic_error("Diffusivity of air in the gas phase "
333 "is constraint by sum of diffusive fluxes = 0 !\n");
334 }
335 else if (phaseIdx==waterPhaseIdx){
336 const LhsEval& diffACl = 1.e-9; // BinaryCoeff::Air_Mesitylene::liquidDiffCoeff(temperature, pressure);
337 const LhsEval& diffWCl = 1.e-9; // BinaryCoeff::H2O_Mesitylene::liquidDiffCoeff(temperature, pressure);
338 const LhsEval& diffAWl = 1.e-9; // BinaryCoeff::H2O_Air::liquidDiffCoeff(temperature, pressure);
339
340 const LhsEval& xwa = decay<LhsEval>(fluidState.moleFraction(waterPhaseIdx, airIdx));
341 const LhsEval& xww = decay<LhsEval>(fluidState.moleFraction(waterPhaseIdx, H2OIdx));
342 const LhsEval& xwc = decay<LhsEval>(fluidState.moleFraction(waterPhaseIdx, NAPLIdx));
343
344 switch (compIdx) {
345 case NAPLIdx:
346 diffCont = (1.- xww)/(xwa/diffAWl + xwc/diffWCl);
347 return diffCont;
348 case airIdx:
349 diffCont = (1.- xwc)/(xww/diffWCl + xwa/diffACl);
350 return diffCont;
351 case H2OIdx:
352 throw std::logic_error("Diffusivity of water in the water phase "
353 "is constraint by sum of diffusive fluxes = 0 !\n");
354 };
355 }
356 else if (phaseIdx==naplPhaseIdx) {
357 throw std::logic_error("Diffusion coefficients of "
358 "substances in liquid phase are undefined!\n");
359 }
360 return 0;
361#endif
362 }
363
365 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
366 static LhsEval fugacityCoefficient(const FluidState& fluidState,
367 const ParameterCache<ParamCacheEval>& /*paramCache*/,
368 unsigned phaseIdx,
369 unsigned compIdx)
370 {
371 assert(phaseIdx < numPhases);
372 assert(compIdx < numComponents);
373
374 const LhsEval& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
375 const LhsEval& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
378
379 if (phaseIdx == waterPhaseIdx) {
380 if (compIdx == H2OIdx)
381 return H2O::vaporPressure(T)/p;
382 else if (compIdx == airIdx)
383 return BinaryCoeff::H2O_N2::henry(T)/p;
384 else if (compIdx == NAPLIdx)
386 assert(false);
387 }
388 // for the NAPL phase, we assume currently that nothing is
389 // dissolved. this means that the affinity of the NAPL
390 // component to the NAPL phase is much higher than for the
391 // other components, i.e. the fugacity cofficient is much
392 // smaller.
393 else if (phaseIdx == naplPhaseIdx) {
394 const LhsEval& phiNapl = NAPL::vaporPressure(T)/p;
395 if (compIdx == NAPLIdx)
396 return phiNapl;
397 else if (compIdx == airIdx)
398 return 1e6*phiNapl;
399 else if (compIdx == H2OIdx)
400 return 1e6*phiNapl;
401 assert(false);
402 }
403
404 // for the gas phase, assume an ideal gas when it comes to
405 // fugacity (-> fugacity == partial pressure)
406 assert(phaseIdx == gasPhaseIdx);
407 return 1.0;
408 }
409
410
412 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
413 static LhsEval enthalpy(const FluidState& fluidState,
414 const ParameterCache<ParamCacheEval>& /*paramCache*/,
415 unsigned phaseIdx)
416 {
417 const LhsEval& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
418 const LhsEval& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
419
420 if (phaseIdx == waterPhaseIdx) {
421 return H2O::liquidEnthalpy(T, p);
422 }
423 else if (phaseIdx == naplPhaseIdx) {
424 return NAPL::liquidEnthalpy(T, p);
425 }
426 else if (phaseIdx == gasPhaseIdx) {
427 // gas phase enthalpy depends strongly on composition
428 LhsEval result = 0;
429 result += H2O::gasEnthalpy(T, p) * decay<LhsEval>(fluidState.massFraction(gasPhaseIdx, H2OIdx));
430 result += NAPL::gasEnthalpy(T, p) * decay<LhsEval>(fluidState.massFraction(gasPhaseIdx, airIdx));
431 result += Air::gasEnthalpy(T, p) * decay<LhsEval>(fluidState.massFraction(gasPhaseIdx, NAPLIdx));
432
433 return result;
434 }
435 throw std::logic_error("Invalid phase index "+std::to_string(phaseIdx));
436 }
437
439 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
440 static LhsEval thermalConductivity(const FluidState& fluidState,
441 const ParameterCache<ParamCacheEval>& /*paramCache*/,
442 unsigned phaseIdx)
443 {
444 assert(phaseIdx < numPhases);
445
446 if (phaseIdx == waterPhaseIdx){ // water phase
447 const LhsEval& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
448 const LhsEval& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
449
451 }
452 else if (phaseIdx == gasPhaseIdx) { // gas phase
453 const LhsEval& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
454 const LhsEval& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
455
456 return Air::gasThermalConductivity(T, p);
457 }
458
459 assert(phaseIdx == naplPhaseIdx);
460
461 // Taken from:
462 //
463 // D. K. H. Briggs: "Thermal Conductivity of Liquids",
464 // Ind. Eng. Chem., 1957, 49 (3), pp 418–421
465 //
466 // Convertion to SI units:
467 // 344e-6 cal/(s cm K) = 0.0143964 J/(s m K)
468 return 0.0143964;
469 }
470};
471} // namespace Opm
472
473#endif
static Evaluation gasThermalConductivity(const Evaluation &, const Evaluation &)
Specific heat conductivity of steam .
Definition: Air.hpp:217
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
static Evaluation gasViscosity(const Evaluation &temperature, const Evaluation &)
The dynamic viscosity of at a given pressure and temperature.
Definition: Air.hpp:137
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 mesitylene. I used the method according to Wilke and ...
Definition: Air_Mesitylene.hpp:59
static Evaluation gasDiffCoeff(const Evaluation &temperature, const Evaluation &pressure)
Binary diffusion coefficent for molecular water and air.
Definition: H2O_Air.hpp:70
static Evaluation henry(const Evaluation &)
Henry coefficent for mesitylene in liquid water.
Definition: H2O_Mesitylene.hpp:52
static Evaluation gasDiffCoeff(Evaluation temperature, Evaluation pressure)
Binary diffusion coefficent [m^2/s] for molecular water and mesitylene.
Definition: H2O_Mesitylene.hpp:67
static Evaluation henry(const Evaluation &temperature)
Henry coefficent for molecular nitrogen in liquid water.
Definition: H2O_N2.hpp:52
A fluid system with water, gas and NAPL as phases and water, air and mesitylene (DNAPL) as components...
Definition: H2OAirMesityleneFluidSystem.hpp:57
static const int numComponents
Number of chemical species in the fluid system.
Definition: H2OAirMesityleneFluidSystem.hpp:83
static const int airIdx
The index of the air pseudo-component.
Definition: H2OAirMesityleneFluidSystem.hpp:97
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: H2OAirMesityleneFluidSystem.hpp:366
static const int H2OIdx
The index of the water component.
Definition: H2OAirMesityleneFluidSystem.hpp:93
static const int waterPhaseIdx
The index of the water phase.
Definition: H2OAirMesityleneFluidSystem.hpp:86
static const int naplPhaseIdx
The index of the NAPL phase.
Definition: H2OAirMesityleneFluidSystem.hpp:88
Mesitylene< Scalar > NAPL
The type of the mesithylene/napl component.
Definition: H2OAirMesityleneFluidSystem.hpp:70
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: H2OAirMesityleneFluidSystem.hpp:121
static bool isLiquid(unsigned phaseIdx)
Return whether a phase is liquid.
Definition: H2OAirMesityleneFluidSystem.hpp:131
static LhsEval thermalConductivity(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Thermal conductivity of a fluid phase [W/(m K)].
Definition: H2OAirMesityleneFluidSystem.hpp:440
TabulatedH2O H2O
The type of the water component.
Definition: H2OAirMesityleneFluidSystem.hpp:77
static bool isCompressible(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: H2OAirMesityleneFluidSystem.hpp:142
static const int NAPLIdx
The index of the NAPL component.
Definition: H2OAirMesityleneFluidSystem.hpp:95
static Scalar molarMass(unsigned compIdx)
Return the molar mass of a component in [kg/mol].
Definition: H2OAirMesityleneFluidSystem.hpp:186
static const char * componentName(unsigned compIdx)
Return the human readable name of a component.
Definition: H2OAirMesityleneFluidSystem.hpp:175
static LhsEval density(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Calculate the density [kg/m^3] of a fluid phase.
Definition: H2OAirMesityleneFluidSystem.hpp:201
static const char * phaseName(unsigned phaseIdx)
Return the human readable name of a fluid phase.
Definition: H2OAirMesityleneFluidSystem.hpp:164
static LhsEval viscosity(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Calculate the dynamic viscosity of a fluid phase [Pa*s].
Definition: H2OAirMesityleneFluidSystem.hpp:246
static LhsEval diffusionCoefficient(const FluidState &, const ParameterCache< ParamCacheEval > &, unsigned, unsigned)
Calculate the binary molecular diffusion coefficient for a component in a fluid phase [mol^2 * s / (k...
Definition: H2OAirMesityleneFluidSystem.hpp:308
static const int numPhases
Number of fluid phases in the fluid system.
Definition: H2OAirMesityleneFluidSystem.hpp:81
static const int gasPhaseIdx
The index of the gas phase.
Definition: H2OAirMesityleneFluidSystem.hpp:90
static void init()
Initialize the fluid system's static parameters.
Definition: H2OAirMesityleneFluidSystem.hpp:100
static bool isIdealMixture(unsigned)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: H2OAirMesityleneFluidSystem.hpp:154
::Opm::Air< Scalar > Air
The type of the air component.
Definition: H2OAirMesityleneFluidSystem.hpp:73
static bool isIdealGas(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: H2OAirMesityleneFluidSystem.hpp:138
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: H2OAirMesityleneFluidSystem.hpp:413
Material properties of pure water .
Definition: H2O.hpp:62
Component for Mesitylene.
Definition: Mesitylene.hpp:45
static bool liquidIsCompressible()
Returns true iff the liquid phase is assumed to be compressible.
Definition: Mesitylene.hpp:218
static Evaluation vaporPressure(const Evaluation &temperature)
The saturation vapor pressure in of pure mesitylene at a given temperature according to Antoine afte...
Definition: Mesitylene.hpp:99
static Evaluation liquidEnthalpy(const Evaluation &temperature, const Evaluation &pressure)
Specific enthalpy of liquid mesitylene .
Definition: Mesitylene.hpp:118
static Scalar molarMass()
The molar mass in of mesitylene.
Definition: Mesitylene.hpp:58
static Evaluation gasDensity(const Evaluation &temperature, const Evaluation &pressure)
The density of pure mesitylene vapor at a given pressure and temperature .
Definition: Mesitylene.hpp:190
static const char * name()
A human readable name for the mesitylene.
Definition: Mesitylene.hpp:52
static bool gasIsIdeal()
Returns true iff the gas phase is assumed to be ideal.
Definition: Mesitylene.hpp:212
static Evaluation gasEnthalpy(const Evaluation &temperature, const Evaluation &pressure)
Specific enthalpy of mesitylene vapor .
Definition: Mesitylene.hpp:178
static Evaluation liquidDensity(const Evaluation &temperature, const Evaluation &)
The density of pure mesitylene at a given pressure and temperature .
Definition: Mesitylene.hpp:200
static Evaluation liquidViscosity(Evaluation temperature, const Evaluation &)
The dynamic viscosity of pure mesitylene.
Definition: Mesitylene.hpp:255
static Evaluation gasViscosity(Evaluation temperature, const Evaluation &, bool=true)
The dynamic viscosity of mesitylene vapor.
Definition: Mesitylene.hpp:229
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 gasEnthalpy(const Evaluation &temperature, const Evaluation &pressure)
Specific enthalpy of the gas .
Definition: TabulatedComponent.hpp:282
static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp, Scalar pressMin, Scalar pressMax, unsigned nPress)
Initialize the tables.
Definition: TabulatedComponent.hpp:72
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 bool gasIsIdeal()
Returns true iff the gas phase is assumed to be ideal.
Definition: TabulatedComponent.hpp:414
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
Definition: H2OAirMesityleneFluidSystem.hpp:67
Definition: MathToolbox.hpp:50