3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
h2oairxylene.hh
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program 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 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program 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 this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_H2O_AIR_XYLENE_FLUID_SYSTEM_HH
25#define DUMUX_H2O_AIR_XYLENE_FLUID_SYSTEM_HH
26
32
36
38
39#include <dumux/io/name.hh>
40
41namespace Dumux {
42namespace FluidSystems {
43
52template <class Scalar,
53 class H2OType = Components::TabulatedComponent<Components::H2O<Scalar> > >
55 : public Base<Scalar, H2OAirXylene<Scalar, H2OType> >
56{
59
60public:
61 using H2O = H2OType;
64
65 static const int numPhases = 3;
66 static const int numComponents = 3;
67
68 static const int wPhaseIdx = 0; // index of the water phase
69 static const int nPhaseIdx = 1; // index of the NAPL phase
70 static const int gPhaseIdx = 2; // index of the gas phase
71
72 static const int H2OIdx = 0;
73 static const int NAPLIdx = 1;
74 static const int AirIdx = 2;
75
76 // export component indices to indicate the main component
77 // of the corresponding phase at atmospheric pressure 1 bar
78 // and room temperature 20°C:
79 static const int wCompIdx = H2OIdx;
80 static const int nCompIdx = NAPLIdx;
81 static const int gCompIdx = AirIdx;
82
89 static void init()
90 {
91 init(/*tempMin=*/273.15,
92 /*tempMax=*/623.15,
93 /*numTemp=*/100,
94 /*pMin=*/0.0,
95 /*pMax=*/20e6,
96 /*numP=*/200);
97 }
98
110 static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp,
111 Scalar pressMin, Scalar pressMax, unsigned nPress)
112 {
113 if (H2O::isTabulated)
114 {
115 H2O::init(tempMin, tempMax, nTemp,
116 pressMin, pressMax, nPress);
117 }
118 }
119
123 static constexpr bool isMiscible()
124 { return true; }
125
131 static constexpr bool isGas(int phaseIdx)
132 {
133 assert(0 <= phaseIdx && phaseIdx < numPhases);
134 return phaseIdx == gPhaseIdx;
135 }
136
143 static bool isIdealGas(int phaseIdx)
144 { return phaseIdx == gPhaseIdx && H2O::gasIsIdeal() && Air::gasIsIdeal() && NAPL::gasIsIdeal(); }
145
160 static bool isIdealMixture(int phaseIdx)
161 {
162 assert(0 <= phaseIdx && phaseIdx < numPhases);
163 // we assume Henry's and Raoult's laws for the water phase and
164 // and no interaction between gas molecules of different
165 // components, so all phases are ideal mixtures!
166 return true;
167 }
168
178 static constexpr bool isCompressible(int phaseIdx)
179 {
180 assert(0 <= phaseIdx && phaseIdx < numPhases);
181 // gases are always compressible
182 if (phaseIdx == gPhaseIdx)
183 return true;
184 else if (phaseIdx == wPhaseIdx)
185 // the water component decides for the water phase...
186 return H2O::liquidIsCompressible();
187
188 // the NAPL component decides for the napl phase...
190 }
191
196 static std::string phaseName(int phaseIdx)
197 {
198 assert(0 <= phaseIdx && phaseIdx < numPhases);
199 switch (phaseIdx)
200 {
201 case wPhaseIdx: return IOName::aqueousPhase();
202 case nPhaseIdx: return IOName::naplPhase();
203 case gPhaseIdx: return IOName::gaseousPhase();
204 }
205 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
206 }
207
212 static std::string componentName(int compIdx)
213 {
214 switch (compIdx) {
215 case H2OIdx: return H2O::name();
216 case AirIdx: return Air::name();
217 case NAPLIdx: return NAPL::name();
218 }
219 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx);
220 }
221
226 static Scalar molarMass(int compIdx)
227 {
228 switch (compIdx) {
229 case H2OIdx: return H2O::molarMass();
230 case AirIdx: return Air::molarMass();
231 case NAPLIdx: return NAPL::molarMass();
232 }
233 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx);
234 }
235
236 using Base::density;
249 template <class FluidState>
250 static Scalar density(const FluidState &fluidState, int phaseIdx)
251 {
252 if (phaseIdx == wPhaseIdx) {
253 // This assumes each gas molecule displaces exactly one
254 // molecule in the liquid.
255 return H2O::liquidMolarDensity(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx))
256 * (H2O::molarMass()*fluidState.moleFraction(wPhaseIdx, H2OIdx)
257 + Air::molarMass()*fluidState.moleFraction(wPhaseIdx, AirIdx)
258 + NAPL::molarMass()*fluidState.moleFraction(wPhaseIdx, NAPLIdx));
259 }
260 else if (phaseIdx == nPhaseIdx) {
261 // assume pure NAPL for the NAPL phase
262 Scalar pressure = NAPL::liquidIsCompressible()?fluidState.pressure(phaseIdx):1e100;
263 return NAPL::liquidDensity(fluidState.temperature(phaseIdx), pressure);
264 }
265
266 assert (phaseIdx == gPhaseIdx);
267 Scalar pH2O =
268 fluidState.moleFraction(gPhaseIdx, H2OIdx) *
269 fluidState.pressure(gPhaseIdx);
270 Scalar pAir =
271 fluidState.moleFraction(gPhaseIdx, AirIdx) *
272 fluidState.pressure(gPhaseIdx);
273 Scalar pNAPL =
274 fluidState.moleFraction(gPhaseIdx, NAPLIdx) *
275 fluidState.pressure(gPhaseIdx);
276 return H2O::gasDensity(fluidState.temperature(phaseIdx), pH2O)
277 + Air::gasDensity(fluidState.temperature(phaseIdx), pAir)
278 + NAPL::gasDensity(fluidState.temperature(phaseIdx), pNAPL);
279 }
280
281 using Base::molarDensity;
291 template <class FluidState>
292 static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
293 {
294 Scalar temperature = fluidState.temperature(phaseIdx);
295 Scalar pressure = fluidState.pressure(phaseIdx);
296 if (phaseIdx == nPhaseIdx)
297 {
299 }
300 else if (phaseIdx == wPhaseIdx)
301 { // This assumes each gas molecule displaces exactly one
302 // molecule in the liquid.
303 return H2O::liquidMolarDensity(temperature, pressure);
304 }
305 else
306 {
307 return H2O::gasMolarDensity(temperature, fluidState.partialPressure(gPhaseIdx, H2OIdx))
308 + NAPL::gasMolarDensity(temperature, fluidState.partialPressure(gPhaseIdx, NAPLIdx))
309 + Air::gasMolarDensity(temperature, fluidState.partialPressure(gPhaseIdx, AirIdx));
310 }
311 }
312
313 using Base::viscosity;
321 template <class FluidState>
322 static Scalar viscosity(const FluidState &fluidState,
323 int phaseIdx)
324 {
325 if (phaseIdx == wPhaseIdx) {
326 // assume pure water viscosity
327 return H2O::liquidViscosity(fluidState.temperature(phaseIdx),
328 fluidState.pressure(phaseIdx));
329 }
330 else if (phaseIdx == nPhaseIdx) {
331 // assume pure NAPL viscosity
332 return NAPL::liquidViscosity(fluidState.temperature(phaseIdx),
333 fluidState.pressure(phaseIdx));
334 }
335
336 assert (phaseIdx == gPhaseIdx);
337
338 /* Wilke method. See:
339 *
340 * See: R. Reid, et al.: The Properties of Gases and Liquids,
341 * 4th edition, McGraw-Hill, 1987, 407-410
342 * 5th edition, McGraw-Hill, 20001, p. 9.21/22
343 *
344 * in this case, we use a simplified version in order to avoid
345 * computationally costly evaluation of sqrt and pow functions and
346 * divisions
347 * -- compare e.g. with Promo Class p. 32/33
348 */
349 Scalar muResult;
350 const Scalar mu[numComponents] = {
351 h2oGasViscosityInMixture(fluidState.temperature(phaseIdx),fluidState.pressure(phaseIdx)),
352 Air::gasViscosity(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx)),
353 NAPL::gasViscosity(fluidState.temperature(phaseIdx), NAPL::vaporPressure(fluidState.temperature(phaseIdx)))
354 };
355 // molar masses
356 const Scalar M[numComponents] = {
357 H2O::molarMass(),
360 };
361
362 Scalar muAW = mu[AirIdx]*fluidState.moleFraction(gPhaseIdx, AirIdx)
363 + mu[H2OIdx]*fluidState.moleFraction(gPhaseIdx, H2OIdx)
364 / (fluidState.moleFraction(gPhaseIdx, AirIdx)
365 + fluidState.moleFraction(gPhaseIdx, H2OIdx));
366 Scalar xAW = fluidState.moleFraction(gPhaseIdx, AirIdx)
367 + fluidState.moleFraction(gPhaseIdx, H2OIdx);
368
369 Scalar MAW = (fluidState.moleFraction(gPhaseIdx, AirIdx)*Air::molarMass()
370 + fluidState.moleFraction(gPhaseIdx, H2OIdx)*H2O::molarMass())
371 / xAW;
372
373 Scalar phiCAW = 0.3; // simplification for this particular system
374 /* actually like this
375 * using std::sqrt;
376 * using std::pow;
377 * Scalar phiCAW = pow(1.+sqrt(mu[NAPLIdx]/muAW)*pow(MAW/M[NAPLIdx],0.25),2)
378 * / sqrt(8.*(1.+M[NAPLIdx]/MAW));
379 */
380 Scalar phiAWC = phiCAW * muAW*M[NAPLIdx]/(mu[NAPLIdx]*MAW);
381
382 muResult = (xAW*muAW)/(xAW+fluidState.moleFraction(gPhaseIdx, NAPLIdx)*phiAWC)
383 + (fluidState.moleFraction(gPhaseIdx, NAPLIdx) * mu[NAPLIdx])
384 / (fluidState.moleFraction(gPhaseIdx, NAPLIdx) + xAW*phiCAW);
385 return muResult;
386 }
387
388
397 template <class FluidState>
398 static Scalar diffusionCoefficient(const FluidState &fluidState,
399 int phaseIdx,
400 int compIdx)
401 {
402 Scalar diffCont;
403
404 Scalar temperature = fluidState.temperature(phaseIdx);
405 Scalar pressure = fluidState.pressure(phaseIdx);
406 if (phaseIdx==gPhaseIdx) {
410
411 const Scalar xga = fluidState.moleFraction(gPhaseIdx, AirIdx);
412 const Scalar xgw = fluidState.moleFraction(gPhaseIdx, H2OIdx);
413 const Scalar xgc = fluidState.moleFraction(gPhaseIdx, NAPLIdx);
414
415 if (compIdx==NAPLIdx)
416 return (1.- xgw)/(xga/diffAW + xgc/diffWC);
417 else if (compIdx==H2OIdx)
418 return (1.- xgc)/(xgw/diffWC + xga/diffAC);
419 else if (compIdx==AirIdx)
420 DUNE_THROW(Dune::InvalidStateException,
421 "Diffusivity of Air in the gas phase "
422 "is constraint by sum of diffusive fluxes = 0 !\n");
423 } else if (phaseIdx==wPhaseIdx){
427
428 Scalar xwa = fluidState.moleFraction(wPhaseIdx, AirIdx);
429 Scalar xww = fluidState.moleFraction(wPhaseIdx, H2OIdx);
430 Scalar xwc = fluidState.moleFraction(wPhaseIdx, NAPLIdx);
431
432 switch (compIdx) {
433 case NAPLIdx:
434 diffCont = (1.- xww)/(xwa/diffAWl + xwc/diffWCl);
435 return diffCont;
436 case AirIdx:
437 diffCont = (1.- xwc)/(xww/diffWCl + xwa/diffACl);
438 return diffCont;
439 case H2OIdx:
440 DUNE_THROW(Dune::InvalidStateException,
441 "Diffusivity of water in the water phase "
442 "is constraint by sum of diffusive fluxes = 0 !\n");
443 }
444 } else if (phaseIdx==nPhaseIdx) {
445 return 0.0;
446 }
447 return 0;
448 }
449
451 template <class FluidState>
452 static Scalar binaryDiffusionCoefficient(const FluidState &fluidState,
453 int phaseIdx,
454 int compIIdx,
455 int compJIdx)
456 {
457 DUNE_THROW(Dune::NotImplemented, "FluidSystems::H2OAirXylene::binaryDiffusionCoefficient()");
458 }
459
474 template <class FluidState>
475 static Scalar fugacityCoefficient(const FluidState &fluidState,
476 int phaseIdx,
477 int compIdx)
478 {
479 assert(0 <= phaseIdx && phaseIdx < numPhases);
480 assert(0 <= compIdx && compIdx < numComponents);
481
482 Scalar T = fluidState.temperature(phaseIdx);
483 Scalar p = fluidState.pressure(phaseIdx);
484
485 if (phaseIdx == wPhaseIdx) {
486 if (compIdx == H2OIdx)
487 return H2O::vaporPressure(T)/p;
488 else if (compIdx == AirIdx)
489 return BinaryCoeff::H2O_Air::henry(T)/p;
490 else if (compIdx == NAPLIdx)
492 }
493
494 // for the NAPL phase, we assume currently that nothing is
495 // dissolved. this means that the affinity of the NAPL
496 // component to the NAPL phase is much higher than for the
497 // other components, i.e. the fugacity coefficient is much
498 // smaller.
499 if (phaseIdx == nPhaseIdx) {
500 Scalar phiNapl = NAPL::vaporPressure(T)/p;
501 if (compIdx == NAPLIdx)
502 return phiNapl;
503 else if (compIdx == AirIdx)
504 return 1e6*phiNapl;
505 else if (compIdx == H2OIdx)
506 return 1e6*phiNapl;
507 }
508
509 // for the gas phase, assume an ideal gas when it comes to
510 // fugacity (-> fugacity == partial pressure)
511 assert(phaseIdx == gPhaseIdx);
512 return 1.0;
513 }
514
515 template <class FluidState>
516 static Scalar kelvinVaporPressure(const FluidState &fluidState,
517 const int phaseIdx,
518 const int compIdx)
519 {
520 DUNE_THROW(Dune::NotImplemented, "FluidSystems::H2OAirXylene::kelvinVaporPressure()");
521 }
522
523 using Base::enthalpy;
533 template <class FluidState>
534 static Scalar enthalpy(const FluidState &fluidState,
535 int phaseIdx)
536 {
537 if (phaseIdx == wPhaseIdx) {
538 return H2O::liquidEnthalpy(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
539 }
540 else if (phaseIdx == nPhaseIdx) {
541 return NAPL::liquidEnthalpy(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
542 }
543 else if (phaseIdx == gPhaseIdx) { // gas phase enthalpy depends strongly on composition
544 Scalar hgc = NAPL::gasEnthalpy(fluidState.temperature(phaseIdx),
545 fluidState.pressure(phaseIdx));
546 Scalar hgw = H2O::gasEnthalpy(fluidState.temperature(phaseIdx),
547 fluidState.pressure(phaseIdx));
548 // pressure is only a dummy here (not dependent on pressure, just temperature)
549 Scalar hga = Air::gasEnthalpy(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
550
551 Scalar result = 0;
552 result += hgw * fluidState.massFraction(gPhaseIdx, H2OIdx);
553 result += hga * fluidState.massFraction(gPhaseIdx, AirIdx);
554 result += hgc * fluidState.massFraction(gPhaseIdx, NAPLIdx);
555
556 return result;
557 }
558 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
559 }
560
567 template <class FluidState>
568 static Scalar componentEnthalpy(const FluidState& fluidState, int phaseIdx, int componentIdx)
569 {
570 const Scalar T = fluidState.temperature(phaseIdx);
571 const Scalar p = fluidState.pressure(phaseIdx);
572
573 if (phaseIdx == wPhaseIdx)
574 {
575 if (componentIdx == H2OIdx)
576 return H2O::liquidEnthalpy(T, p);
577 else if (componentIdx == NAPLIdx)
578 DUNE_THROW(Dune::NotImplemented, "The component enthalpy for NAPL in water is not implemented.");
579 else if (componentIdx == AirIdx)
580 DUNE_THROW(Dune::NotImplemented, "The component enthalpy for Air in water is not implemented.");
581 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx);
582 }
583 else if (phaseIdx == nPhaseIdx)
584 {
585 if (componentIdx == H2OIdx)
586 DUNE_THROW(Dune::NotImplemented, "The component enthalpy for water in NAPL is not implemented.");
587 else if (componentIdx == NAPLIdx)
588 return NAPL::liquidEnthalpy(T, p);
589 else if (componentIdx == AirIdx)
590 DUNE_THROW(Dune::NotImplemented, "The component enthalpy for air in NAPL is not implemented.");
591 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx);
592 }
593 else if (phaseIdx == gPhaseIdx)
594 {
595 if (componentIdx == H2OIdx)
596 return H2O::gasEnthalpy(T, p);
597 else if (componentIdx == NAPLIdx)
598 return NAPL::gasEnthalpy(T, p);
599 else if (componentIdx == AirIdx)
600 return Air::gasEnthalpy(T,p);
601 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx);
602 }
603 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
604 }
605
606 using Base::heatCapacity;
612 template <class FluidState>
613 static Scalar heatCapacity(const FluidState &fluidState,
614 int phaseIdx)
615 {
616 DUNE_THROW(Dune::NotImplemented, "FluidSystems::H2OAirXylene::heatCapacity()");
617 }
618
625 template <class FluidState>
626 static Scalar thermalConductivity(const FluidState &fluidState,
627 int phaseIdx)
628 {
629 const Scalar temperature = fluidState.temperature(phaseIdx) ;
630 const Scalar pressure = fluidState.pressure(phaseIdx);
631 if (phaseIdx == wPhaseIdx)
632 {
633 return H2O::liquidThermalConductivity(temperature, pressure);
634 }
635 else if (phaseIdx == nPhaseIdx)
636 {
638 }
639 else if (phaseIdx == gPhaseIdx)
640 {
642 }
643 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
644 }
645
646private:
647 static Scalar waterPhaseDensity_(Scalar T, Scalar pw, Scalar xww, Scalar xwa, Scalar xwc)
648 {
649 Scalar rholH2O = H2O::liquidDensity(T, pw);
650 Scalar clH2O = rholH2O/H2O::molarMass();
651
652 // this assumes each dissolved molecule displaces exactly one
653 // water molecule in the liquid
654 return
655 clH2O*(xww*H2O::molarMass() + xwa*Air::molarMass() + xwc*NAPL::molarMass());
656 }
657
658 static Scalar gasPhaseDensity_(Scalar T, Scalar pg, Scalar xgw, Scalar xga, Scalar xgc)
659 {
660 return H2O::gasDensity(T, pg*xgw) + Air::gasDensity(T, pg*xga) + NAPL::gasDensity(T, pg*xgc);
661 }
662
663 static Scalar NAPLPhaseDensity_(Scalar T, Scalar pn)
664 {
665 return
666 NAPL::liquidDensity(T, pn);
667 }
668
669};
670
671} // end namespace FluidSystems
672} // end namespace Dumux
673
674#endif
A collection of input/output field names for common physical quantities.
Binary coefficients for water and air.
Binary coefficients for air and xylene.
Binary coefficients for water and xylene.
Material properties of pure water .
A simple class for the air fluid properties.
Tabulates all thermodynamic properties of a given untabulated chemical species.
Properties of xylene.
Relations valid for an ideal gas.
Definition: adapt.hh:29
std::string temperature() noexcept
I/O name of temperature for equilibrium models.
Definition: name.hh:51
std::string gaseousPhase() noexcept
I/O name of gaseous phase.
Definition: name.hh:123
std::string naplPhase() noexcept
I/O name of napl phase.
Definition: name.hh:131
std::string aqueousPhase() noexcept
I/O name of aqueous phase.
Definition: name.hh:127
std::string pressure(int phaseIdx) noexcept
I/O name of pressure for multiphase systems.
Definition: name.hh:34
Scalar h2oGasViscosityInMixture(Scalar temperature, Scalar pressure)
The dynamic viscosity of steam in a gas mixture.
Definition: h2o.hh:977
static Scalar gasDiffCoeff(Scalar temperature, Scalar pressure)
Binary diffusion coefficient for air and xylene. method according to Wilke and Lee see W....
Definition: air_xylene.hh:62
static Scalar liquidDiffCoeff(Scalar temperature, Scalar pressure)
Diffusion coefficient for air and xylene in liquid water.
Definition: air_xylene.hh:108
static Scalar henry(Scalar temperature)
Henry coefficient for air in liquid water.
Definition: h2o_air.hh:48
static Scalar gasDiffCoeff(Scalar temperature, Scalar pressure)
Binary diffusion coefficient for molecular water and air.
Definition: h2o_air.hh:68
static Scalar liquidDiffCoeff(Scalar temperature, Scalar pressure)
Diffusion coefficient for molecular nitrogen in liquid water.
Definition: h2o_air.hh:101
static Scalar liquidDiffCoeff(Scalar temperature, Scalar pressure)
Diffusion coefficient for xylene in liquid water.
Definition: h2o_xylene.hh:117
static Scalar henry(Scalar temperature)
Henry coefficient for xylene in liquid water.
Definition: h2o_xylene.hh:52
static Scalar gasDiffCoeff(Scalar temperature, Scalar pressure)
Binary diffusion coefficient for molecular water and xylene.
Definition: h2o_xylene.hh:69
A class for the air fluid properties.
Definition: air.hh:47
static Scalar gasDensity(Scalar temperature, Scalar pressure)
The density of Air at a given pressure and temperature.
Definition: air.hh:85
static constexpr Scalar molarMass()
The molar mass in of Air.
Definition: air.hh:62
static Scalar gasViscosity(Scalar temperature, Scalar pressure)
The dynamic viscosity of Air at a given pressure and temperature.
Definition: air.hh:193
static Scalar gasThermalConductivity(Scalar temperature, Scalar pressure)
Thermal conductivity of air.
Definition: air.hh:351
static constexpr bool gasIsIdeal()
Returns true, the gas phase is assumed to be ideal.
Definition: air.hh:109
static Scalar gasEnthalpy(Scalar temperature, Scalar pressure)
Specific enthalpy of Air with 273.15 as basis.
Definition: air.hh:276
static std::string name()
A human readable name for Air.
Definition: air.hh:54
static Scalar gasMolarDensity(Scalar temperature, Scalar pressure)
The molar density of air in , depending on pressure and temperature.
Definition: air.hh:97
Properties of xylene.
Definition: xylene.hh:49
static Scalar gasViscosity(Scalar temp, Scalar pressure)
The dynamic viscosity of xylene vapor.
Definition: xylene.hh:328
static std::string name()
A human readable name for the xylene.
Definition: xylene.hh:57
static Scalar liquidDensity(Scalar temperature, Scalar pressure)
The density of pure xylene at a given pressure and temperature .
Definition: xylene.hh:299
static Scalar liquidEnthalpy(const Scalar temperature, const Scalar pressure)
Specific enthalpy of liquid xylene .
Definition: xylene.hh:167
static Scalar gasMolarDensity(Scalar temperature, Scalar pressure)
The molar gas density of xylene gas at a given pressure and temperature.
Definition: xylene.hh:262
static constexpr bool gasIsIdeal()
Returns true if the gas phase is assumed to be ideal.
Definition: xylene.hh:313
static Scalar vaporPressure(Scalar temperature)
The saturation vapor pressure in of pure xylene at a given temperature according to Antoine after Be...
Definition: xylene.hh:106
static Scalar liquidThermalConductivity(Scalar temperature, Scalar pressure)
Thermal conductivity of xylene.
Definition: xylene.hh:384
static constexpr bool liquidIsCompressible()
Returns true if the liquid phase is assumed to be compressible.
Definition: xylene.hh:319
static Scalar gasEnthalpy(Scalar temperature, Scalar pressure)
Specific enthalpy of xylene vapor .
Definition: xylene.hh:238
static constexpr Scalar molarMass()
The molar mass in of xylene.
Definition: xylene.hh:63
static Scalar gasDensity(Scalar temperature, Scalar pressure)
The density of xylene gas at a given pressure and temperature.
Definition: xylene.hh:249
static Scalar liquidMolarDensity(Scalar temp, Scalar pressure)
The molar liquid density of pure xylene at a given pressure and temperature .
Definition: xylene.hh:274
static Scalar liquidViscosity(Scalar temp, Scalar pressure)
The dynamic viscosity of pure xylene.
Definition: xylene.hh:356
Fluid system base class.
Definition: fluidsystems/base.hh:45
static Scalar density(const FluidState &fluidState, int phaseIdx)
Calculate the density of a fluid phase.
Definition: fluidsystems/base.hh:134
static Scalar thermalConductivity(const FluidState &fluidState, int phaseIdx)
Thermal conductivity of a fluid phase .
Definition: fluidsystems/base.hh:390
static Scalar fugacityCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Calculate the fugacity coefficient of an individual component in a fluid phase.
Definition: fluidsystems/base.hh:197
static Scalar diffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Calculate the binary molecular diffusion coefficient for a component in a fluid phase .
Definition: fluidsystems/base.hh:278
static Scalar binaryDiffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIIdx, int compJIdx)
Given a phase's composition, temperature and pressure, return the binary diffusion coefficient for c...
Definition: fluidsystems/base.hh:326
static Scalar enthalpy(const FluidState &fluidState, int phaseIdx)
Given a phase's composition, temperature, pressure and density, calculate its specific enthalpy .
Definition: fluidsystems/base.hh:363
static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
Calculate the molar density of a fluid phase.
Definition: fluidsystems/base.hh:160
static Scalar viscosity(const FluidState &fluidState, int phaseIdx)
Calculate the dynamic viscosity of a fluid phase .
Definition: fluidsystems/base.hh:236
static Scalar heatCapacity(const FluidState &fluidState, int phaseIdx)
Specific isobaric heat capacity of a fluid phase .
Definition: fluidsystems/base.hh:424
A three-phase fluid system featuring gas, NAPL and water as phases and distilled water and air (Pseu...
Definition: h2oairxylene.hh:56
static Scalar diffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Given all mole fractions, return the diffusion coefficient of a component in a phase.
Definition: h2oairxylene.hh:398
static Scalar density(const FluidState &fluidState, int phaseIdx)
Given a phase's composition, temperature, pressure, and the partial pressures of all components,...
Definition: h2oairxylene.hh:250
static constexpr bool isGas(int phaseIdx)
Return whether a phase is gaseous.
Definition: h2oairxylene.hh:131
static const int gCompIdx
Definition: h2oairxylene.hh:81
static Scalar viscosity(const FluidState &fluidState, int phaseIdx)
Return the viscosity of a phase .
Definition: h2oairxylene.hh:322
static void init()
Initialize the fluid system's static parameters generically.
Definition: h2oairxylene.hh:89
static Scalar heatCapacity(const FluidState &fluidState, int phaseIdx)
Return the heat capacity in .
Definition: h2oairxylene.hh:613
static const int nCompIdx
Definition: h2oairxylene.hh:80
static Scalar molarMass(int compIdx)
Return the molar mass of a component in .
Definition: h2oairxylene.hh:226
static const int AirIdx
Definition: h2oairxylene.hh:74
static const int gPhaseIdx
Definition: h2oairxylene.hh:70
static std::string componentName(int compIdx)
Return the human readable name of a component (used in indices)
Definition: h2oairxylene.hh:212
static std::string phaseName(int phaseIdx)
Return the human readable name of a phase (used in indices)
Definition: h2oairxylene.hh:196
static const int numComponents
Definition: h2oairxylene.hh:66
static Scalar binaryDiffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIIdx, int compJIdx)
Definition: h2oairxylene.hh:452
static const int H2OIdx
Definition: h2oairxylene.hh:72
static const int wCompIdx
Definition: h2oairxylene.hh:79
static constexpr bool isCompressible(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: h2oairxylene.hh:178
static constexpr bool isMiscible()
Returns whether the fluids are miscible.
Definition: h2oairxylene.hh:123
static bool isIdealGas(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: h2oairxylene.hh:143
H2OType H2O
Definition: h2oairxylene.hh:61
static const int nPhaseIdx
Definition: h2oairxylene.hh:69
static const int numPhases
Definition: h2oairxylene.hh:65
static const int NAPLIdx
Definition: h2oairxylene.hh:73
static Scalar fugacityCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Returns the fugacity coefficient of a component in a phase.
Definition: h2oairxylene.hh:475
static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
The molar density of a fluid phase in .
Definition: h2oairxylene.hh:292
static Scalar componentEnthalpy(const FluidState &fluidState, int phaseIdx, int componentIdx)
Returns the specific enthalpy of a component in a specific phase.
Definition: h2oairxylene.hh:568
static Scalar enthalpy(const FluidState &fluidState, int phaseIdx)
Given all mole fractions in a phase, return the specific phase enthalpy .
Definition: h2oairxylene.hh:534
static const int wPhaseIdx
Definition: h2oairxylene.hh:68
static Scalar thermalConductivity(const FluidState &fluidState, int phaseIdx)
Return the thermal conductivity .
Definition: h2oairxylene.hh:626
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: h2oairxylene.hh:110
static bool isIdealMixture(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: h2oairxylene.hh:160
static Scalar kelvinVaporPressure(const FluidState &fluidState, const int phaseIdx, const int compIdx)
Definition: h2oairxylene.hh:516
Fluid system base class.