version 3.8
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_H2O_AIR_XYLENE_FLUID_SYSTEM_HH
13#define DUMUX_H2O_AIR_XYLENE_FLUID_SYSTEM_HH
14
20
24
26
27#include <dumux/io/name.hh>
28
29namespace Dumux {
30namespace FluidSystems {
31
40template <class Scalar,
41 class H2OType = Components::TabulatedComponent<Components::H2O<Scalar> > >
43 : public Base<Scalar, H2OAirXylene<Scalar, H2OType> >
44{
46
47public:
48 using H2O = H2OType;
51
52 static const int numPhases = 3;
53 static const int numComponents = 3;
54
55 static const int wPhaseIdx = 0; // index of the water phase
56 static const int nPhaseIdx = 1; // index of the NAPL phase
57 static const int gPhaseIdx = 2; // index of the gas phase
58
59 static const int H2OIdx = 0;
60 static const int NAPLIdx = 1;
61 static const int AirIdx = 2;
62
63 // export component indices to indicate the main component
64 // of the corresponding phase at atmospheric pressure 1 bar
65 // and room temperature 20°C:
66 static const int wCompIdx = H2OIdx;
67 static const int nCompIdx = NAPLIdx;
68 static const int gCompIdx = AirIdx;
69
76 static void init()
77 {
78 init(/*tempMin=*/273.15,
79 /*tempMax=*/623.15,
80 /*numTemp=*/100,
81 /*pMin=*/0.0,
82 /*pMax=*/20e6,
83 /*numP=*/200);
84 }
85
97 static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp,
98 Scalar pressMin, Scalar pressMax, unsigned nPress)
99 {
100 if (H2O::isTabulated)
101 {
102 H2O::init(tempMin, tempMax, nTemp,
103 pressMin, pressMax, nPress);
104 }
105 }
106
110 static constexpr bool isMiscible()
111 { return true; }
112
118 static constexpr bool isGas(int phaseIdx)
119 {
120 assert(0 <= phaseIdx && phaseIdx < numPhases);
121 return phaseIdx == gPhaseIdx;
122 }
123
130 static bool isIdealGas(int phaseIdx)
131 { return phaseIdx == gPhaseIdx && H2O::gasIsIdeal() && Air::gasIsIdeal() && NAPL::gasIsIdeal(); }
132
147 static bool isIdealMixture(int phaseIdx)
148 {
149 assert(0 <= phaseIdx && phaseIdx < numPhases);
150 // we assume Henry's and Raoult's laws for the water phase and
151 // and no interaction between gas molecules of different
152 // components, so all phases are ideal mixtures!
153 return true;
154 }
155
165 static constexpr bool isCompressible(int phaseIdx)
166 {
167 assert(0 <= phaseIdx && phaseIdx < numPhases);
168 // gases are always compressible
169 if (phaseIdx == gPhaseIdx)
170 return true;
171 else if (phaseIdx == wPhaseIdx)
172 // the water component decides for the water phase...
173 return H2O::liquidIsCompressible();
174
175 // the NAPL component decides for the napl phase...
177 }
178
183 static std::string phaseName(int phaseIdx)
184 {
185 assert(0 <= phaseIdx && phaseIdx < numPhases);
186 switch (phaseIdx)
187 {
188 case wPhaseIdx: return IOName::aqueousPhase();
189 case nPhaseIdx: return IOName::naplPhase();
190 case gPhaseIdx: return IOName::gaseousPhase();
191 }
192 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
193 }
194
199 static std::string componentName(int compIdx)
200 {
201 switch (compIdx) {
202 case H2OIdx: return H2O::name();
203 case AirIdx: return Air::name();
204 case NAPLIdx: return NAPL::name();
205 }
206 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx);
207 }
208
213 static Scalar molarMass(int compIdx)
214 {
215 switch (compIdx) {
216 case H2OIdx: return H2O::molarMass();
217 case AirIdx: return Air::molarMass();
218 case NAPLIdx: return NAPL::molarMass();
219 }
220 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx);
221 }
222
236 template <class FluidState>
237 static Scalar density(const FluidState &fluidState, int phaseIdx)
238 {
239 if (phaseIdx == wPhaseIdx) {
240 // This assumes each gas molecule displaces exactly one
241 // molecule in the liquid.
242 return H2O::liquidMolarDensity(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx))
243 * (H2O::molarMass()*fluidState.moleFraction(wPhaseIdx, H2OIdx)
244 + Air::molarMass()*fluidState.moleFraction(wPhaseIdx, AirIdx)
245 + NAPL::molarMass()*fluidState.moleFraction(wPhaseIdx, NAPLIdx));
246 }
247 else if (phaseIdx == nPhaseIdx) {
248 // assume pure NAPL for the NAPL phase
249 Scalar pressure = NAPL::liquidIsCompressible()?fluidState.pressure(phaseIdx):1e100;
250 return NAPL::liquidDensity(fluidState.temperature(phaseIdx), pressure);
251 }
252
253 assert (phaseIdx == gPhaseIdx);
254 Scalar pH2O =
255 fluidState.moleFraction(gPhaseIdx, H2OIdx) *
256 fluidState.pressure(gPhaseIdx);
257 Scalar pAir =
258 fluidState.moleFraction(gPhaseIdx, AirIdx) *
259 fluidState.pressure(gPhaseIdx);
260 Scalar pNAPL =
261 fluidState.moleFraction(gPhaseIdx, NAPLIdx) *
262 fluidState.pressure(gPhaseIdx);
263 return H2O::gasDensity(fluidState.temperature(phaseIdx), pH2O)
264 + Air::gasDensity(fluidState.temperature(phaseIdx), pAir)
265 + NAPL::gasDensity(fluidState.temperature(phaseIdx), pNAPL);
266 }
267
270 template <class FluidState>
271 static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
272 {
273 Scalar temperature = fluidState.temperature(phaseIdx);
274 Scalar pressure = fluidState.pressure(phaseIdx);
275 if (phaseIdx == nPhaseIdx)
276 {
278 }
279 else if (phaseIdx == wPhaseIdx)
280 { // This assumes each gas molecule displaces exactly one
281 // molecule in the liquid.
282 return H2O::liquidMolarDensity(temperature, pressure);
283 }
284 else
285 {
286 return H2O::gasMolarDensity(temperature, fluidState.partialPressure(gPhaseIdx, H2OIdx))
287 + NAPL::gasMolarDensity(temperature, fluidState.partialPressure(gPhaseIdx, NAPLIdx))
288 + Air::gasMolarDensity(temperature, fluidState.partialPressure(gPhaseIdx, AirIdx));
289 }
290 }
291
300 template <class FluidState>
301 static Scalar viscosity(const FluidState &fluidState,
302 int phaseIdx)
303 {
304 if (phaseIdx == wPhaseIdx) {
305 // assume pure water viscosity
306 return H2O::liquidViscosity(fluidState.temperature(phaseIdx),
307 fluidState.pressure(phaseIdx));
308 }
309 else if (phaseIdx == nPhaseIdx) {
310 // assume pure NAPL viscosity
311 return NAPL::liquidViscosity(fluidState.temperature(phaseIdx),
312 fluidState.pressure(phaseIdx));
313 }
314
315 assert (phaseIdx == gPhaseIdx);
316
317 /* Wilke method. See:
318 *
319 * See: R. Reid, et al.: The Properties of Gases and Liquids,
320 * 4th edition, McGraw-Hill, 1987, 407-410
321 * 5th edition, McGraw-Hill, 20001, p. 9.21/22
322 *
323 * in this case, we use a simplified version in order to avoid
324 * computationally costly evaluation of sqrt and pow functions and
325 * divisions
326 * -- compare e.g. with Promo Class p. 32/33
327 */
328 Scalar muResult;
329 const Scalar mu[numComponents] = {
330 h2oGasViscosityInMixture(fluidState.temperature(phaseIdx),fluidState.pressure(phaseIdx)),
331 Air::gasViscosity(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx)),
332 NAPL::gasViscosity(fluidState.temperature(phaseIdx), NAPL::vaporPressure(fluidState.temperature(phaseIdx)))
333 };
334 // molar masses
335 const Scalar M[numComponents] = {
336 H2O::molarMass(),
339 };
340
341 Scalar muAW = mu[AirIdx]*fluidState.moleFraction(gPhaseIdx, AirIdx)
342 + mu[H2OIdx]*fluidState.moleFraction(gPhaseIdx, H2OIdx)
343 / (fluidState.moleFraction(gPhaseIdx, AirIdx)
344 + fluidState.moleFraction(gPhaseIdx, H2OIdx));
345 Scalar xAW = fluidState.moleFraction(gPhaseIdx, AirIdx)
346 + fluidState.moleFraction(gPhaseIdx, H2OIdx);
347
348 Scalar MAW = (fluidState.moleFraction(gPhaseIdx, AirIdx)*Air::molarMass()
349 + fluidState.moleFraction(gPhaseIdx, H2OIdx)*H2O::molarMass())
350 / xAW;
351
352 Scalar phiCAW = 0.3; // simplification for this particular system
353 /* actually like this
354 * using std::sqrt;
355 * using std::pow;
356 * Scalar phiCAW = pow(1.+sqrt(mu[NAPLIdx]/muAW)*pow(MAW/M[NAPLIdx],0.25),2)
357 * / sqrt(8.*(1.+M[NAPLIdx]/MAW));
358 */
359 Scalar phiAWC = phiCAW * muAW*M[NAPLIdx]/(mu[NAPLIdx]*MAW);
360
361 muResult = (xAW*muAW)/(xAW+fluidState.moleFraction(gPhaseIdx, NAPLIdx)*phiAWC)
362 + (fluidState.moleFraction(gPhaseIdx, NAPLIdx) * mu[NAPLIdx])
363 / (fluidState.moleFraction(gPhaseIdx, NAPLIdx) + xAW*phiCAW);
364 return muResult;
365 }
366
367
370 template <class FluidState>
371 static Scalar diffusionCoefficient(const FluidState &fluidState,
372 int phaseIdx,
373 int compIdx)
374 {
375 Scalar diffCont;
376
377 Scalar temperature = fluidState.temperature(phaseIdx);
378 Scalar pressure = fluidState.pressure(phaseIdx);
379 if (phaseIdx==gPhaseIdx) {
383
384 const Scalar xga = fluidState.moleFraction(gPhaseIdx, AirIdx);
385 const Scalar xgw = fluidState.moleFraction(gPhaseIdx, H2OIdx);
386 const Scalar xgc = fluidState.moleFraction(gPhaseIdx, NAPLIdx);
387
388 if (compIdx==NAPLIdx)
389 return (1.- xgw)/(xga/diffAW + xgc/diffWC);
390 else if (compIdx==H2OIdx)
391 return (1.- xgc)/(xgw/diffWC + xga/diffAC);
392 else if (compIdx==AirIdx)
393 DUNE_THROW(Dune::InvalidStateException,
394 "Diffusivity of Air in the gas phase "
395 "is constraint by sum of diffusive fluxes = 0 !\n");
396 } else if (phaseIdx==wPhaseIdx){
400
401 Scalar xwa = fluidState.moleFraction(wPhaseIdx, AirIdx);
402 Scalar xww = fluidState.moleFraction(wPhaseIdx, H2OIdx);
403 Scalar xwc = fluidState.moleFraction(wPhaseIdx, NAPLIdx);
404
405 switch (compIdx) {
406 case NAPLIdx:
407 diffCont = (1.- xww)/(xwa/diffAWl + xwc/diffWCl);
408 return diffCont;
409 case AirIdx:
410 diffCont = (1.- xwc)/(xww/diffWCl + xwa/diffACl);
411 return diffCont;
412 case H2OIdx:
413 DUNE_THROW(Dune::InvalidStateException,
414 "Diffusivity of water in the water phase "
415 "is constraint by sum of diffusive fluxes = 0 !\n");
416 }
417 } else if (phaseIdx==nPhaseIdx) {
418 return 0.0;
419 }
420 return 0;
421 }
422
425 template <class FluidState>
426 static Scalar binaryDiffusionCoefficient(const FluidState &fluidState,
427 int phaseIdx,
428 int compIIdx,
429 int compJIdx)
430 {
431 DUNE_THROW(Dune::NotImplemented, "FluidSystems::H2OAirXylene::binaryDiffusionCoefficient()");
432 }
433
448 template <class FluidState>
449 static Scalar fugacityCoefficient(const FluidState &fluidState,
450 int phaseIdx,
451 int compIdx)
452 {
453 assert(0 <= phaseIdx && phaseIdx < numPhases);
454 assert(0 <= compIdx && compIdx < numComponents);
455
456 Scalar T = fluidState.temperature(phaseIdx);
457 Scalar p = fluidState.pressure(phaseIdx);
458
459 if (phaseIdx == wPhaseIdx) {
460 if (compIdx == H2OIdx)
461 return H2O::vaporPressure(T)/p;
462 else if (compIdx == AirIdx)
463 return BinaryCoeff::H2O_Air::henry(T)/p;
464 else if (compIdx == NAPLIdx)
466 }
467
468 // for the NAPL phase, we assume currently that nothing is
469 // dissolved. this means that the affinity of the NAPL
470 // component to the NAPL phase is much higher than for the
471 // other components, i.e. the fugacity coefficient is much
472 // smaller.
473 if (phaseIdx == nPhaseIdx) {
474 Scalar phiNapl = NAPL::vaporPressure(T)/p;
475 if (compIdx == NAPLIdx)
476 return phiNapl;
477 else if (compIdx == AirIdx)
478 return 1e6*phiNapl;
479 else if (compIdx == H2OIdx)
480 return 1e6*phiNapl;
481 }
482
483 // for the gas phase, assume an ideal gas when it comes to
484 // fugacity (-> fugacity == partial pressure)
485 assert(phaseIdx == gPhaseIdx);
486 return 1.0;
487 }
488
489 template <class FluidState>
490 static Scalar kelvinVaporPressure(const FluidState &fluidState,
491 const int phaseIdx,
492 const int compIdx)
493 {
494 DUNE_THROW(Dune::NotImplemented, "FluidSystems::H2OAirXylene::kelvinVaporPressure()");
495 }
496
507 template <class FluidState>
508 static Scalar enthalpy(const FluidState &fluidState,
509 int phaseIdx)
510 {
511 if (phaseIdx == wPhaseIdx) {
512 return H2O::liquidEnthalpy(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
513 }
514 else if (phaseIdx == nPhaseIdx) {
515 return NAPL::liquidEnthalpy(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
516 }
517 else if (phaseIdx == gPhaseIdx) { // gas phase enthalpy depends strongly on composition
518 Scalar hgc = NAPL::gasEnthalpy(fluidState.temperature(phaseIdx),
519 fluidState.pressure(phaseIdx));
520 Scalar hgw = H2O::gasEnthalpy(fluidState.temperature(phaseIdx),
521 fluidState.pressure(phaseIdx));
522 // pressure is only a dummy here (not dependent on pressure, just temperature)
523 Scalar hga = Air::gasEnthalpy(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
524
525 Scalar result = 0;
526 result += hgw * fluidState.massFraction(gPhaseIdx, H2OIdx);
527 result += hga * fluidState.massFraction(gPhaseIdx, AirIdx);
528 result += hgc * fluidState.massFraction(gPhaseIdx, NAPLIdx);
529
530 return result;
531 }
532 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
533 }
534
541 template <class FluidState>
542 static Scalar componentEnthalpy(const FluidState& fluidState, int phaseIdx, int componentIdx)
543 {
544 const Scalar T = fluidState.temperature(phaseIdx);
545 const Scalar p = fluidState.pressure(phaseIdx);
546
547 if (phaseIdx == wPhaseIdx)
548 {
549 if (componentIdx == H2OIdx)
550 return H2O::liquidEnthalpy(T, p);
551 else if (componentIdx == NAPLIdx)
552 DUNE_THROW(Dune::NotImplemented, "The component enthalpy for NAPL in water is not implemented.");
553 else if (componentIdx == AirIdx)
554 DUNE_THROW(Dune::NotImplemented, "The component enthalpy for Air in water is not implemented.");
555 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx);
556 }
557 else if (phaseIdx == nPhaseIdx)
558 {
559 if (componentIdx == H2OIdx)
560 DUNE_THROW(Dune::NotImplemented, "The component enthalpy for water in NAPL is not implemented.");
561 else if (componentIdx == NAPLIdx)
562 return NAPL::liquidEnthalpy(T, p);
563 else if (componentIdx == AirIdx)
564 DUNE_THROW(Dune::NotImplemented, "The component enthalpy for air in NAPL is not implemented.");
565 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx);
566 }
567 else if (phaseIdx == gPhaseIdx)
568 {
569 if (componentIdx == H2OIdx)
570 return H2O::gasEnthalpy(T, p);
571 else if (componentIdx == NAPLIdx)
572 return NAPL::gasEnthalpy(T, p);
573 else if (componentIdx == AirIdx)
574 return Air::gasEnthalpy(T,p);
575 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx);
576 }
577 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
578 }
579
582 template <class FluidState>
583 static Scalar heatCapacity(const FluidState &fluidState,
584 int phaseIdx)
585 {
586 DUNE_THROW(Dune::NotImplemented, "FluidSystems::H2OAirXylene::heatCapacity()");
587 }
588
591 template <class FluidState>
592 static Scalar thermalConductivity(const FluidState &fluidState,
593 int phaseIdx)
594 {
595 const Scalar temperature = fluidState.temperature(phaseIdx) ;
596 const Scalar pressure = fluidState.pressure(phaseIdx);
597 if (phaseIdx == wPhaseIdx)
598 {
599 return H2O::liquidThermalConductivity(temperature, pressure);
600 }
601 else if (phaseIdx == nPhaseIdx)
602 {
604 }
605 else if (phaseIdx == gPhaseIdx)
606 {
608 }
609 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
610 }
611
612private:
613 static Scalar waterPhaseDensity_(Scalar T, Scalar pw, Scalar xww, Scalar xwa, Scalar xwc)
614 {
615 Scalar rholH2O = H2O::liquidDensity(T, pw);
616 Scalar clH2O = rholH2O/H2O::molarMass();
617
618 // this assumes each dissolved molecule displaces exactly one
619 // water molecule in the liquid
620 return
621 clH2O*(xww*H2O::molarMass() + xwa*Air::molarMass() + xwc*NAPL::molarMass());
622 }
623
624 static Scalar gasPhaseDensity_(Scalar T, Scalar pg, Scalar xgw, Scalar xga, Scalar xgc)
625 {
626 return H2O::gasDensity(T, pg*xgw) + Air::gasDensity(T, pg*xga) + NAPL::gasDensity(T, pg*xgc);
627 }
628
629 static Scalar NAPLPhaseDensity_(Scalar T, Scalar pn)
630 {
631 return
632 NAPL::liquidDensity(T, pn);
633 }
634
635};
636
637} // end namespace FluidSystems
638} // end namespace Dumux
639
640#endif
A simple class for the air fluid properties.
Binary coefficients for air and xylene.
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:50
static Scalar liquidDiffCoeff(Scalar temperature, Scalar pressure)
Diffusion coefficient for air and xylene in liquid water.
Definition: air_xylene.hh:96
static Scalar henry(Scalar temperature)
Henry coefficient for air in liquid water.
Definition: h2o_air.hh:36
static Scalar gasDiffCoeff(Scalar temperature, Scalar pressure)
Binary diffusion coefficient for molecular water and air.
Definition: h2o_air.hh:56
static Scalar liquidDiffCoeff(Scalar temperature, Scalar pressure)
Diffusion coefficient for molecular nitrogen in liquid water.
Definition: h2o_air.hh:89
static Scalar liquidDiffCoeff(Scalar temperature, Scalar pressure)
Diffusion coefficient for xylene in liquid water.
Definition: h2o_xylene.hh:105
static Scalar henry(Scalar temperature)
Henry coefficient for xylene in liquid water.
Definition: h2o_xylene.hh:40
static Scalar gasDiffCoeff(Scalar temperature, Scalar pressure)
Binary diffusion coefficient for molecular water and xylene.
Definition: h2o_xylene.hh:57
A class for the air fluid properties.
Definition: air.hh:35
static Scalar gasDensity(Scalar temperature, Scalar pressure)
The density of Air at a given pressure and temperature.
Definition: air.hh:73
static constexpr Scalar molarMass()
The molar mass in of Air.
Definition: air.hh:50
static Scalar gasViscosity(Scalar temperature, Scalar pressure)
The dynamic viscosity of Air at a given pressure and temperature.
Definition: air.hh:181
static Scalar gasThermalConductivity(Scalar temperature, Scalar pressure)
Thermal conductivity of air.
Definition: air.hh:339
static constexpr bool gasIsIdeal()
Returns true, the gas phase is assumed to be ideal.
Definition: air.hh:97
static Scalar gasEnthalpy(Scalar temperature, Scalar pressure)
Specific enthalpy of Air with 273.15 as basis.
Definition: air.hh:264
static std::string name()
A human readable name for Air.
Definition: air.hh:42
static Scalar gasMolarDensity(Scalar temperature, Scalar pressure)
The molar density of air in , depending on pressure and temperature.
Definition: air.hh:85
Properties of xylene.
Definition: xylene.hh:37
static Scalar gasViscosity(Scalar temp, Scalar pressure)
The dynamic viscosity of xylene vapor.
Definition: xylene.hh:316
static std::string name()
A human readable name for the xylene.
Definition: xylene.hh:45
static Scalar liquidDensity(Scalar temperature, Scalar pressure)
The density of pure xylene at a given pressure and temperature .
Definition: xylene.hh:287
static Scalar liquidEnthalpy(const Scalar temperature, const Scalar pressure)
Specific enthalpy of liquid xylene .
Definition: xylene.hh:155
static Scalar gasMolarDensity(Scalar temperature, Scalar pressure)
The molar gas density of xylene gas at a given pressure and temperature.
Definition: xylene.hh:250
static constexpr bool gasIsIdeal()
Returns true if the gas phase is assumed to be ideal.
Definition: xylene.hh:301
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:94
static Scalar liquidThermalConductivity(Scalar temperature, Scalar pressure)
Thermal conductivity of xylene.
Definition: xylene.hh:372
static constexpr bool liquidIsCompressible()
Returns true if the liquid phase is assumed to be compressible.
Definition: xylene.hh:307
static Scalar gasEnthalpy(Scalar temperature, Scalar pressure)
Specific enthalpy of xylene vapor .
Definition: xylene.hh:226
static constexpr Scalar molarMass()
The molar mass in of xylene.
Definition: xylene.hh:51
static Scalar gasDensity(Scalar temperature, Scalar pressure)
The density of xylene gas at a given pressure and temperature.
Definition: xylene.hh:237
static Scalar liquidMolarDensity(Scalar temp, Scalar pressure)
The molar liquid density of pure xylene at a given pressure and temperature .
Definition: xylene.hh:262
static Scalar liquidViscosity(Scalar temp, Scalar pressure)
The dynamic viscosity of pure xylene.
Definition: xylene.hh:344
Fluid system base class.
Definition: fluidsystems/base.hh:33
A three-phase fluid system featuring gas, NAPL and water as phases and distilled water and air (Pseu...
Definition: h2oairxylene.hh:44
static Scalar diffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Calculate the binary molecular diffusion coefficient for a component in a fluid phase .
Definition: h2oairxylene.hh:371
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:237
static constexpr bool isGas(int phaseIdx)
Return whether a phase is gaseous.
Definition: h2oairxylene.hh:118
static const int gCompIdx
Definition: h2oairxylene.hh:68
static Scalar viscosity(const FluidState &fluidState, int phaseIdx)
Return the viscosity of a phase .
Definition: h2oairxylene.hh:301
static void init()
Initialize the fluid system's static parameters generically.
Definition: h2oairxylene.hh:76
static Scalar heatCapacity(const FluidState &fluidState, int phaseIdx)
Specific isobaric heat capacity of a fluid phase .
Definition: h2oairxylene.hh:583
static const int nCompIdx
Definition: h2oairxylene.hh:67
static Scalar molarMass(int compIdx)
Return the molar mass of a component in .
Definition: h2oairxylene.hh:213
static const int AirIdx
Definition: h2oairxylene.hh:61
static const int gPhaseIdx
Definition: h2oairxylene.hh:57
static std::string componentName(int compIdx)
Return the human readable name of a component (used in indices)
Definition: h2oairxylene.hh:199
static std::string phaseName(int phaseIdx)
Return the human readable name of a phase (used in indices)
Definition: h2oairxylene.hh:183
static const int numComponents
Definition: h2oairxylene.hh:53
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: h2oairxylene.hh:426
static const int H2OIdx
Definition: h2oairxylene.hh:59
static const int wCompIdx
Definition: h2oairxylene.hh:66
static constexpr bool isCompressible(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: h2oairxylene.hh:165
static constexpr bool isMiscible()
Returns whether the fluids are miscible.
Definition: h2oairxylene.hh:110
static bool isIdealGas(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: h2oairxylene.hh:130
H2OType H2O
Definition: h2oairxylene.hh:48
static const int nPhaseIdx
Definition: h2oairxylene.hh:56
static const int numPhases
Definition: h2oairxylene.hh:52
static const int NAPLIdx
Definition: h2oairxylene.hh:60
static Scalar fugacityCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Returns the fugacity coefficient of a component in a phase.
Definition: h2oairxylene.hh:449
static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
Calculate the molar density of a fluid phase.
Definition: h2oairxylene.hh:271
static Scalar componentEnthalpy(const FluidState &fluidState, int phaseIdx, int componentIdx)
Returns the specific enthalpy of a component in a specific phase.
Definition: h2oairxylene.hh:542
static Scalar enthalpy(const FluidState &fluidState, int phaseIdx)
Given all mole fractions in a phase, return the specific phase enthalpy .
Definition: h2oairxylene.hh:508
static const int wPhaseIdx
Definition: h2oairxylene.hh:55
static Scalar thermalConductivity(const FluidState &fluidState, int phaseIdx)
Thermal conductivity of a fluid phase .
Definition: h2oairxylene.hh:592
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:97
static bool isIdealMixture(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: h2oairxylene.hh:147
static Scalar kelvinVaporPressure(const FluidState &fluidState, const int phaseIdx, const int compIdx)
Definition: h2oairxylene.hh:490
Fluid system base class.
Material properties of pure water .
Binary coefficients for water and air.
Binary coefficients for water and xylene.
Relations valid for an ideal gas.
A collection of input/output field names for common physical quantities.
Scalar h2oGasViscosityInMixture(Scalar temperature, Scalar pressure)
The dynamic viscosity of steam in a gas mixture.
Definition: h2o.hh:961
std::string temperature() noexcept
I/O name of temperature for equilibrium models.
Definition: name.hh:39
std::string gaseousPhase() noexcept
I/O name of gaseous phase.
Definition: name.hh:111
std::string naplPhase() noexcept
I/O name of napl phase.
Definition: name.hh:119
std::string aqueousPhase() noexcept
I/O name of aqueous phase.
Definition: name.hh:115
std::string pressure(int phaseIdx) noexcept
I/O name of pressure for multiphase systems.
Definition: name.hh:22
Definition: adapt.hh:17
Tabulates all thermodynamic properties of a given untabulated chemical species.
Properties of xylene.