3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
h2on2.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_N2_FLUID_SYSTEM_HH
25#define DUMUX_H2O_N2_FLUID_SYSTEM_HH
26
27#include <cassert>
28#include <iomanip>
29
31
33
38
39#include <dumux/io/name.hh>
40
41#include "base.hh"
42
43namespace Dumux {
44namespace FluidSystems {
49template<bool fastButSimplifiedRelations = false>
51{
52 static constexpr bool useH2ODensityAsLiquidMixtureDensity() { return fastButSimplifiedRelations; }
53 static constexpr bool useIdealGasDensity() { return fastButSimplifiedRelations; }
54 static constexpr bool useN2ViscosityAsGasMixtureViscosity() { return fastButSimplifiedRelations; }
55 static constexpr bool useN2HeatConductivityAsGasMixtureHeatConductivity() { return fastButSimplifiedRelations; }
56 static constexpr bool useIdealGasHeatCapacities() { return fastButSimplifiedRelations; }
57};
58
65template <class Scalar, class Policy = H2ON2DefaultPolicy<>>
66class H2ON2
67 : public Base<Scalar, H2ON2<Scalar, Policy> >
68{
71
72 // convenience using declarations
76
77public:
78 using H2O = TabulatedH2O;
79 using N2 = SimpleN2;
80
81 static constexpr int numPhases = 2;
82 static constexpr int numComponents = 2;
83
84 static constexpr int liquidPhaseIdx = 0;
85 static constexpr int gasPhaseIdx = 1;
86 static constexpr int phase0Idx = liquidPhaseIdx;
87 static constexpr int phase1Idx = gasPhaseIdx;
88
89 static constexpr int H2OIdx = 0;
90 static constexpr int N2Idx = 1;
91 static constexpr int comp0Idx = H2OIdx;
92 static constexpr int comp1Idx = N2Idx;
93 static constexpr int liquidCompIdx = H2OIdx;
94 static constexpr int gasCompIdx = N2Idx;
95
96 /****************************************
97 * Fluid phase related static parameters
98 ****************************************/
104 static std::string phaseName(int phaseIdx)
105 {
106 assert(0 <= phaseIdx && phaseIdx < numPhases);
107 switch (phaseIdx)
108 {
109 case liquidPhaseIdx: return IOName::liquidPhase();
110 case gasPhaseIdx: return IOName::gaseousPhase();
111 }
112 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
113 }
114
118 static constexpr bool isMiscible()
119 { return true; }
120
126 static constexpr bool isGas(int phaseIdx)
127 {
128 assert(0 <= phaseIdx && phaseIdx < numPhases);
129 return phaseIdx == gasPhaseIdx;
130 }
131
146 static bool isIdealMixture(int phaseIdx)
147 {
148 assert(0 <= phaseIdx && phaseIdx < numPhases);
149 // we assume Henry's and Raoult's laws for the water phase and
150 // and no interaction between gas molecules of different
151 // components, so all phases are ideal mixtures!
152 return true;
153 }
154
164 static constexpr bool isCompressible(int phaseIdx)
165 {
166 assert(0 <= phaseIdx && phaseIdx < numPhases);
167 // gases are always compressible
168 if (phaseIdx == gasPhaseIdx)
169 return true;
170 // the water component decides for the liquid phase...
172 }
173
180 static bool isIdealGas(int phaseIdx)
181 {
182 assert(0 <= phaseIdx && phaseIdx < numPhases);
183
184 if (phaseIdx == gasPhaseIdx)
185 // let the components decide
186 return H2O::gasIsIdeal() && N2::gasIsIdeal();
187 return false; // not a gas
188 }
189
190 /****************************************
191 * Component related static parameters
192 ****************************************/
198 static std::string componentName(int compIdx)
199 {
200 switch (compIdx)
201 {
202 case H2OIdx: return H2O::name();
203 case N2Idx: return N2::name();
204 }
205
206 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx);
207 }
208
214 static Scalar molarMass(int compIdx)
215 {
216 static const Scalar M[] = {
219 };
220
221 assert(0 <= compIdx && compIdx < numComponents);
222 return M[compIdx];
223 }
224
230 static Scalar criticalTemperature(int compIdx)
231 {
232 static const Scalar Tcrit[] = {
235 };
236
237 assert(0 <= compIdx && compIdx < numComponents);
238 return Tcrit[compIdx];
239 }
240
246 static Scalar criticalPressure(int compIdx)
247 {
248 static const Scalar pcrit[] = {
251 };
252
253 assert(0 <= compIdx && compIdx < numComponents);
254 return pcrit[compIdx];
255 }
256
266 template <class FluidState>
267 static Scalar kelvinVaporPressure(const FluidState &fluidState,
268 const int phaseIdx,
269 const int compIdx)
270 {
271 assert(compIdx == H2OIdx && phaseIdx == liquidPhaseIdx);
272
273 using std::exp;
274 return fugacityCoefficient(fluidState, phaseIdx, compIdx)
275 * fluidState.pressure(phaseIdx)
276 * exp(-(fluidState.pressure(gasPhaseIdx)-fluidState.pressure(liquidPhaseIdx))
277 / density(fluidState, phaseIdx)
279 / fluidState.temperature());
280 }
281
287 static Scalar criticalMolarVolume(int compIdx)
288 {
289 DUNE_THROW(Dune::NotImplemented,
290 "H2ON2FluidSystem::criticalMolarVolume()");
291 }
292
298 static Scalar acentricFactor(int compIdx)
299 {
300 static const Scalar accFac[] = {
301 H2O::acentricFactor(),
302 N2::acentricFactor()
303 };
304
305 assert(0 <= compIdx && compIdx < numComponents);
306 return accFac[compIdx];
307 }
308
309 /****************************************
310 * thermodynamic relations
311 ****************************************/
312
319 static void init()
320 {
321 init(/*tempMin=*/273.15,
322 /*tempMax=*/623.15,
323 /*numTemp=*/100,
324 /*pMin=*/0.0,
325 /*pMax=*/20e6,
326 /*numP=*/200);
327 }
328
340 static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp,
341 Scalar pressMin, Scalar pressMax, unsigned nPress)
342 {
343 std::cout << "The H2O-N2 fluid system was configured with the following policy:\n";
344 std::cout << " - use H2O density as liquid mixture density: " << std::boolalpha << Policy::useH2ODensityAsLiquidMixtureDensity() << "\n";
345 std::cout << " - use ideal gas density: " << std::boolalpha << Policy::useIdealGasDensity() << "\n";
346 std::cout << " - use N2 viscosity as gas mixture viscosity: " << std::boolalpha << Policy::useN2ViscosityAsGasMixtureViscosity() << "\n";
347 std::cout << " - use N2 heat conductivity as gas mixture heat conductivity: " << std::boolalpha << Policy::useN2HeatConductivityAsGasMixtureHeatConductivity() << "\n";
348 std::cout << " - use ideal gas heat capacities: " << std::boolalpha << Policy::useIdealGasHeatCapacities() << std::endl;
349
351 {
352 TabulatedH2O::init(tempMin, tempMax, nTemp,
353 pressMin, pressMax, nPress);
354 }
355 }
356
357 using Base::density;
370 template <class FluidState>
371 static Scalar density(const FluidState &fluidState,
372 int phaseIdx)
373 {
374 assert(0 <= phaseIdx && phaseIdx < numPhases);
375
376 Scalar T = fluidState.temperature(phaseIdx);
377 Scalar p = fluidState.pressure(phaseIdx);
378
379 // liquid phase
380 if (phaseIdx == liquidPhaseIdx) {
381 if (Policy::useH2ODensityAsLiquidMixtureDensity())
382 // assume pure water
383 return H2O::liquidDensity(T, p);
384 else
385 {
386 // See: Eq. (7) in Class et al. (2002a)
387 // This assumes each gas molecule displaces exactly one
388 // molecule in the liquid.
389 return H2O::liquidMolarDensity(T, p)
390 * (H2O::molarMass()*fluidState.moleFraction(liquidPhaseIdx, H2OIdx)
391 + N2::molarMass()*fluidState.moleFraction(liquidPhaseIdx, N2Idx));
392 }
393 }
394
395 // gas phase
396 using std::max;
397 if (Policy::useIdealGasDensity())
398 // for the gas phase assume an ideal gas
399 {
400 const Scalar averageMolarMass = fluidState.averageMolarMass(gasPhaseIdx);
401 return IdealGas::density(averageMolarMass, T, p);
402 }
403
404 // assume ideal mixture: steam and nitrogen don't "see" each other
405 Scalar rho_gH2O = H2O::gasDensity(T, fluidState.partialPressure(gasPhaseIdx, H2OIdx));
406 Scalar rho_gN2 = N2::gasDensity(T, fluidState.partialPressure(gasPhaseIdx, N2Idx));
407 return (rho_gH2O + rho_gN2);
408 }
409
410 using Base::molarDensity;
423 template <class FluidState>
424 static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
425 {
426 assert(0 <= phaseIdx && phaseIdx < numPhases);
427
428 Scalar T = fluidState.temperature(phaseIdx);
429 Scalar p = fluidState.pressure(phaseIdx);
430
431 // liquid phase
432 if (phaseIdx == liquidPhaseIdx)
433 {
434 // assume pure water or that each gas molecule displaces exactly one
435 // molecule in the liquid.
436 return H2O::liquidMolarDensity(T, p);
437 }
438
439 // gas phase
440 using std::max;
441 if (Policy::useIdealGasDensity())
442 // for the gas phase assume an ideal gas
443 {
444 return IdealGas::molarDensity(T, p);
445 }
446
447 // assume ideal mixture: steam and nitrogen don't "see" each other
448 Scalar rho_gH2O = H2O::gasMolarDensity(T, fluidState.partialPressure(gasPhaseIdx, H2OIdx));
449 Scalar rho_gN2 = N2::gasMolarDensity(T, fluidState.partialPressure(gasPhaseIdx, N2Idx));
450 return rho_gH2O + rho_gN2;
451 }
452
453 using Base::viscosity;
466 template <class FluidState>
467 static Scalar viscosity(const FluidState &fluidState,
468 int phaseIdx)
469 {
470 assert(0 <= phaseIdx && phaseIdx < numPhases);
471
472 Scalar T = fluidState.temperature(phaseIdx);
473 Scalar p = fluidState.pressure(phaseIdx);
474
475 // liquid phase
476 if (phaseIdx == liquidPhaseIdx) {
477 // assume pure water for the liquid phase
478 return H2O::liquidViscosity(T, p);
479 }
480
481 // gas phase
482 if (Policy::useN2ViscosityAsGasMixtureViscosity())
483 {
484 // assume pure nitrogen for the gas phase
485 return N2::gasViscosity(T, p);
486 }
487 else
488 {
489 // Wilke method (Reid et al.):
490 Scalar muResult = 0;
491 const Scalar mu[numComponents] = {
493 N2::gasViscosity(T, p)
494 };
495
496 Scalar sumx = 0.0;
497 using std::max;
498 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
499 sumx += fluidState.moleFraction(phaseIdx, compIdx);
500 sumx = max(1e-10, sumx);
501
502 for (int i = 0; i < numComponents; ++i) {
503 Scalar divisor = 0;
504// using std::sqrt;
505// using std::pow;
506 for (int j = 0; j < numComponents; ++j) {
507 Scalar phiIJ = 1 + sqrt(mu[i]/mu[j]) * pow(molarMass(j)/molarMass(i), 1/4.0);
508 phiIJ *= phiIJ;
509 phiIJ /= sqrt(8*(1 + molarMass(i)/molarMass(j)));
510 divisor += fluidState.moleFraction(phaseIdx, j)/sumx * phiIJ;
511 }
512 muResult += fluidState.moleFraction(phaseIdx, i)/sumx * mu[i] / divisor;
513 }
514
515 return muResult;
516 }
517 }
518
547 template <class FluidState>
548 static Scalar fugacityCoefficient(const FluidState &fluidState,
549 int phaseIdx,
550 int compIdx)
551 {
552 assert(0 <= phaseIdx && phaseIdx < numPhases);
553 assert(0 <= compIdx && compIdx < numComponents);
554
555 Scalar T = fluidState.temperature(phaseIdx);
556 Scalar p = fluidState.pressure(phaseIdx);
557
558 // liquid phase
559 if (phaseIdx == liquidPhaseIdx) {
560 if (compIdx == H2OIdx)
561 return H2O::vaporPressure(T)/p;
562 return BinaryCoeff::H2O_N2::henry(T)/p;
563 }
564
565 // for the gas phase, assume an ideal gas when it comes to
566 // fugacity (-> fugacity == partial pressure)
567 return 1.0;
568 }
569
594 template <class FluidState>
595 static Scalar diffusionCoefficient(const FluidState &fluidState,
596 int phaseIdx,
597 int compIdx)
598 {
599 DUNE_THROW(Dune::NotImplemented, "Diffusion coefficients");
600 }
601
613 template <class FluidState>
614 static Scalar binaryDiffusionCoefficient(const FluidState &fluidState,
615 int phaseIdx,
616 int compIIdx,
617 int compJIdx)
618
619 {
620 if (compIIdx > compJIdx)
621 {
622 using std::swap;
623 swap(compIIdx, compJIdx);
624 }
625
626 const Scalar T = fluidState.temperature(phaseIdx);
627 const Scalar p = fluidState.pressure(phaseIdx);
628
629 if (phaseIdx == liquidPhaseIdx && compIIdx == H2OIdx && compJIdx == N2Idx)
631
632 else if (phaseIdx == gasPhaseIdx && compIIdx == H2OIdx && compJIdx == N2Idx)
634
635 else
636 DUNE_THROW(Dune::InvalidStateException,
637 "Binary diffusion coefficient of components "
638 << compIIdx << " and " << compJIdx
639 << " in phase " << phaseIdx << " is unavailable!\n");
640 }
641
642 using Base::enthalpy;
655 template <class FluidState>
656 static Scalar enthalpy(const FluidState &fluidState,
657 int phaseIdx)
658 {
659 const Scalar T = fluidState.temperature(phaseIdx);
660 const Scalar p = fluidState.pressure(phaseIdx);
661
662 // liquid phase
663 if (phaseIdx == liquidPhaseIdx) {
664 return H2O::liquidEnthalpy(T, p);
665 }
666 // gas phase
667 else {
668 // assume ideal mixture: which means
669 // that the total specific enthalpy is the sum of the
670 // "partial specific enthalpies" of the components.
671 Scalar hH2O =
672 fluidState.massFraction(gasPhaseIdx, H2OIdx)
673 * H2O::gasEnthalpy(T, p);
674 Scalar hN2 =
675 fluidState.massFraction(gasPhaseIdx, N2Idx)
676 * N2::gasEnthalpy(T, p);
677 return hH2O + hN2;
678 }
679 }
680
687 template <class FluidState>
688 static Scalar componentEnthalpy(const FluidState &fluidState,
689 int phaseIdx,
690 int componentIdx)
691 {
692 const Scalar T = fluidState.temperature(phaseIdx);
693 const Scalar p = fluidState.pressure(phaseIdx);
694
695 if (phaseIdx == liquidPhaseIdx)
696 {
697 if (componentIdx == H2OIdx)
698 return H2O::liquidEnthalpy(T, p);
699 else if (componentIdx == N2Idx)
700 DUNE_THROW(Dune::NotImplemented, "Component enthalpy of nitrogen in liquid phase");
701 else
702 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx);
703 }
704 else if (phaseIdx == gasPhaseIdx)
705 {
706 if (componentIdx == H2OIdx)
707 return H2O::gasEnthalpy(T, p);
708 else if (componentIdx == N2Idx)
709 return N2::gasEnthalpy(T, p);
710 else
711 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx);
712 }
713 else
714 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
715 }
716
726 template <class FluidState>
727 static Scalar thermalConductivity(const FluidState &fluidState,
728 const int phaseIdx)
729 {
730 assert(0 <= phaseIdx && phaseIdx < numPhases);
731
732 const Scalar temperature = fluidState.temperature(phaseIdx) ;
733 const Scalar pressure = fluidState.pressure(phaseIdx);
734 if (phaseIdx == liquidPhaseIdx)
735 {
737 }
738 else
739 {
741 if (!Policy::useN2HeatConductivityAsGasMixtureHeatConductivity())
742 {
743 Scalar xN2 = fluidState.moleFraction(phaseIdx, N2Idx);
744 Scalar xH2O = fluidState.moleFraction(phaseIdx, H2OIdx);
745 Scalar lambdaN2 = xN2 * lambdaPureN2;
746 Scalar partialPressure = pressure * xH2O;
747 Scalar lambdaH2O = xH2O * H2O::gasThermalConductivity(temperature, partialPressure);
748 return lambdaN2 + lambdaH2O;
749 }
750 else
751 return lambdaPureN2;
752 }
753 }
754
755 using Base::heatCapacity;
763 template <class FluidState>
764 static Scalar heatCapacity(const FluidState &fluidState,
765 int phaseIdx)
766 {
767 if (phaseIdx == liquidPhaseIdx) {
768 return H2O::liquidHeatCapacity(fluidState.temperature(phaseIdx),
769 fluidState.pressure(phaseIdx));
770 }
771
772 // for the gas phase, assume ideal mixture
773 Scalar c_pN2;
774 Scalar c_pH2O;
775 // let the water and nitrogen components do things their own way
776 if (!Policy::useIdealGasHeatCapacities()) {
777 c_pN2 = N2::gasHeatCapacity(fluidState.temperature(phaseIdx),
778 fluidState.pressure(phaseIdx)
779 * fluidState.moleFraction(phaseIdx, N2Idx));
780
781 c_pH2O = H2O::gasHeatCapacity(fluidState.temperature(phaseIdx),
782 fluidState.pressure(phaseIdx)
783 * fluidState.moleFraction(phaseIdx, H2OIdx));
784 }
785 else {
786 // assume an ideal gas for both components. See:
787 // http://en.wikipedia.org/wiki/Heat_capacity
788 Scalar c_vN2molar = Constants<Scalar>::R*2.39;
789 Scalar c_pN2molar = Constants<Scalar>::R + c_vN2molar;
790
791 Scalar c_vH2Omolar = Constants<Scalar>::R*3.37; // <- correct??
792 Scalar c_pH2Omolar = Constants<Scalar>::R + c_vH2Omolar;
793
794 c_pN2 = c_pN2molar/molarMass(N2Idx);
795 c_pH2O = c_pH2Omolar/molarMass(H2OIdx);
796 }
797
798 // mangle both components together
799 return c_pH2O*fluidState.massFraction(gasPhaseIdx, H2OIdx)
800 + c_pN2*fluidState.massFraction(gasPhaseIdx, N2Idx);
801 }
802};
803
804} // end namespace FluidSystems
805
806} // end namespace Dumux
807
808#endif
Some exceptions thrown in DuMux
A collection of input/output field names for common physical quantities.
Relations valid for an ideal gas.
Material properties of pure water .
Tabulates all thermodynamic properties of a given untabulated chemical species.
Properties of pure molecular nitrogen .
Binary coefficients for water and nitrogen.
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 liquidPhase() noexcept
I/O name of liquid phase.
Definition: name.hh:119
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 liquidDiffCoeff(Scalar temperature, Scalar pressure)
Diffusion coefficient for molecular nitrogen in liquid water.
Definition: h2o_n2.hh:97
static Scalar henry(Scalar temperature)
Henry coefficient for molecular nitrogen in liquid water.
Definition: h2o_n2.hh:48
static Scalar gasDiffCoeff(Scalar temperature, Scalar pressure)
Binary diffusion coefficient for molecular water and nitrogen.
Definition: h2o_n2.hh:66
Properties of pure molecular nitrogen .
Definition: n2.hh:47
static Scalar criticalTemperature()
Returns the critical temperature of molecular nitrogen.
Definition: n2.hh:66
static Scalar criticalPressure()
Returns the critical pressure of molecular nitrogen.
Definition: n2.hh:72
static Scalar gasViscosity(Scalar temperature, Scalar pressure)
The dynamic viscosity of at a given pressure and temperature.
Definition: n2.hh:243
static Scalar gasMolarDensity(Scalar temperature, Scalar pressure)
The molar density of gas in at a given pressure and temperature.
Definition: n2.hh:143
static std::string name()
A human readable name for nitrogen.
Definition: n2.hh:54
static Scalar gasDensity(Scalar temperature, Scalar pressure)
The density of gas at a given pressure and temperature.
Definition: n2.hh:130
static constexpr Scalar molarMass()
The molar mass in of molecular nitrogen.
Definition: n2.hh:60
static Scalar gasThermalConductivity(Scalar temperature, Scalar pressure)
Thermal conductivity of nitrogen.
Definition: n2.hh:281
static constexpr bool gasIsIdeal()
Returns true if the gas phase is assumed to be ideal.
Definition: n2.hh:155
static const Scalar gasEnthalpy(Scalar temperature, Scalar pressure)
Specific enthalpy of pure nitrogen gas.
Definition: n2.hh:176
static const Scalar gasHeatCapacity(Scalar T, Scalar pressure)
Specific isobaric heat capacity of pure nitrogen gas.
Definition: n2.hh:213
Tabulates all thermodynamic properties of a given untabulated chemical species.
Definition: tabulatedcomponent.hh:82
static const Scalar gasEnthalpy(Scalar temperature, Scalar pressure)
Specific enthalpy of the gas .
Definition: tabulatedcomponent.hh:238
static Scalar criticalTemperature()
Returns the critical temperature in of the component.
Definition: tabulatedcomponent.hh:184
static const Scalar gasHeatCapacity(Scalar temperature, Scalar pressure)
Specific isobaric heat capacity of the gas .
Definition: tabulatedcomponent.hh:292
static std::string name()
A human readable name for the component.
Definition: tabulatedcomponent.hh:172
static Scalar liquidThermalConductivity(Scalar temperature, Scalar pressure)
The thermal conductivity of liquid water .
Definition: tabulatedcomponent.hh:619
static Scalar gasDensity(Scalar temperature, Scalar pressure)
The density of gas at a given pressure and temperature .
Definition: tabulatedcomponent.hh:456
static const Scalar liquidEnthalpy(Scalar temperature, Scalar pressure)
Specific enthalpy of the liquid .
Definition: tabulatedcomponent.hh:265
static Scalar criticalPressure()
Returns the critical pressure in of the component.
Definition: tabulatedcomponent.hh:190
static Scalar liquidMolarDensity(Scalar temperature, Scalar pressure)
The molar density of liquid in at a given pressure and temperature.
Definition: tabulatedcomponent.hh:529
static constexpr Scalar molarMass()
The molar mass in of the component.
Definition: tabulatedcomponent.hh:178
static Scalar liquidViscosity(Scalar temperature, Scalar pressure)
The dynamic viscosity of liquid.
Definition: tabulatedcomponent.hh:565
static constexpr bool liquidIsCompressible()
Returns true if the liquid phase is assumed to be compressible.
Definition: tabulatedcomponent.hh:439
static Scalar liquidDensity(Scalar temperature, Scalar pressure)
The density of liquid at a given pressure and temperature .
Definition: tabulatedcomponent.hh:495
static Scalar gasThermalConductivity(Scalar temperature, Scalar pressure)
The thermal conductivity of gaseous water .
Definition: tabulatedcomponent.hh:592
static Scalar gasMolarDensity(Scalar temperature, Scalar pressure)
The molar density of gas in at a given pressure and temperature.
Definition: tabulatedcomponent.hh:485
static constexpr bool isTabulated
state that we are tabulated
Definition: tabulatedcomponent.hh:88
static Scalar vaporPressure(Scalar T)
The vapor pressure in of the component at a given temperature.
Definition: tabulatedcomponent.hh:211
static void init(Scalar tempMin, Scalar tempMax, std::size_t nTemp, Scalar pressMin, Scalar pressMax, std::size_t nPress)
Initialize the tables.
Definition: tabulatedcomponent.hh:100
static constexpr bool gasIsIdeal()
Returns true if the gas phase is assumed to be ideal.
Definition: tabulatedcomponent.hh:445
static const Scalar liquidHeatCapacity(Scalar temperature, Scalar pressure)
Specific isobaric heat capacity of the liquid .
Definition: tabulatedcomponent.hh:319
A central place for various physical constants occuring in some equations.
Definition: constants.hh:39
Fluid system base class.
Definition: fluidsystems/base.hh:45
Scalar Scalar
export the scalar type
Definition: fluidsystems/base.hh:48
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
Policy for the H2O-N2 fluid system.
Definition: h2on2.hh:51
static constexpr bool useIdealGasHeatCapacities()
Definition: h2on2.hh:56
static constexpr bool useN2HeatConductivityAsGasMixtureHeatConductivity()
Definition: h2on2.hh:55
static constexpr bool useIdealGasDensity()
Definition: h2on2.hh:53
static constexpr bool useH2ODensityAsLiquidMixtureDensity()
Definition: h2on2.hh:52
static constexpr bool useN2ViscosityAsGasMixtureViscosity()
Definition: h2on2.hh:54
A two-phase fluid system with two components water Nitrogen for non-equilibrium models.
Definition: h2on2.hh:68
static constexpr bool isGas(int phaseIdx)
Return whether a phase is gaseous.
Definition: h2on2.hh:126
static constexpr int liquidPhaseIdx
index of the liquid phase
Definition: h2on2.hh:84
static bool isIdealMixture(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: h2on2.hh:146
static Scalar molarMass(int compIdx)
Return the molar mass of a component in .
Definition: h2on2.hh:214
static Scalar heatCapacity(const FluidState &fluidState, int phaseIdx)
Specific isobaric heat capacity of a fluid phase. .
Definition: h2on2.hh:764
static Scalar thermalConductivity(const FluidState &fluidState, const int phaseIdx)
Thermal conductivity of a fluid phase .
Definition: h2on2.hh:727
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: h2on2.hh:614
static constexpr int phase1Idx
index of the second phase
Definition: h2on2.hh:87
static constexpr int liquidCompIdx
index of the liquid component
Definition: h2on2.hh:93
static Scalar criticalMolarVolume(int compIdx)
Molar volume of a component at the critical point .
Definition: h2on2.hh:287
static constexpr int N2Idx
Definition: h2on2.hh:90
static constexpr bool isCompressible(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: h2on2.hh:164
static constexpr int comp1Idx
index of the second component
Definition: h2on2.hh:92
static constexpr int gasPhaseIdx
index of the gas phase
Definition: h2on2.hh:85
static std::string componentName(int compIdx)
Return the human readable name of a component.
Definition: h2on2.hh:198
static constexpr int H2OIdx
Definition: h2on2.hh:89
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: h2on2.hh:340
static void init()
Initialize the fluid system's static parameters generically.
Definition: h2on2.hh:319
static constexpr int comp0Idx
index of the first component
Definition: h2on2.hh:91
static Scalar componentEnthalpy(const FluidState &fluidState, int phaseIdx, int componentIdx)
Returns the specific enthalpy of a component in the specified phase.
Definition: h2on2.hh:688
static Scalar acentricFactor(int compIdx)
The acentric factor of a component .
Definition: h2on2.hh:298
static Scalar fugacityCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Calculate the fugacity coefficient of an individual component in a fluid phase.
Definition: h2on2.hh:548
static constexpr int phase0Idx
index of the first phase
Definition: h2on2.hh:86
static Scalar density(const FluidState &fluidState, int phaseIdx)
Given a phase's composition, temperature, pressure, and the partial pressures of all components,...
Definition: h2on2.hh:371
static Scalar kelvinVaporPressure(const FluidState &fluidState, const int phaseIdx, const int compIdx)
Vapor pressure including the Kelvin equation in .
Definition: h2on2.hh:267
static Scalar criticalPressure(int compIdx)
Critical pressure of a component .
Definition: h2on2.hh:246
static constexpr int numComponents
Number of components in the fluid system.
Definition: h2on2.hh:82
static constexpr int numPhases
Number of phases in the fluid system.
Definition: h2on2.hh:81
static constexpr bool isMiscible()
Returns whether the fluids are miscible.
Definition: h2on2.hh:118
static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
The molar density of a fluid phase in .
Definition: h2on2.hh:424
static Scalar enthalpy(const FluidState &fluidState, int phaseIdx)
Given a phase's composition, temperature, pressure and density, calculate its specific enthalpy .
Definition: h2on2.hh:656
static constexpr int gasCompIdx
index of the gas component
Definition: h2on2.hh:94
static Scalar diffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Calculate the molecular diffusion coefficient for a component in a fluid phase .
Definition: h2on2.hh:595
static Scalar viscosity(const FluidState &fluidState, int phaseIdx)
Calculate the dynamic viscosity of a fluid phase .
Definition: h2on2.hh:467
static bool isIdealGas(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: h2on2.hh:180
static Scalar criticalTemperature(int compIdx)
Critical temperature of a component .
Definition: h2on2.hh:230
static std::string phaseName(int phaseIdx)
Return the human readable name of a fluid phase.
Definition: h2on2.hh:104
Relations valid for an ideal gas.
Definition: idealgas.hh:37
static constexpr Scalar density(Scalar avgMolarMass, Scalar temperature, Scalar pressure)
The density of the gas in , depending on pressure, temperature and average molar mass of the gas.
Definition: idealgas.hh:49
static constexpr Scalar molarDensity(Scalar temperature, Scalar pressure)
The molar density of the gas , depending on pressure and temperature.
Definition: idealgas.hh:70
Fluid system base class.