3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * 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_MESITYLENE_FLUID_SYSTEM_HH
25#define DUMUX_H2O_AIR_MESITYLENE_FLUID_SYSTEM_HH
26
33
37
39
40#include <dumux/io/name.hh>
41
42namespace Dumux {
43namespace FluidSystems {
44
53template <class Scalar,
54 class H2OType = Components::TabulatedComponent<Components::H2O<Scalar> > >
56 : public Base<Scalar, H2OAirMesitylene<Scalar, H2OType> >
57{
60
61public:
64 using H2O = H2OType;
65
66
67 static const int numPhases = 3;
68 static const int numComponents = 3;
69
70 static const int wPhaseIdx = 0; // index of the water phase
71 static const int nPhaseIdx = 1; // index of the NAPL phase
72 static const int gPhaseIdx = 2; // index of the gas phase
73
74 static const int H2OIdx = 0;
75 static const int NAPLIdx = 1;
76 static const int AirIdx = 2;
77
78 // export component indices to indicate the main component
79 // of the corresponding phase at atmospheric pressure 1 bar
80 // and room temperature 20°C:
81 static const int wCompIdx = H2OIdx;
82 static const int nCompIdx = NAPLIdx;
83 static const int gCompIdx = AirIdx;
84
91 static void init()
92 {
93 init(/*tempMin=*/273.15,
94 /*tempMax=*/623.15,
95 /*numTemp=*/100,
96 /*pMin=*/0.0,
97 /*pMax=*/20e6,
98 /*numP=*/200);
99 }
100
112 static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp,
113 Scalar pressMin, Scalar pressMax, unsigned nPress)
114 {
115 if (H2O::isTabulated)
116 {
117 H2O::init(tempMin, tempMax, nTemp,
118 pressMin, pressMax, nPress);
119 }
120 }
121
125 static constexpr bool isMiscible()
126 { return true; }
127
133 static constexpr bool isGas(int phaseIdx)
134 {
135 assert(0 <= phaseIdx && phaseIdx < numPhases);
136 return phaseIdx == gPhaseIdx;
137 }
138
145 static bool isIdealGas(int phaseIdx)
146 { return phaseIdx == gPhaseIdx && H2O::gasIsIdeal() && Air::gasIsIdeal() && NAPL::gasIsIdeal(); }
147
162 static bool isIdealMixture(int phaseIdx)
163 {
164 assert(0 <= phaseIdx && phaseIdx < numPhases);
165 // we assume Henry's and Raoult's laws for the water phase and
166 // and no interaction between gas molecules of different
167 // components, so all phases are ideal mixtures!
168 return true;
169 }
170
180 static constexpr bool isCompressible(int phaseIdx)
181 {
182 assert(0 <= phaseIdx && phaseIdx < numPhases);
183 // gases are always compressible
184 if (phaseIdx == gPhaseIdx)
185 return true;
186 else if (phaseIdx == wPhaseIdx)
187 // the water component decides for the water phase...
188 return H2O::liquidIsCompressible();
189
190 // the NAPL component decides for the napl phase...
192 }
193
197 static std::string phaseName(int phaseIdx)
198 {
199 assert(0 <= phaseIdx && phaseIdx < numPhases);
200 switch (phaseIdx)
201 {
202 case wPhaseIdx: return IOName::aqueousPhase();
203 case nPhaseIdx: return IOName::naplPhase();
204 case gPhaseIdx: return IOName::gaseousPhase();
205 }
206 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
207 }
208
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 // See: Eq. (7) in Class et al. (2002a)
254 // this assumes each dissolved molecule displaces exactly one
255 // water molecule in the liquid
256 return H2O::liquidMolarDensity(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx))
257 * (H2O::molarMass()*fluidState.moleFraction(wPhaseIdx, H2OIdx)
258 + Air::molarMass()*fluidState.moleFraction(wPhaseIdx, AirIdx)
259 + NAPL::molarMass()*fluidState.moleFraction(wPhaseIdx, NAPLIdx));
260 }
261 else if (phaseIdx == nPhaseIdx) {
262 // assume pure NAPL for the NAPL phase
263 Scalar pressure = NAPL::liquidIsCompressible()?fluidState.pressure(phaseIdx):1e100;
264 return NAPL::liquidDensity(fluidState.temperature(phaseIdx), pressure);
265 }
266
267 assert (phaseIdx == gPhaseIdx);
268 Scalar pH2O =
269 fluidState.moleFraction(gPhaseIdx, H2OIdx) *
270 fluidState.pressure(gPhaseIdx);
271 Scalar pAir =
272 fluidState.moleFraction(gPhaseIdx, AirIdx) *
273 fluidState.pressure(gPhaseIdx);
274 Scalar pNAPL =
275 fluidState.moleFraction(gPhaseIdx, NAPLIdx) *
276 fluidState.pressure(gPhaseIdx);
277 return H2O::gasDensity(fluidState.temperature(phaseIdx), pH2O)
278 + Air::gasDensity(fluidState.temperature(phaseIdx), pAir)
279 + NAPL::gasDensity(fluidState.temperature(phaseIdx), pNAPL);
280 }
281
282 using Base::molarDensity;
292 template <class FluidState>
293 static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
294 {
295 Scalar temperature = fluidState.temperature(phaseIdx);
296 Scalar pressure = fluidState.pressure(phaseIdx);
297
298 if (phaseIdx == nPhaseIdx)
299 {
300 // assume pure NAPL for the NAPL phase
302 }
303 else if (phaseIdx == wPhaseIdx)
304 {
305 return H2O::liquidMolarDensity(temperature, pressure);
306 }
307 else
308 {
309 return H2O::gasMolarDensity(temperature, fluidState.partialPressure(gPhaseIdx, H2OIdx))
310 + NAPL::gasMolarDensity(temperature, fluidState.partialPressure(gPhaseIdx, NAPLIdx))
311 + Air::gasMolarDensity(temperature, fluidState.partialPressure(gPhaseIdx, AirIdx));
312 }
313 }
314
315 using Base::viscosity;
322 template <class FluidState>
323 static Scalar viscosity(const FluidState &fluidState,
324 int phaseIdx)
325 {
326 if (phaseIdx == wPhaseIdx) {
327 // assume pure water viscosity
328 return H2O::liquidViscosity(fluidState.temperature(phaseIdx),
329 fluidState.pressure(phaseIdx));
330 }
331 else if (phaseIdx == nPhaseIdx) {
332 // assume pure NAPL viscosity
333 return NAPL::liquidViscosity(fluidState.temperature(phaseIdx),
334 fluidState.pressure(phaseIdx));
335 }
336
337 assert (phaseIdx == gPhaseIdx);
338
339 /* Wilke method. See:
340 *
341 * See: R. Reid, et al.: The Properties of Gases and Liquids,
342 * 4th edition, McGraw-Hill, 1987, 407-410
343 * 5th edition, McGraw-Hill, 2001, p. 9.21/22
344 *
345 * in this case, we use a simplified version in order to avoid
346 * computationally costly evaluation of sqrt and pow functions and
347 * divisions
348 * -- compare e.g. with Promo Class p. 32/33
349 */
350 Scalar muResult;
351 const Scalar mu[numComponents] = {
352 H2O::gasViscosity(fluidState.temperature(phaseIdx), H2O::vaporPressure(fluidState.temperature(phaseIdx))),
353 Air::gasViscosity(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx)),
354 NAPL::gasViscosity(fluidState.temperature(phaseIdx), NAPL::vaporPressure(fluidState.temperature(phaseIdx)))
355 };
356 // molar masses
357 const Scalar M[numComponents] = {
358 H2O::molarMass(),
361 };
362
363 Scalar muAW = mu[AirIdx]*fluidState.moleFraction(gPhaseIdx, AirIdx)
364 + mu[H2OIdx]*fluidState.moleFraction(gPhaseIdx, H2OIdx)
365 / (fluidState.moleFraction(gPhaseIdx, AirIdx)
366 + fluidState.moleFraction(gPhaseIdx, H2OIdx));
367 Scalar xAW = fluidState.moleFraction(gPhaseIdx, AirIdx)
368 + fluidState.moleFraction(gPhaseIdx, H2OIdx);
369
370 Scalar MAW = (fluidState.moleFraction(gPhaseIdx, AirIdx)*Air::molarMass()
371 + fluidState.moleFraction(gPhaseIdx, H2OIdx)*H2O::molarMass())
372 / xAW;
373
374 Scalar phiCAW = 0.3; // simplification for this particular system
375 /* actually like this
376 * using std::sqrt;
377 * using std::pow;
378 * Scalar phiCAW = pow(1.+sqrt(mu[NAPLIdx]/muAW)*pow(MAW/M[NAPLIdx],0.25),2)
379 * / sqrt(8.*(1.+M[NAPLIdx]/MAW));
380 */
381 Scalar phiAWC = phiCAW * muAW*M[NAPLIdx]/(mu[NAPLIdx]*MAW);
382
383 muResult = (xAW*muAW)/(xAW+fluidState.moleFraction(gPhaseIdx, NAPLIdx)*phiAWC)
384 + (fluidState.moleFraction(gPhaseIdx, NAPLIdx) * mu[NAPLIdx])
385 / (fluidState.moleFraction(gPhaseIdx, NAPLIdx) + xAW*phiCAW);
386 return muResult;
387 }
388
389
398 template <class FluidState>
399 static Scalar diffusionCoefficient(const FluidState &fluidState,
400 int phaseIdx,
401 int compIdx)
402 {
403 switch (phaseIdx)
404 {
405 case gPhaseIdx:
406 {
407 switch (compIdx)
408 {
409 case NAPLIdx:
410 {
411 Scalar diffWC = BinaryCoeff::H2O_Mesitylene::gasDiffCoeff(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
412 Scalar diffAW = BinaryCoeff::H2O_Air::gasDiffCoeff(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
413 const Scalar xga = fluidState.moleFraction(gPhaseIdx, AirIdx);
414 const Scalar xgw = fluidState.moleFraction(gPhaseIdx, H2OIdx);
415 const Scalar xgc = fluidState.moleFraction(gPhaseIdx, NAPLIdx);
416 return (1.- xgw)/(xga/diffAW + xgc/diffWC);
417 }
418 case H2OIdx:
419 {
420 Scalar diffAC = BinaryCoeff::Air_Mesitylene::gasDiffCoeff(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
421 Scalar diffWC = BinaryCoeff::H2O_Mesitylene::gasDiffCoeff(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
422 const Scalar xga = fluidState.moleFraction(gPhaseIdx, AirIdx);
423 const Scalar xgw = fluidState.moleFraction(gPhaseIdx, H2OIdx);
424 const Scalar xgc = fluidState.moleFraction(gPhaseIdx, NAPLIdx);
425 return (1.- xgc)/(xgw/diffWC + xga/diffAC);
426 }
427 case AirIdx:
428 DUNE_THROW(Dune::InvalidStateException, "Diffusivity of Air in the gas phase is constraint by sum of diffusive fluxes = 0 !");
429 }
430 }
431 case wPhaseIdx:
432 {
433 Scalar diffACl = BinaryCoeff::Air_Mesitylene::liquidDiffCoeff(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
434 Scalar diffWCl = BinaryCoeff::H2O_Mesitylene::liquidDiffCoeff(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
435 Scalar diffAWl = BinaryCoeff::H2O_Air::liquidDiffCoeff(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
436
437 Scalar xwa = fluidState.moleFraction(wPhaseIdx, AirIdx);
438 Scalar xww = fluidState.moleFraction(wPhaseIdx, H2OIdx);
439 Scalar xwc = fluidState.moleFraction(wPhaseIdx, NAPLIdx);
440
441 switch (compIdx)
442 {
443 case NAPLIdx:
444 return (1.- xww)/(xwa/diffAWl + xwc/diffWCl);
445 case AirIdx:
446 return (1.- xwc)/(xww/diffWCl + xwa/diffACl);
447 case H2OIdx:
448 DUNE_THROW(Dune::InvalidStateException,
449 "Diffusivity of water in the water phase "
450 "is constraint by sum of diffusive fluxes = 0 !\n");
451 }
452 }
453 case nPhaseIdx:
454 DUNE_THROW(Dune::InvalidStateException, "Diffusion coefficients of substances in non-wetting liquid phase are undefined!\n");
455 }
456 return 0;
457 }
458
460 template <class FluidState>
461 static Scalar binaryDiffusionCoefficient(const FluidState &fluidState,
462 int phaseIdx,
463 int compIIdx,
464 int compJIdx)
465 {
466 DUNE_THROW(Dune::NotImplemented, "FluidSystems::H2OAirMesitylene::binaryDiffusionCoefficient()");
467 }
468
483 template <class FluidState>
484 static Scalar fugacityCoefficient(const FluidState &fluidState,
485 int phaseIdx,
486 int compIdx)
487 {
488 assert(0 <= phaseIdx && phaseIdx < numPhases);
489 assert(0 <= compIdx && compIdx < numComponents);
490
491 Scalar T = fluidState.temperature(phaseIdx);
492 Scalar p = fluidState.pressure(phaseIdx);
493
494 if (phaseIdx == wPhaseIdx) {
495 if (compIdx == H2OIdx)
496 return H2O::vaporPressure(T)/p;
497 else if (compIdx == AirIdx)
498 return BinaryCoeff::H2O_Air::henry(T)/p;
499 else if (compIdx == NAPLIdx)
501 }
502
503 // for the NAPL phase, we assume currently that nothing is
504 // dissolved. this means that the affinity of the NAPL
505 // component to the NAPL phase is much higher than for the
506 // other components, i.e. the fugacity coefficient is much
507 // smaller.
508 if (phaseIdx == nPhaseIdx) {
509 Scalar phiNapl = NAPL::vaporPressure(T)/p;
510 if (compIdx == NAPLIdx)
511 return phiNapl;
512 else if (compIdx == AirIdx)
513 return 1e6*phiNapl;
514 else if (compIdx == H2OIdx)
515 return 1e6*phiNapl;
516 }
517
518 // for the gas phase, assume an ideal gas when it comes to
519 // fugacity (-> fugacity == partial pressure)
520 assert(phaseIdx == gPhaseIdx);
521 return 1.0;
522 }
523
524 template <class FluidState>
525 static Scalar kelvinVaporPressure(const FluidState &fluidState,
526 const int phaseIdx,
527 const int compIdx)
528 {
529 DUNE_THROW(Dune::NotImplemented, "FluidSystems::H2OAirMesitylene::kelvinVaporPressure()");
530 }
531
532 using Base::enthalpy;
542 template <class FluidState>
543 static Scalar enthalpy(const FluidState &fluidState,
544 int phaseIdx)
545 {
546 if (phaseIdx == wPhaseIdx) {
547 return H2O::liquidEnthalpy(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
548 }
549 else if (phaseIdx == nPhaseIdx) {
550 return NAPL::liquidEnthalpy(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
551 }
552 else if (phaseIdx == gPhaseIdx) { // gas phase enthalpy depends strongly on composition
553 Scalar hgc = NAPL::gasEnthalpy(fluidState.temperature(phaseIdx),
554 fluidState.pressure(phaseIdx));
555 Scalar hgw = H2O::gasEnthalpy(fluidState.temperature(phaseIdx),
556 fluidState.pressure(phaseIdx));
557 // pressure is only a dummy here (not dependent on pressure, just temperature)
558 Scalar hga = Air::gasEnthalpy(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
559
560 Scalar result = 0;
561 result += hgw * fluidState.massFraction(gPhaseIdx, H2OIdx);
562 result += hga * fluidState.massFraction(gPhaseIdx, AirIdx);
563 result += hgc * fluidState.massFraction(gPhaseIdx, NAPLIdx);
564
565 return result;
566 }
567 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
568 }
569
576 template <class FluidState>
577 static Scalar componentEnthalpy(const FluidState& fluidState, int phaseIdx, int componentIdx)
578 {
579 const Scalar T = fluidState.temperature(phaseIdx);
580 const Scalar p = fluidState.pressure(phaseIdx);
581
582 if (phaseIdx == wPhaseIdx)
583 {
584 if (componentIdx == H2OIdx)
585 return H2O::liquidEnthalpy(T, p);
586 else if (componentIdx == NAPLIdx)
587 DUNE_THROW(Dune::NotImplemented, "The component enthalpy for NAPL in water is not implemented.");
588 else if (componentIdx == AirIdx)
589 DUNE_THROW(Dune::NotImplemented, "The component enthalpy for Air in water is not implemented.");
590 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx);
591 }
592 else if (phaseIdx == nPhaseIdx)
593 {
594 if (componentIdx == H2OIdx)
595 DUNE_THROW(Dune::NotImplemented, "The component enthalpy for water in NAPL is not implemented.");
596 else if (componentIdx == NAPLIdx)
597 return NAPL::liquidEnthalpy(T, p);
598 else if (componentIdx == AirIdx)
599 DUNE_THROW(Dune::NotImplemented, "The component enthalpy for air in NAPL is not implemented.");
600 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx);
601 }
602 else if (phaseIdx == gPhaseIdx)
603 {
604 if (componentIdx == H2OIdx)
605 return H2O::gasEnthalpy(T, p);
606 else if (componentIdx == NAPLIdx)
607 return NAPL::gasEnthalpy(T, p);
608 else if (componentIdx == AirIdx)
609 return Air::gasEnthalpy(T,p);
610 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx);
611 }
612 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
613 }
614
615 using Base::heatCapacity;
621 template <class FluidState>
622 static Scalar heatCapacity(const FluidState &fluidState,
623 int phaseIdx)
624 {
625 DUNE_THROW(Dune::NotImplemented, "FluidSystems::H2OAirMesitylene::heatCapacity()");
626 }
627
634 template <class FluidState>
635 static Scalar thermalConductivity(const FluidState &fluidState,
636 int phaseIdx)
637 {
638 const Scalar temperature = fluidState.temperature(phaseIdx) ;
639 const Scalar pressure = fluidState.pressure(phaseIdx);
640 if (phaseIdx == wPhaseIdx)
641 {
642 return H2O::liquidThermalConductivity(temperature, pressure);
643 }
644 else if (phaseIdx == nPhaseIdx)
645 {
647 }
648 else if (phaseIdx == gPhaseIdx)
649 {
651 }
652 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
653 }
654
655private:
656 static Scalar waterPhaseDensity_(Scalar T, Scalar pw, Scalar xww, Scalar xwa, Scalar xwc)
657 {
658 Scalar rholH2O = H2O::liquidDensity(T, pw);
659 Scalar clH2O = rholH2O/H2O::molarMass();
660
661 // this assumes each dissolved molecule displaces exactly one
662 // water molecule in the liquid
663 return clH2O*(xww*H2O::molarMass() + xwa*Air::molarMass() + xwc*NAPL::molarMass());
664 }
665
666 static Scalar gasPhaseDensity_(Scalar T, Scalar pg, Scalar xgw, Scalar xga, Scalar xgc)
667 {
668 return H2O::gasDensity(T, pg*xgw) + Air::gasDensity(T, pg*xga) + NAPL::gasDensity(T, pg*xgc);
669 }
670
671 static Scalar NAPLPhaseDensity_(Scalar T, Scalar pn)
672 {
673 return NAPL::liquidDensity(T, pn);
674 }
675
676};
677
678} // end namespace FluidSystems
679} // end namespace Dumux
680
681#endif
A collection of input/output field names for common physical quantities.
Binary coefficients for air and mesitylene.
Binary coefficients for water and air.
Binary coefficients for water and mesitylene.
A simple class for the air fluid properties.
Material properties of pure water .
Properties of mesitylene.
Tabulates all thermodynamic properties of a given untabulated chemical species.
Relations valid for an ideal gas.
make the local view function available whenever we use the grid geometry
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
static Scalar liquidDiffCoeff(Scalar temperature, Scalar pressure)
Diffusion coefficient for air and mesitylene in liquid water.
Definition: air_mesitylene.hh:108
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:59
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 mesitylene in liquid water.
Definition: h2o_mesitylene.hh:114
static Scalar henry(Scalar temperature)
Henry coefficient for mesitylene in liquid water.
Definition: h2o_mesitylene.hh:47
static Scalar gasDiffCoeff(Scalar temperature, Scalar pressure)
Binary diffusion coefficient for molecular water and mesitylene.
Definition: h2o_mesitylene.hh:63
A class for the air fluid properties.
Definition: air.hh:46
static Scalar gasDensity(Scalar temperature, Scalar pressure)
The density of Air at a given pressure and temperature.
Definition: air.hh:84
static constexpr Scalar molarMass()
The molar mass in of Air.
Definition: air.hh:61
static Scalar gasViscosity(Scalar temperature, Scalar pressure)
The dynamic viscosity of Air at a given pressure and temperature.
Definition: air.hh:186
static Scalar gasThermalConductivity(Scalar temperature, Scalar pressure)
Thermal conductivity of air.
Definition: air.hh:342
static constexpr bool gasIsIdeal()
Returns true, the gas phase is assumed to be ideal.
Definition: air.hh:108
static Scalar gasEnthalpy(Scalar temperature, Scalar pressure)
Specific enthalpy of Air with 273.15 as basis.
Definition: air.hh:268
static std::string name()
A human readable name for Air.
Definition: air.hh:53
static Scalar gasMolarDensity(Scalar temperature, Scalar pressure)
The molar density of air in , depending on pressure and temperature.
Definition: air.hh:96
mesitylene
Definition: mesitylene.hh:48
static Scalar gasMolarDensity(Scalar temperature, Scalar pressure)
The molar density of mesitylene in , depending on pressure and temperature.
Definition: mesitylene.hh:219
static Scalar gasEnthalpy(Scalar temperature, Scalar pressure)
Specific enthalpy of mesitylene vapor .
Definition: mesitylene.hh:195
static std::string name()
A human readable name for the mesitylene.
Definition: mesitylene.hh:55
static Scalar liquidMolarDensity(Scalar temperature, Scalar pressure)
The molar density of pure mesitylene at a given pressure and temperature .
Definition: mesitylene.hh:242
static constexpr Scalar molarMass()
The molar mass in of mesitylene.
Definition: mesitylene.hh:61
static Scalar liquidThermalConductivity(Scalar temperature, Scalar pressure)
Thermal conductivity of mesitylene.
Definition: mesitylene.hh:379
static Scalar vaporPressure(Scalar temperature)
The saturation vapor pressure in of pure mesitylene at a given temperature according to Antoine afte...
Definition: mesitylene.hh:105
static constexpr bool liquidIsCompressible()
Returns true if the liquid phase is assumed to be compressible.
Definition: mesitylene.hh:273
static Scalar gasDensity(Scalar temperature, Scalar pressure)
The density of mesitylene at a given pressure and temperature .
Definition: mesitylene.hh:206
static Scalar liquidEnthalpy(const Scalar temperature, const Scalar pressure)
Specific enthalpy of liquid mesitylene .
Definition: mesitylene.hh:124
static Scalar liquidDensity(Scalar temperature, Scalar pressure)
The density of pure mesitylene at a given pressure and temperature .
Definition: mesitylene.hh:228
static Scalar gasViscosity(Scalar temperature, Scalar pressure, bool regularize=true)
The dynamic viscosity of mesitylene vapor.
Definition: mesitylene.hh:283
static constexpr bool gasIsIdeal()
Returns true if the gas phase is assumed to be ideal.
Definition: mesitylene.hh:267
static Scalar liquidViscosity(Scalar temperature, Scalar pressure)
The dynamic viscosity of pure mesitylene.
Definition: mesitylene.hh:313
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: h2oairmesitylene.hh:57
static const int numComponents
Definition: h2oairmesitylene.hh:68
static Scalar thermalConductivity(const FluidState &fluidState, int phaseIdx)
Return the thermal conductivity .
Definition: h2oairmesitylene.hh:635
static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
The molar density of a fluid phase in .
Definition: h2oairmesitylene.hh:293
static void init()
Initialize the fluid system's static parameters generically.
Definition: h2oairmesitylene.hh:91
static const int numPhases
Definition: h2oairmesitylene.hh:67
static const int nPhaseIdx
Definition: h2oairmesitylene.hh:71
static Scalar componentEnthalpy(const FluidState &fluidState, int phaseIdx, int componentIdx)
Returns the specific enthalpy of a component in a specific phase.
Definition: h2oairmesitylene.hh:577
static const int H2OIdx
Definition: h2oairmesitylene.hh:74
static const int AirIdx
Definition: h2oairmesitylene.hh:76
static Scalar molarMass(int compIdx)
Return the molar mass of a component in .
Definition: h2oairmesitylene.hh:226
static Scalar viscosity(const FluidState &fluidState, int phaseIdx)
Return the viscosity of a phase .
Definition: h2oairmesitylene.hh:323
static Scalar enthalpy(const FluidState &fluidState, int phaseIdx)
Given all mole fractions in a phase, return the specific phase enthalpy .
Definition: h2oairmesitylene.hh:543
static bool isIdealGas(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: h2oairmesitylene.hh:145
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:399
static const int wPhaseIdx
Definition: h2oairmesitylene.hh:70
static constexpr bool isGas(int phaseIdx)
Return whether a phase is gaseous.
Definition: h2oairmesitylene.hh:133
static constexpr bool isCompressible(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: h2oairmesitylene.hh:180
static const int wCompIdx
Definition: h2oairmesitylene.hh:81
static Scalar fugacityCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Returns the fugacity coefficient of a component in a phase.
Definition: h2oairmesitylene.hh:484
static const int gPhaseIdx
Definition: h2oairmesitylene.hh:72
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:112
static constexpr bool isMiscible()
Returns whether the fluids are miscible.
Definition: h2oairmesitylene.hh:125
static Scalar heatCapacity(const FluidState &fluidState, int phaseIdx)
Return the heat capacity in .
Definition: h2oairmesitylene.hh:622
static bool isIdealMixture(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: h2oairmesitylene.hh:162
static std::string componentName(int compIdx)
Return the human readable name of a component (used in indices)
Definition: h2oairmesitylene.hh:212
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:250
static const int NAPLIdx
Definition: h2oairmesitylene.hh:75
H2OType H2O
Definition: h2oairmesitylene.hh:64
static const int gCompIdx
Definition: h2oairmesitylene.hh:83
static const int nCompIdx
Definition: h2oairmesitylene.hh:82
static std::string phaseName(int phaseIdx)
Return the human readable name of a phase (used in indices)
Definition: h2oairmesitylene.hh:197
static Scalar kelvinVaporPressure(const FluidState &fluidState, const int phaseIdx, const int compIdx)
Definition: h2oairmesitylene.hh:525
static Scalar binaryDiffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIIdx, int compJIdx)
Definition: h2oairmesitylene.hh:461
Fluid system base class.