version 3.8
h2oairmesitylene.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_MESITYLENE_FLUID_SYSTEM_HH
13#define DUMUX_H2O_AIR_MESITYLENE_FLUID_SYSTEM_HH
14
21
25
27
28#include <dumux/io/name.hh>
29
30namespace Dumux {
31namespace FluidSystems {
32
41template <class Scalar,
42 class H2OType = Components::TabulatedComponent<Components::H2O<Scalar> > >
44 : public Base<Scalar, H2OAirMesitylene<Scalar, H2OType> >
45{
47
48public:
51 using H2O = H2OType;
52
53
54 static const int numPhases = 3;
55 static const int numComponents = 3;
56
57 static const int wPhaseIdx = 0; // index of the water phase
58 static const int nPhaseIdx = 1; // index of the NAPL phase
59 static const int gPhaseIdx = 2; // index of the gas phase
60
61 static const int H2OIdx = 0;
62 static const int NAPLIdx = 1;
63 static const int AirIdx = 2;
64
65 // export component indices to indicate the main component
66 // of the corresponding phase at atmospheric pressure 1 bar
67 // and room temperature 20°C:
68 static const int wCompIdx = H2OIdx;
69 static const int nCompIdx = NAPLIdx;
70 static const int gCompIdx = AirIdx;
71
78 static void init()
79 {
80 init(/*tempMin=*/273.15,
81 /*tempMax=*/623.15,
82 /*numTemp=*/100,
83 /*pMin=*/0.0,
84 /*pMax=*/20e6,
85 /*numP=*/200);
86 }
87
99 static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp,
100 Scalar pressMin, Scalar pressMax, unsigned nPress)
101 {
102 if (H2O::isTabulated)
103 {
104 H2O::init(tempMin, tempMax, nTemp,
105 pressMin, pressMax, nPress);
106 }
107 }
108
112 static constexpr bool isMiscible()
113 { return true; }
114
120 static constexpr bool isGas(int phaseIdx)
121 {
122 assert(0 <= phaseIdx && phaseIdx < numPhases);
123 return phaseIdx == gPhaseIdx;
124 }
125
132 static bool isIdealGas(int phaseIdx)
133 { return phaseIdx == gPhaseIdx && H2O::gasIsIdeal() && Air::gasIsIdeal() && NAPL::gasIsIdeal(); }
134
149 static bool isIdealMixture(int phaseIdx)
150 {
151 assert(0 <= phaseIdx && phaseIdx < numPhases);
152 // we assume Henry's and Raoult's laws for the water phase and
153 // and no interaction between gas molecules of different
154 // components, so all phases are ideal mixtures!
155 return true;
156 }
157
167 static constexpr bool isCompressible(int phaseIdx)
168 {
169 assert(0 <= phaseIdx && phaseIdx < numPhases);
170 // gases are always compressible
171 if (phaseIdx == gPhaseIdx)
172 return true;
173 else if (phaseIdx == wPhaseIdx)
174 // the water component decides for the water phase...
175 return H2O::liquidIsCompressible();
176
177 // the NAPL component decides for the napl phase...
179 }
180
184 static std::string phaseName(int phaseIdx)
185 {
186 assert(0 <= phaseIdx && phaseIdx < numPhases);
187 switch (phaseIdx)
188 {
189 case wPhaseIdx: return IOName::aqueousPhase();
190 case nPhaseIdx: return IOName::naplPhase();
191 case gPhaseIdx: return IOName::gaseousPhase();
192 }
193 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
194 }
195
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 // See: Eq. (7) in Class et al. (2002a)
241 // this assumes each dissolved molecule displaces exactly one
242 // water molecule in the liquid
243 return H2O::liquidMolarDensity(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx))
244 * (H2O::molarMass()*fluidState.moleFraction(wPhaseIdx, H2OIdx)
245 + Air::molarMass()*fluidState.moleFraction(wPhaseIdx, AirIdx)
246 + NAPL::molarMass()*fluidState.moleFraction(wPhaseIdx, NAPLIdx));
247 }
248 else if (phaseIdx == nPhaseIdx) {
249 // assume pure NAPL for the NAPL phase
250 Scalar pressure = NAPL::liquidIsCompressible()?fluidState.pressure(phaseIdx):1e100;
251 return NAPL::liquidDensity(fluidState.temperature(phaseIdx), pressure);
252 }
253
254 assert (phaseIdx == gPhaseIdx);
255 Scalar pH2O =
256 fluidState.moleFraction(gPhaseIdx, H2OIdx) *
257 fluidState.pressure(gPhaseIdx);
258 Scalar pAir =
259 fluidState.moleFraction(gPhaseIdx, AirIdx) *
260 fluidState.pressure(gPhaseIdx);
261 Scalar pNAPL =
262 fluidState.moleFraction(gPhaseIdx, NAPLIdx) *
263 fluidState.pressure(gPhaseIdx);
264 return H2O::gasDensity(fluidState.temperature(phaseIdx), pH2O)
265 + Air::gasDensity(fluidState.temperature(phaseIdx), pAir)
266 + NAPL::gasDensity(fluidState.temperature(phaseIdx), pNAPL);
267 }
268
271 template <class FluidState>
272 static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
273 {
274 Scalar temperature = fluidState.temperature(phaseIdx);
275 Scalar pressure = fluidState.pressure(phaseIdx);
276
277 if (phaseIdx == nPhaseIdx)
278 {
279 // assume pure NAPL for the NAPL phase
281 }
282 else if (phaseIdx == wPhaseIdx)
283 {
284 return H2O::liquidMolarDensity(temperature, pressure);
285 }
286 else
287 {
288 return H2O::gasMolarDensity(temperature, fluidState.partialPressure(gPhaseIdx, H2OIdx))
289 + NAPL::gasMolarDensity(temperature, fluidState.partialPressure(gPhaseIdx, NAPLIdx))
290 + Air::gasMolarDensity(temperature, fluidState.partialPressure(gPhaseIdx, AirIdx));
291 }
292 }
293
301 template <class FluidState>
302 static Scalar viscosity(const FluidState &fluidState,
303 int phaseIdx)
304 {
305 if (phaseIdx == wPhaseIdx) {
306 // assume pure water viscosity
307 return H2O::liquidViscosity(fluidState.temperature(phaseIdx),
308 fluidState.pressure(phaseIdx));
309 }
310 else if (phaseIdx == nPhaseIdx) {
311 // assume pure NAPL viscosity
312 return NAPL::liquidViscosity(fluidState.temperature(phaseIdx),
313 fluidState.pressure(phaseIdx));
314 }
315
316 assert (phaseIdx == gPhaseIdx);
317
318 /* Wilke method. See:
319 *
320 * See: R. Reid, et al.: The Properties of Gases and Liquids,
321 * 4th edition, McGraw-Hill, 1987, 407-410
322 * 5th edition, McGraw-Hill, 2001, p. 9.21/22
323 *
324 * in this case, we use a simplified version in order to avoid
325 * computationally costly evaluation of sqrt and pow functions and
326 * divisions
327 * -- compare e.g. with Promo Class p. 32/33
328 */
329 Scalar muResult;
330 const Scalar mu[numComponents] = {
331 h2oGasViscosityInMixture(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx)),
332 Air::gasViscosity(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx)),
333 NAPL::gasViscosity(fluidState.temperature(phaseIdx), NAPL::vaporPressure(fluidState.temperature(phaseIdx)))
334 };
335 // molar masses
336 const Scalar M[numComponents] = {
337 H2O::molarMass(),
340 };
341
342 Scalar muAW = mu[AirIdx]*fluidState.moleFraction(gPhaseIdx, AirIdx)
343 + mu[H2OIdx]*fluidState.moleFraction(gPhaseIdx, H2OIdx)
344 / (fluidState.moleFraction(gPhaseIdx, AirIdx)
345 + fluidState.moleFraction(gPhaseIdx, H2OIdx));
346 Scalar xAW = fluidState.moleFraction(gPhaseIdx, AirIdx)
347 + fluidState.moleFraction(gPhaseIdx, H2OIdx);
348
349 Scalar MAW = (fluidState.moleFraction(gPhaseIdx, AirIdx)*Air::molarMass()
350 + fluidState.moleFraction(gPhaseIdx, H2OIdx)*H2O::molarMass())
351 / xAW;
352
353 Scalar phiCAW = 0.3; // simplification for this particular system
354 /* actually like this
355 * using std::sqrt;
356 * using std::pow;
357 * Scalar phiCAW = pow(1.+sqrt(mu[NAPLIdx]/muAW)*pow(MAW/M[NAPLIdx],0.25),2)
358 * / sqrt(8.*(1.+M[NAPLIdx]/MAW));
359 */
360 Scalar phiAWC = phiCAW * muAW*M[NAPLIdx]/(mu[NAPLIdx]*MAW);
361
362 muResult = (xAW*muAW)/(xAW+fluidState.moleFraction(gPhaseIdx, NAPLIdx)*phiAWC)
363 + (fluidState.moleFraction(gPhaseIdx, NAPLIdx) * mu[NAPLIdx])
364 / (fluidState.moleFraction(gPhaseIdx, NAPLIdx) + xAW*phiCAW);
365 return muResult;
366 }
367
368
377 template <class FluidState>
378 static Scalar diffusionCoefficient(const FluidState &fluidState,
379 int phaseIdx,
380 int compIdx)
381 {
382 switch (phaseIdx)
383 {
384 case gPhaseIdx:
385 {
386 switch (compIdx)
387 {
388 case NAPLIdx:
389 {
390 Scalar diffWC = BinaryCoeff::H2O_Mesitylene::gasDiffCoeff(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
391 Scalar diffAW = BinaryCoeff::H2O_Air::gasDiffCoeff(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
392 const Scalar xga = fluidState.moleFraction(gPhaseIdx, AirIdx);
393 const Scalar xgw = fluidState.moleFraction(gPhaseIdx, H2OIdx);
394 const Scalar xgc = fluidState.moleFraction(gPhaseIdx, NAPLIdx);
395 return (1.- xgw)/(xga/diffAW + xgc/diffWC);
396 }
397 case H2OIdx:
398 {
399 Scalar diffAC = BinaryCoeff::Air_Mesitylene::gasDiffCoeff(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
400 Scalar diffWC = BinaryCoeff::H2O_Mesitylene::gasDiffCoeff(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
401 const Scalar xga = fluidState.moleFraction(gPhaseIdx, AirIdx);
402 const Scalar xgw = fluidState.moleFraction(gPhaseIdx, H2OIdx);
403 const Scalar xgc = fluidState.moleFraction(gPhaseIdx, NAPLIdx);
404 return (1.- xgc)/(xgw/diffWC + xga/diffAC);
405 }
406 case AirIdx:
407 DUNE_THROW(Dune::InvalidStateException, "Diffusivity of Air in the gas phase is constraint by sum of diffusive fluxes = 0 !");
408 }
409 }
410 case wPhaseIdx:
411 {
412 Scalar diffACl = BinaryCoeff::Air_Mesitylene::liquidDiffCoeff(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
413 Scalar diffWCl = BinaryCoeff::H2O_Mesitylene::liquidDiffCoeff(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
414 Scalar diffAWl = BinaryCoeff::H2O_Air::liquidDiffCoeff(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
415
416 Scalar xwa = fluidState.moleFraction(wPhaseIdx, AirIdx);
417 Scalar xww = fluidState.moleFraction(wPhaseIdx, H2OIdx);
418 Scalar xwc = fluidState.moleFraction(wPhaseIdx, NAPLIdx);
419
420 switch (compIdx)
421 {
422 case NAPLIdx:
423 return (1.- xww)/(xwa/diffAWl + xwc/diffWCl);
424 case AirIdx:
425 return (1.- xwc)/(xww/diffWCl + xwa/diffACl);
426 case H2OIdx:
427 DUNE_THROW(Dune::InvalidStateException,
428 "Diffusivity of water in the water phase "
429 "is constraint by sum of diffusive fluxes = 0 !\n");
430 }
431 }
432 case nPhaseIdx:
433 {
434 return 0;
435 }
436 }
437 return 0;
438 }
439
442 template <class FluidState>
443 static Scalar binaryDiffusionCoefficient(const FluidState &fluidState,
444 int phaseIdx,
445 int compIIdx,
446 int compJIdx)
447 {
448 DUNE_THROW(Dune::NotImplemented, "FluidSystems::H2OAirMesitylene::binaryDiffusionCoefficient()");
449 }
450
465 template <class FluidState>
466 static Scalar fugacityCoefficient(const FluidState &fluidState,
467 int phaseIdx,
468 int compIdx)
469 {
470 assert(0 <= phaseIdx && phaseIdx < numPhases);
471 assert(0 <= compIdx && compIdx < numComponents);
472
473 Scalar T = fluidState.temperature(phaseIdx);
474 Scalar p = fluidState.pressure(phaseIdx);
475
476 if (phaseIdx == wPhaseIdx) {
477 if (compIdx == H2OIdx)
478 return H2O::vaporPressure(T)/p;
479 else if (compIdx == AirIdx)
480 return BinaryCoeff::H2O_Air::henry(T)/p;
481 else if (compIdx == NAPLIdx)
483 }
484
485 // for the NAPL phase, we assume currently that nothing is
486 // dissolved. this means that the affinity of the NAPL
487 // component to the NAPL phase is much higher than for the
488 // other components, i.e. the fugacity coefficient is much
489 // smaller.
490 if (phaseIdx == nPhaseIdx) {
491 Scalar phiNapl = NAPL::vaporPressure(T)/p;
492 if (compIdx == NAPLIdx)
493 return phiNapl;
494 else if (compIdx == AirIdx)
495 return 1e6*phiNapl;
496 else if (compIdx == H2OIdx)
497 return 1e6*phiNapl;
498 }
499
500 // for the gas phase, assume an ideal gas when it comes to
501 // fugacity (-> fugacity == partial pressure)
502 assert(phaseIdx == gPhaseIdx);
503 return 1.0;
504 }
505
506 template <class FluidState>
507 static Scalar kelvinVaporPressure(const FluidState &fluidState,
508 const int phaseIdx,
509 const int compIdx)
510 {
511 DUNE_THROW(Dune::NotImplemented, "FluidSystems::H2OAirMesitylene::kelvinVaporPressure()");
512 }
513
524 template <class FluidState>
525 static Scalar enthalpy(const FluidState &fluidState,
526 int phaseIdx)
527 {
528 if (phaseIdx == wPhaseIdx) {
529 return H2O::liquidEnthalpy(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
530 }
531 else if (phaseIdx == nPhaseIdx) {
532 return NAPL::liquidEnthalpy(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
533 }
534 else if (phaseIdx == gPhaseIdx) { // gas phase enthalpy depends strongly on composition
535 Scalar hgc = NAPL::gasEnthalpy(fluidState.temperature(phaseIdx),
536 fluidState.pressure(phaseIdx));
537 Scalar hgw = H2O::gasEnthalpy(fluidState.temperature(phaseIdx),
538 fluidState.pressure(phaseIdx));
539 // pressure is only a dummy here (not dependent on pressure, just temperature)
540 Scalar hga = Air::gasEnthalpy(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
541
542 Scalar result = 0;
543 result += hgw * fluidState.massFraction(gPhaseIdx, H2OIdx);
544 result += hga * fluidState.massFraction(gPhaseIdx, AirIdx);
545 result += hgc * fluidState.massFraction(gPhaseIdx, NAPLIdx);
546
547 return result;
548 }
549 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
550 }
551
558 template <class FluidState>
559 static Scalar componentEnthalpy(const FluidState& fluidState, int phaseIdx, int componentIdx)
560 {
561 const Scalar T = fluidState.temperature(phaseIdx);
562 const Scalar p = fluidState.pressure(phaseIdx);
563
564 if (phaseIdx == wPhaseIdx)
565 {
566 if (componentIdx == H2OIdx)
567 return H2O::liquidEnthalpy(T, p);
568 else if (componentIdx == NAPLIdx)
569 DUNE_THROW(Dune::NotImplemented, "The component enthalpy for NAPL in water is not implemented.");
570 else if (componentIdx == AirIdx)
571 DUNE_THROW(Dune::NotImplemented, "The component enthalpy for Air in water is not implemented.");
572 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx);
573 }
574 else if (phaseIdx == nPhaseIdx)
575 {
576 if (componentIdx == H2OIdx)
577 DUNE_THROW(Dune::NotImplemented, "The component enthalpy for water in NAPL is not implemented.");
578 else if (componentIdx == NAPLIdx)
579 return NAPL::liquidEnthalpy(T, p);
580 else if (componentIdx == AirIdx)
581 DUNE_THROW(Dune::NotImplemented, "The component enthalpy for air in NAPL is not implemented.");
582 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx);
583 }
584 else if (phaseIdx == gPhaseIdx)
585 {
586 if (componentIdx == H2OIdx)
587 return H2O::gasEnthalpy(T, p);
588 else if (componentIdx == NAPLIdx)
589 return NAPL::gasEnthalpy(T, p);
590 else if (componentIdx == AirIdx)
591 return Air::gasEnthalpy(T,p);
592 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx);
593 }
594 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
595 }
596
599 template <class FluidState>
600 static Scalar heatCapacity(const FluidState &fluidState,
601 int phaseIdx)
602 {
603 DUNE_THROW(Dune::NotImplemented, "FluidSystems::H2OAirMesitylene::heatCapacity()");
604 }
605
608 template <class FluidState>
609 static Scalar thermalConductivity(const FluidState &fluidState,
610 int phaseIdx)
611 {
612 const Scalar temperature = fluidState.temperature(phaseIdx) ;
613 const Scalar pressure = fluidState.pressure(phaseIdx);
614 if (phaseIdx == wPhaseIdx)
615 {
616 return H2O::liquidThermalConductivity(temperature, pressure);
617 }
618 else if (phaseIdx == nPhaseIdx)
619 {
621 }
622 else if (phaseIdx == gPhaseIdx)
623 {
625 }
626 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
627 }
628
629private:
630 static Scalar waterPhaseDensity_(Scalar T, Scalar pw, Scalar xww, Scalar xwa, Scalar xwc)
631 {
632 Scalar rholH2O = H2O::liquidDensity(T, pw);
633 Scalar clH2O = rholH2O/H2O::molarMass();
634
635 // this assumes each dissolved molecule displaces exactly one
636 // water molecule in the liquid
637 return clH2O*(xww*H2O::molarMass() + xwa*Air::molarMass() + xwc*NAPL::molarMass());
638 }
639
640 static Scalar gasPhaseDensity_(Scalar T, Scalar pg, Scalar xgw, Scalar xga, Scalar xgc)
641 {
642 return H2O::gasDensity(T, pg*xgw) + Air::gasDensity(T, pg*xga) + NAPL::gasDensity(T, pg*xgc);
643 }
644
645 static Scalar NAPLPhaseDensity_(Scalar T, Scalar pn)
646 {
647 return NAPL::liquidDensity(T, pn);
648 }
649
650};
651
652} // end namespace FluidSystems
653} // end namespace Dumux
654
655#endif
A simple class for the air fluid properties.
Binary coefficients for air and mesitylene.
static Scalar liquidDiffCoeff(Scalar temperature, Scalar pressure)
Diffusion coefficient for air and mesitylene in liquid water.
Definition: air_mesitylene.hh:97
static Scalar gasDiffCoeff(Scalar temperature, Scalar pressure)
Binary diffusion coefficient for air and mesitylene. I used the method according to Wilke and Lee se...
Definition: air_mesitylene.hh:47
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 mesitylene in liquid water.
Definition: h2o_mesitylene.hh:104
static Scalar henry(Scalar temperature)
Henry coefficient for mesitylene in liquid water.
Definition: h2o_mesitylene.hh:38
static Scalar gasDiffCoeff(Scalar temperature, Scalar pressure)
Binary diffusion coefficient for molecular water and mesitylene.
Definition: h2o_mesitylene.hh:54
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
mesitylene
Definition: mesitylene.hh:36
static Scalar gasMolarDensity(Scalar temperature, Scalar pressure)
The molar density of mesitylene in , depending on pressure and temperature.
Definition: mesitylene.hh:207
static Scalar gasEnthalpy(Scalar temperature, Scalar pressure)
Specific enthalpy of mesitylene vapor .
Definition: mesitylene.hh:183
static std::string name()
A human readable name for the mesitylene.
Definition: mesitylene.hh:43
static Scalar liquidMolarDensity(Scalar temperature, Scalar pressure)
The molar density of pure mesitylene at a given pressure and temperature .
Definition: mesitylene.hh:230
static constexpr Scalar molarMass()
The molar mass in of mesitylene.
Definition: mesitylene.hh:49
static Scalar liquidThermalConductivity(Scalar temperature, Scalar pressure)
Thermal conductivity of mesitylene.
Definition: mesitylene.hh:366
static Scalar gasViscosity(Scalar temperature, Scalar pressure)
The dynamic viscosity of mesitylene vapor.
Definition: mesitylene.hh:270
static Scalar vaporPressure(Scalar temperature)
The saturation vapor pressure in of pure mesitylene at a given temperature according to Antoine afte...
Definition: mesitylene.hh:93
static constexpr bool liquidIsCompressible()
Returns true if the liquid phase is assumed to be compressible.
Definition: mesitylene.hh:261
static Scalar gasDensity(Scalar temperature, Scalar pressure)
The density of mesitylene at a given pressure and temperature .
Definition: mesitylene.hh:194
static Scalar liquidEnthalpy(const Scalar temperature, const Scalar pressure)
Specific enthalpy of liquid mesitylene .
Definition: mesitylene.hh:112
static Scalar liquidDensity(Scalar temperature, Scalar pressure)
The density of pure mesitylene at a given pressure and temperature .
Definition: mesitylene.hh:216
static constexpr bool gasIsIdeal()
Returns true if the gas phase is assumed to be ideal.
Definition: mesitylene.hh:255
static Scalar liquidViscosity(Scalar temperature, Scalar pressure)
The dynamic viscosity of pure mesitylene.
Definition: mesitylene.hh:300
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: h2oairmesitylene.hh:45
static const int numComponents
Definition: h2oairmesitylene.hh:55
static Scalar thermalConductivity(const FluidState &fluidState, int phaseIdx)
Thermal conductivity of a fluid phase .
Definition: h2oairmesitylene.hh:609
static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
Calculate the molar density of a fluid phase.
Definition: h2oairmesitylene.hh:272
static void init()
Initialize the fluid system's static parameters generically.
Definition: h2oairmesitylene.hh:78
static const int numPhases
Definition: h2oairmesitylene.hh:54
static const int nPhaseIdx
Definition: h2oairmesitylene.hh:58
static Scalar componentEnthalpy(const FluidState &fluidState, int phaseIdx, int componentIdx)
Returns the specific enthalpy of a component in a specific phase.
Definition: h2oairmesitylene.hh:559
static const int H2OIdx
Definition: h2oairmesitylene.hh:61
static const int AirIdx
Definition: h2oairmesitylene.hh:63
static Scalar molarMass(int compIdx)
Return the molar mass of a component in .
Definition: h2oairmesitylene.hh:213
static Scalar viscosity(const FluidState &fluidState, int phaseIdx)
Return the viscosity of a phase .
Definition: h2oairmesitylene.hh:302
static Scalar enthalpy(const FluidState &fluidState, int phaseIdx)
Given all mole fractions in a phase, return the specific phase enthalpy .
Definition: h2oairmesitylene.hh:525
static bool isIdealGas(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: h2oairmesitylene.hh:132
static Scalar diffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Given all mole fractions, return the diffusion coefficient in of a component in a phase.
Definition: h2oairmesitylene.hh:378
static const int wPhaseIdx
Definition: h2oairmesitylene.hh:57
static constexpr bool isGas(int phaseIdx)
Return whether a phase is gaseous.
Definition: h2oairmesitylene.hh:120
static constexpr bool isCompressible(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: h2oairmesitylene.hh:167
static const int wCompIdx
Definition: h2oairmesitylene.hh:68
static Scalar fugacityCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Returns the fugacity coefficient of a component in a phase.
Definition: h2oairmesitylene.hh:466
static const int gPhaseIdx
Definition: h2oairmesitylene.hh:59
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: h2oairmesitylene.hh:99
static constexpr bool isMiscible()
Returns whether the fluids are miscible.
Definition: h2oairmesitylene.hh:112
static Scalar heatCapacity(const FluidState &fluidState, int phaseIdx)
Specific isobaric heat capacity of a fluid phase .
Definition: h2oairmesitylene.hh:600
static bool isIdealMixture(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: h2oairmesitylene.hh:149
static std::string componentName(int compIdx)
Return the human readable name of a component (used in indices)
Definition: h2oairmesitylene.hh:199
static Scalar density(const FluidState &fluidState, int phaseIdx)
Given a phase's composition, temperature, pressure, and the partial pressures of all components,...
Definition: h2oairmesitylene.hh:237
static const int NAPLIdx
Definition: h2oairmesitylene.hh:62
H2OType H2O
Definition: h2oairmesitylene.hh:51
static const int gCompIdx
Definition: h2oairmesitylene.hh:70
static const int nCompIdx
Definition: h2oairmesitylene.hh:69
static std::string phaseName(int phaseIdx)
Return the human readable name of a phase (used in indices)
Definition: h2oairmesitylene.hh:184
static Scalar kelvinVaporPressure(const FluidState &fluidState, const int phaseIdx, const int compIdx)
Definition: h2oairmesitylene.hh:507
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: h2oairmesitylene.hh:443
Fluid system base class.
Material properties of pure water .
Binary coefficients for water and air.
Binary coefficients for water and mesitylene.
Relations valid for an ideal gas.
Properties of mesitylene.
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.