3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
h2on2o2.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 *****************************************************************************/
19
25#ifndef DUMUX_H2O_N2_O2_FLUID_SYSTEM_HH
26#define DUMUX_H2O_N2_O2_FLUID_SYSTEM_HH
27
28#include <cassert>
29#include <iomanip>
30
33
37
41
46
47#include <dumux/io/name.hh>
48
49namespace Dumux {
50namespace FluidSystems {
55template<bool fastButSimplifiedRelations = false>
57{
58 static constexpr bool useH2ODensityAsLiquidMixtureDensity() { return fastButSimplifiedRelations; }
59 static constexpr bool useIdealGasDensity() { return fastButSimplifiedRelations; }
60 static constexpr bool useN2ViscosityAsGasMixtureViscosity() { return fastButSimplifiedRelations; }
61 static constexpr bool useN2HeatConductivityAsGasMixtureHeatConductivity() { return fastButSimplifiedRelations; }
62 static constexpr bool useIdealGasHeatCapacities() { return fastButSimplifiedRelations; }
63};
64
75template <class Scalar, class Policy = H2ON2O2DefaultPolicy<>>
77 : public Base<Scalar, H2ON2O2<Scalar, Policy> >
78{
81
87
89 using H2O = TabulatedH2O;
90
92 using N2 = SimpleN2;
93
94public:
95 static constexpr int numPhases = 2;
96 static constexpr int numComponents = 3;
97 static constexpr int numSPhases = 0; // TODO: Remove
98
99 static constexpr int liquidPhaseIdx = 0;
100 static constexpr int gasPhaseIdx = 1;
101 static constexpr int phase0Idx = liquidPhaseIdx;
102 static constexpr int phase1Idx = gasPhaseIdx;
103
104 static constexpr int H2OIdx = 0;
105 static constexpr int N2Idx = 1;
106 static constexpr int O2Idx = 2;
107
108 static constexpr int comp0Idx = H2OIdx; // first major component
109 static constexpr int comp1Idx = N2Idx; // second major component
110 static constexpr int comp2Idx = O2Idx; // secondary component
111
112 // main component at 20°C and 1 bar
113 static constexpr int liquidPhaseMainCompIdx = H2OIdx;
114 static constexpr int gasPhaseMainCompIdx = N2Idx;
115
116 /****************************************
117 * Fluid phase related static parameters
118 ****************************************/
124 static std::string phaseName(int phaseIdx)
125 {
126 assert(0 <= phaseIdx && phaseIdx < numPhases);
127 switch (phaseIdx)
128 {
129 case liquidPhaseIdx: return IOName::liquidPhase();
130 case gasPhaseIdx: return IOName::gaseousPhase();
131 }
132 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
133 }
134
140 static constexpr bool isGas(int phaseIdx)
141 {
142 assert(0 <= phaseIdx && phaseIdx < numPhases);
143 return phaseIdx == gasPhaseIdx;
144 }
145
160 static bool isIdealMixture(int phaseIdx)
161 {
162 assert(0 <= phaseIdx && phaseIdx < numPhases);
163 // we assume Henry's and Raoult's laws for the water phase and
164 // and no interaction between gas molecules of different
165 // components, so all phases are ideal mixtures!
166 return true;
167 }
168
178 static constexpr bool isCompressible(int phaseIdx)
179 {
180 assert(0 <= phaseIdx && phaseIdx < numPhases);
181 // gases are always compressible
182 if (phaseIdx == gasPhaseIdx)
183 return true;
184 // the water component decides for the liquid phase...
186 }
187
194 static bool isIdealGas(int phaseIdx)
195 {
196 assert(0 <= phaseIdx && phaseIdx < numPhases);
197 if (phaseIdx == gasPhaseIdx)
198 // let the components decide
200 return false; // not a gas
201 }
202
206 static constexpr bool isMiscible()
207 { return true; }
208
209 /****************************************
210 * Component related static parameters
211 ****************************************/
217 static std::string componentName(int compIdx)
218 {
219 switch (compIdx)
220 {
221 case H2OIdx: return H2O::name();
222 case N2Idx: return N2::name();
223 case O2Idx: return O2::name();
224 }
225
226 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx);
227 }
228
234 static Scalar molarMass(int compIdx)
235 {
236 static const Scalar M[] = {
240 };
241
242 assert(0 <= compIdx && compIdx < numComponents);
243 return M[compIdx];
244 }
245
251 static Scalar criticalTemperature(int compIdx)
252 {
253 static const Scalar Tcrit[] = {
257 };
258
259 assert(0 <= compIdx && compIdx < numComponents);
260 return Tcrit[compIdx];
261 }
262
268 static Scalar criticalPressure(int compIdx)
269 {
270 static const Scalar pcrit[] = {
274 };
275
276 assert(0 <= compIdx && compIdx < numComponents);
277 return pcrit[compIdx];
278 }
279
285 static Scalar criticalMolarVolume(int compIdx)
286 {
287 DUNE_THROW(Dune::NotImplemented,
288 "H2ON2O2FluidSystem::criticalMolarVolume()");
289 }
290
296 static Scalar acentricFactor(int compIdx)
297 {
298 static const Scalar accFac[] = {
299 H2O::acentricFactor(),
300 N2::acentricFactor(),
301 O2::acentricFactor()
302 };
303
304 assert(0 <= compIdx && compIdx < numComponents);
305 return accFac[compIdx];
306 }
307
319 template <class FluidState>
320 static Scalar kelvinVaporPressure(const FluidState &fluidState,
321 const int phaseIdx,
322 const int compIdx,
323 const Scalar radius)
324 {
325 assert(0 <= phaseIdx && phaseIdx == liquidPhaseIdx);
326 assert(0 <= compIdx && compIdx == liquidPhaseMainCompIdx);
327
328 Scalar T = fluidState.temperature(phaseIdx);
329
330 Scalar vaporPressure = H2O::vaporPressure(T);
331 Scalar exponent = molarMass(compIdx)/(density(fluidState, phaseIdx) * Constants::R * T);
332 exponent *= (2 * surfaceTension(fluidState) / radius);
333 using std::exp;
334 Scalar kelvinVaporPressure = vaporPressure * exp(exponent);
335
336 return kelvinVaporPressure;
337 }
338
348 template <class FluidState>
349 static Scalar kelvinVaporPressure(const FluidState &fluidState,
350 const int phaseIdx,
351 const int compIdx)
352 {
353 assert(compIdx == liquidPhaseMainCompIdx && phaseIdx == liquidPhaseIdx);
354
355 using std::exp;
356 return fugacityCoefficient(fluidState, phaseIdx, compIdx)
357 * fluidState.pressure(phaseIdx)
358 * exp(-(fluidState.pressure(gasPhaseIdx)-fluidState.pressure(liquidPhaseIdx))
359 / density(fluidState, phaseIdx)
361 / fluidState.temperature());
362 }
363
371 template <class FluidState>
372 static Scalar surfaceTension(const FluidState &fluidState)
373 {
374 const Scalar T = fluidState.temperature(); //K
375 const Scalar B = 0.2358 ; // [N/m]
376 const Scalar T_c = H2O::criticalTemperature(); //K
377 const Scalar mu = 1.256;
378 const Scalar b = -0.625;
379 //Equation to calculate surface Tension of Water According to IAPWS Release on Surface Tension from September 1994
380 using std::pow;
381 const Scalar surfaceTension = B*pow((1.-(T/T_c)),mu)*(1.+b*(1.-(T/T_c)));
382 return surfaceTension; //surface Tension [N/m]
383 }
384 /****************************************
385 * thermodynamic relations
386 ****************************************/
387
394 static void init()
395 {
396 init(/*tempMin=*/273.15,
397 /*tempMax=*/623.15,
398 /*numTemp=*/100,
399 /*pMin=*/0.0,
400 /*pMax=*/20e6,
401 /*numP=*/200);
402 }
403
415 static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp,
416 Scalar pressMin, Scalar pressMax, unsigned nPress)
417 {
418 std::cout << "The H2O-N2-O2 fluid system was configured with the following policy:\n";
419 std::cout << " - use H2O density as liquid mixture density: " << std::boolalpha << Policy::useH2ODensityAsLiquidMixtureDensity() << "\n";
420 std::cout << " - use ideal gas density: " << std::boolalpha << Policy::useIdealGasDensity() << "\n";
421 std::cout << " - use N2 viscosity as gas mixture viscosity: " << std::boolalpha << Policy::useN2ViscosityAsGasMixtureViscosity() << "\n";
422 std::cout << " - use N2 heat conductivity as gas mixture heat conductivity: " << std::boolalpha << Policy::useN2HeatConductivityAsGasMixtureHeatConductivity() << "\n";
423 std::cout << " - use ideal gas heat capacities: " << std::boolalpha << Policy::useIdealGasHeatCapacities() << std::endl;
424
426 {
427 TabulatedH2O::init(tempMin, tempMax, nTemp,
428 pressMin, pressMax, nPress);
429 }
430 }
431
432 using Base::density;
445 template <class FluidState>
446 static Scalar density(const FluidState &fluidState,
447 int phaseIdx)
448 {
449 assert(0 <= phaseIdx && phaseIdx < numPhases);
450
451 Scalar T = fluidState.temperature(phaseIdx);
452 Scalar p = fluidState.pressure(phaseIdx);
453
454 // liquid phase
455 if (phaseIdx == liquidPhaseIdx)
456 {
457 // assume pure water
458 if (Policy::useH2ODensityAsLiquidMixtureDensity())
459 return H2O::liquidDensity(T, p);
460
461 // See: Eq. (7) in Class et al. (2002a)
462 // This assumes each gas molecule displaces exactly one
463 // molecule in the liquid.
464 else
465 return H2O::liquidMolarDensity(T, p)
466 * (fluidState.moleFraction(liquidPhaseIdx, H2OIdx)*H2O::molarMass()
467 + fluidState.moleFraction(liquidPhaseIdx, N2Idx)*N2::molarMass()
468 + fluidState.moleFraction(liquidPhaseIdx, O2Idx)*O2::molarMass());
469 }
470
471 // gas phase
472 else if (phaseIdx == gasPhaseIdx)
473 {
474
475 // for the gas phase assume an ideal gas
476 using std::max;
477 if (Policy::useIdealGasDensity())
478 return IdealGas::molarDensity(T, p) * fluidState.averageMolarMass(gasPhaseIdx);
479
480 // assume ideal mixture: steam, nitrogen and oxygen don't "see" each other
481 else
482 return H2O::gasDensity(T, fluidState.partialPressure(gasPhaseIdx, H2OIdx))
483 + N2::gasDensity(T, fluidState.partialPressure(gasPhaseIdx, N2Idx))
484 + O2::gasDensity(T, fluidState.partialPressure(gasPhaseIdx, O2Idx));
485 }
486
487 DUNE_THROW(Dune::InvalidStateException, "Unknown phase index " << phaseIdx);
488 }
489
490 using Base::molarDensity;
503 template <class FluidState>
504 static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
505 {
506 const Scalar T = fluidState.temperature(phaseIdx);
507 const Scalar p = fluidState.pressure(phaseIdx);
508
509 if (phaseIdx == liquidPhaseIdx)
510 {
511 // assume pure water or that each gas molecule displaces exactly one
512 // molecule in the liquid.
513 return H2O::liquidMolarDensity(T, p);
514 }
515 else
516 {
517 if (Policy::useIdealGasDensity())
518 { //assume ideal gas
519 return IdealGas::molarDensity(T,p);
520 }
521
522 return H2O::gasMolarDensity(T, fluidState.partialPressure(gasPhaseIdx, H2OIdx))
523 + N2::gasMolarDensity(T, fluidState.partialPressure(gasPhaseIdx, N2Idx))
524 + O2::gasMolarDensity(T, fluidState.partialPressure(gasPhaseIdx, O2Idx));
525 }
526 }
527
528 using Base::viscosity;
541 template <class FluidState>
542 static Scalar viscosity(const FluidState &fluidState,
543 int phaseIdx)
544 {
545 assert(0 <= phaseIdx && phaseIdx < numPhases);
546
547 Scalar T = fluidState.temperature(phaseIdx);
548 Scalar p = fluidState.pressure(phaseIdx);
549
550 // liquid phase
551 if (phaseIdx == liquidPhaseIdx) {
552 // assume pure water for the liquid phase
553 return H2O::liquidViscosity(T, p);
554 }
555
556 // gas phase
557 if (Policy::useN2ViscosityAsGasMixtureViscosity())
558 {
559 // assume pure nitrogen for the gas phase
560 return N2::gasViscosity(T, p);
561 }
562 else
563 {
564 // Wilke method (Reid et al.):
565 Scalar muResult = 0;
566 const Scalar mu[numComponents] = {
568 N2::gasViscosity(T, p),
569 O2::gasViscosity(T, p)
570 };
571
572 Scalar sumx = 0.0;
573 using std::max;
574 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
575 sumx += fluidState.moleFraction(phaseIdx, compIdx);
576 sumx = max(1e-10, sumx);
577
578 for (int i = 0; i < numComponents; ++i) {
579 Scalar divisor = 0;
580 using std::pow;
581 using std::sqrt;
582 for (int j = 0; j < numComponents; ++j) {
583 Scalar phiIJ = 1 + sqrt(mu[i]/mu[j]) * pow(molarMass(j)/molarMass(i), 1/4.0);
584 phiIJ *= phiIJ;
585 phiIJ /= sqrt(8*(1 + molarMass(i)/molarMass(j)));
586 divisor += fluidState.moleFraction(phaseIdx, j)/sumx * phiIJ;
587 }
588 muResult += fluidState.moleFraction(phaseIdx, i)/sumx * mu[i] / divisor;
589 }
590 return muResult;
591 }
592 }
593
613 template <class FluidState>
614 static Scalar fugacityCoefficient(const FluidState &fluidState,
615 int phaseIdx,
616 int compIdx)
617 {
618 assert(0 <= phaseIdx && phaseIdx < numPhases);
619 assert(0 <= compIdx && compIdx < numComponents);
620
621 Scalar T = fluidState.temperature(phaseIdx);
622 Scalar p = fluidState.pressure(phaseIdx);
623
624 // liquid phase
625 if (phaseIdx == liquidPhaseIdx)
626 {
627 switch(compIdx){
628 case H2OIdx: return H2O::vaporPressure(T)/p;
629 case N2Idx: return BinaryCoeff::H2O_N2::henry(T)/p;
630 case O2Idx: return BinaryCoeff::H2O_O2::henry(T)/p;
631 };
632 }
633
634 // for the gas phase, assume an ideal gas when it comes to
635 // fugacity (-> fugacity == partial pressure)
636 return 1.0;
637 }
638
663 template <class FluidState>
664 static Scalar diffusionCoefficient(const FluidState &fluidState,
665 int phaseIdx,
666 int compIdx)
667 {
668 DUNE_THROW(Dune::NotImplemented, "Diffusion coefficients");
669 }
670
682 template <class FluidState>
683 static Scalar binaryDiffusionCoefficient(const FluidState &fluidState,
684 int phaseIdx,
685 int compIIdx,
686 int compJIdx)
687
688 {
689 if (compIIdx > compJIdx)
690 {
691 using std::swap;
692 swap(compIIdx, compJIdx);
693 }
694
695#ifndef NDEBUG
696 if (compIIdx == compJIdx ||
697 phaseIdx > numPhases - 1 ||
698 compJIdx > numComponents - 1)
699 {
700 DUNE_THROW(Dune::InvalidStateException,
701 "Binary diffusion coefficient of components "
702 << compIIdx << " and " << compJIdx
703 << " in phase " << phaseIdx << " is undefined!\n");
704 }
705#endif
706
707 Scalar T = fluidState.temperature(phaseIdx);
708 Scalar p = fluidState.pressure(phaseIdx);
709
710 // liquid phase
711 if (phaseIdx == liquidPhaseIdx) {
712 if (compIIdx == H2OIdx && compJIdx == N2Idx)
714 if (compIIdx == H2OIdx && compJIdx == O2Idx)
716 DUNE_THROW(Dune::InvalidStateException,
717 "Binary diffusion coefficient of components "
718 << compIIdx << " and " << compJIdx
719 << " in phase " << phaseIdx << " is undefined!\n");
720 }
721 // gas phase
722 if (phaseIdx == gasPhaseIdx) {
723 if (compIIdx == H2OIdx && compJIdx == N2Idx)
725 if (compIIdx == H2OIdx && compJIdx == O2Idx)
727 if(compIIdx == N2Idx && compJIdx == O2Idx)
729 DUNE_THROW(Dune::InvalidStateException,
730 "Binary diffusion coefficient of components "
731 << compIIdx << " and " << compJIdx
732 << " in phase " << phaseIdx << " is undefined!\n");
733 }
734
735 DUNE_THROW(Dune::InvalidStateException,
736 "Binary diffusion coefficient of components "
737 << compIIdx << " and " << compJIdx
738 << " in phase " << phaseIdx << " is undefined!\n");
739 }
740
741 using Base::enthalpy;
754 template <class FluidState>
755 static Scalar enthalpy(const FluidState &fluidState,
756 int phaseIdx)
757 {
758 Scalar T = fluidState.temperature(phaseIdx);
759 Scalar p = fluidState.pressure(phaseIdx);
762
763 // liquid phase
764 if (phaseIdx == liquidPhaseIdx) {
765 return H2O::liquidEnthalpy(T, p);
766 }
767 // gas phase
768 else if (phaseIdx == gasPhaseIdx)
769 {
770 // assume ideal mixture: which means
771 // that the total specific enthalpy is the sum of the
772 // "partial specific enthalpies" of the components.
773 Scalar hH2O =
774 fluidState.massFraction(gasPhaseIdx, H2OIdx)
775 * H2O::gasEnthalpy(T, p);
776 Scalar hN2 =
777 fluidState.massFraction(gasPhaseIdx, N2Idx)
778 * N2::gasEnthalpy(T,p);
779 Scalar hO2 =
780 fluidState.massFraction(gasPhaseIdx, O2Idx)
781 * O2::gasEnthalpy(T,p);
782 return hH2O + hN2 + hO2;
783 }
784 else
785 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
786 }
790 template <class FluidState>
791 static Scalar componentEnthalpy(const FluidState &fluidState,
792 int phaseIdx,
793 int componentIdx)
794 {
795 const Scalar T = fluidState.temperature(phaseIdx);
796 const Scalar p = fluidState.pressure(phaseIdx);
797
798 if (phaseIdx == phase0Idx)
799 {
800 if (componentIdx == H2OIdx)
801 return H2O::liquidEnthalpy(T, p);
802 else if (componentIdx == N2Idx)
803 DUNE_THROW(Dune::NotImplemented, "Component enthalpy of nitrogen in liquid phase");
804 else if (componentIdx == O2Idx)
805 DUNE_THROW(Dune::NotImplemented, "Component enthalpy of oxygen in liquid phase");
806 else
807 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx);
808 }
809 else if (phaseIdx == phase1Idx)
810 {
811 if (componentIdx == H2OIdx)
812 return H2O::gasEnthalpy(T, p);
813 else if (componentIdx == N2Idx)
814 return N2::gasEnthalpy(T, p);
815 else if (componentIdx == O2Idx)
816 return O2::gasEnthalpy(T, p);
817 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx);
818 }
819 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
820 }
821
832 template <class FluidState>
833 static Scalar thermalConductivity(const FluidState &fluidState,
834 const int phaseIdx)
835 {
836 assert(0 <= phaseIdx && phaseIdx < numPhases);
837 Scalar temperature = fluidState.temperature(phaseIdx) ;
838 Scalar pressure = fluidState.pressure(phaseIdx);
839
840 if (phaseIdx == liquidPhaseIdx)
841 {
843 }
844 else
845 {
848 if (!Policy::useN2HeatConductivityAsGasMixtureHeatConductivity())
849 {
850 Scalar xN2 = fluidState.moleFraction(phaseIdx, N2Idx);
851 Scalar xO2 = fluidState.moleFraction(phaseIdx, O2Idx);
852 Scalar xH2O = fluidState.moleFraction(phaseIdx, H2OIdx);
853 Scalar lambdaN2 = xN2 * lambdaPureN2;
854 Scalar lambdaO2 = xO2 * lambdaPureO2;
855 Scalar partialPressure = pressure * xH2O;
856 Scalar lambdaH2O = xH2O * H2O::gasThermalConductivity(temperature, partialPressure);
857 return lambdaN2 + lambdaH2O + lambdaO2;
858 }
859 else
860 return lambdaPureN2;
861 }
862 }
863
864 using Base::heatCapacity;
872 template <class FluidState>
873 static Scalar heatCapacity(const FluidState &fluidState,
874 int phaseIdx)
875 {
876 if (phaseIdx == liquidPhaseIdx) {
877 return H2O::liquidHeatCapacity(fluidState.temperature(phaseIdx),
878 fluidState.pressure(phaseIdx));
879 }
880
881 Scalar c_pN2;
882 Scalar c_pO2;
883 Scalar c_pH2O;
884 // let the water and nitrogen components do things their own way
885 if (!Policy::useIdealGasHeatCapacities()) {
886 c_pN2 = N2::gasHeatCapacity(fluidState.temperature(phaseIdx),
887 fluidState.pressure(phaseIdx)
888 * fluidState.moleFraction(phaseIdx, N2Idx));
889
890 c_pH2O = H2O::gasHeatCapacity(fluidState.temperature(phaseIdx),
891 fluidState.pressure(phaseIdx)
892 * fluidState.moleFraction(phaseIdx, H2OIdx));
893 c_pO2 = O2::gasHeatCapacity(fluidState.temperature(phaseIdx),
894 fluidState.pressure(phaseIdx)
895 * fluidState.moleFraction(phaseIdx, O2Idx));
896 }
897 else {
898 // assume an ideal gas for both components. See:
899 //
900 //http://en.wikipedia.org/wiki/Heat_capacity
901 Scalar c_vN2molar = Constants::R*2.39;
902 Scalar c_pN2molar = Constants::R + c_vN2molar;
903
904 Scalar c_vO2molar = Constants::R*2.43;
905 Scalar c_pO2molar = Constants::R + c_vO2molar;
906
907 Scalar c_vH2Omolar = Constants::R*3.37; // <- correct??
908 Scalar c_pH2Omolar = Constants::R + c_vH2Omolar;
909
910 c_pN2 = c_pN2molar/molarMass(N2Idx);
911 c_pO2 = c_pO2molar/molarMass(O2Idx);
912 c_pH2O = c_pH2Omolar/molarMass(H2OIdx);
913 }
914
915 // mangle all components together
916 return
917 c_pH2O*fluidState.massFraction(gasPhaseIdx, H2OIdx)
918 + c_pN2*fluidState.massFraction(gasPhaseIdx, N2Idx)
919 + c_pO2*fluidState.massFraction(gasPhaseIdx, O2Idx);
920 }
921
922};
923
924} // end namespace FluidSystems
925} // end namespace Dumux
926
927#endif
A collection of input/output field names for common physical quantities.
A central place for various physical constants occuring in some equations.
Material properties of pure water .
Properties of pure molecular nitrogen .
Tabulates all thermodynamic properties of a given untabulated chemical species.
Properties of pure molecular oxygen .
Relations valid for an ideal gas.
Binary coefficients for water and nitrogen.
Binary coefficients for nitrogen and oxygen.
Binary coefficients for water and oxygen.
Some exceptions thrown in DuMux
Some templates to wrap the valgrind macros.
bool CheckDefined(const T &value)
Make valgrind complain if the object occupied by an object is undefined.
Definition: valgrind.hh:72
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 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
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
static Scalar gasDiffCoeff(Scalar temperature, Scalar pressure)
Binary diffusion coefficient for molecular water and oxygen.
Definition: h2o_o2.hh:66
static Scalar henry(Scalar temperature)
Henry coefficient for molecular oxygen in liquid water.
Definition: h2o_o2.hh:48
static Scalar liquidDiffCoeff(Scalar temperature, Scalar pressure)
Diffusion coefficient for molecular oxygen in liquid water.
Definition: h2o_o2.hh:97
static Scalar gasDiffCoeff(Scalar temperature, Scalar pressure)
Binary diffusion coefficient for molecular oxygen in liquid nitrogen.
Definition: n2_o2.hh:61
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
Properties of pure molecular oxygen .
Definition: o2.hh:47
static Scalar gasEnthalpy(Scalar temperature, Scalar pressure)
Specific enthalpy of pure oxygen gas.
Definition: o2.hh:173
static Scalar gasMolarDensity(Scalar temperature, Scalar pressure)
The molar density of pure in , depending on pressure and temperature.
Definition: o2.hh:146
static constexpr Scalar criticalTemperature()
Returns the critical temperature in of molecular oxygen.
Definition: o2.hh:66
static std::string name()
A human readable name for the .
Definition: o2.hh:54
static constexpr Scalar molarMass()
The molar mass in of molecular oxygen.
Definition: o2.hh:60
static Scalar gasHeatCapacity(Scalar T, Scalar pressure)
Specific isobaric heat capacity of pure oxygen gas.
Definition: o2.hh:190
static Scalar gasViscosity(Scalar temperature, Scalar pressure)
The dynamic viscosity of at a given pressure and temperature.
Definition: o2.hh:217
static constexpr Scalar criticalPressure()
Returns the critical pressure in of molecular oxygen.
Definition: o2.hh:72
static constexpr bool gasIsIdeal()
Returns true if the gas phase is assumed to be ideal.
Definition: o2.hh:152
static constexpr Scalar gasDensity(Scalar temperature, Scalar pressure)
The density in of pure at a given pressure and temperature.
Definition: o2.hh:134
static constexpr Scalar gasThermalConductivity(Scalar temperature, Scalar pressure)
Thermal conductivity of nitrogen.
Definition: o2.hh:256
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 Scalar gasViscosity(Scalar temperature, Scalar pressure)
The dynamic viscosity of gas.
Definition: tabulatedcomponent.hh:538
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
static constexpr Scalar R
The ideal gas constant .
Definition: constants.hh:44
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-O2 fluid system.
Definition: h2on2o2.hh:57
static constexpr bool useIdealGasDensity()
Definition: h2on2o2.hh:59
static constexpr bool useIdealGasHeatCapacities()
Definition: h2on2o2.hh:62
static constexpr bool useH2ODensityAsLiquidMixtureDensity()
Definition: h2on2o2.hh:58
static constexpr bool useN2ViscosityAsGasMixtureViscosity()
Definition: h2on2o2.hh:60
static constexpr bool useN2HeatConductivityAsGasMixtureHeatConductivity()
Definition: h2on2o2.hh:61
A two-phase (water and air) fluid system with water, nitrogen and oxygen as components.
Definition: h2on2o2.hh:78
static std::string phaseName(int phaseIdx)
Return the human readable name of a fluid phase.
Definition: h2on2o2.hh:124
static constexpr int comp1Idx
Definition: h2on2o2.hh:109
static Scalar criticalTemperature(int compIdx)
Critical temperature of a component .
Definition: h2on2o2.hh:251
static constexpr int numSPhases
Definition: h2on2o2.hh:97
static constexpr bool isGas(int phaseIdx)
Return whether a phase is gaseous.
Definition: h2on2o2.hh:140
static constexpr int N2Idx
Definition: h2on2o2.hh:105
static constexpr int H2OIdx
Definition: h2on2o2.hh:104
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: h2on2o2.hh:683
static Scalar componentEnthalpy(const FluidState &fluidState, int phaseIdx, int componentIdx)
Returns the specific enthalpy of a component in a specific phase.
Definition: h2on2o2.hh:791
static Scalar molarMass(int compIdx)
Return the molar mass of a component in .
Definition: h2on2o2.hh:234
static bool isIdealGas(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: h2on2o2.hh:194
static constexpr int phase1Idx
index of the second phase
Definition: h2on2o2.hh:102
static constexpr bool isMiscible()
Returns whether the fluids are miscible.
Definition: h2on2o2.hh:206
static std::string componentName(int compIdx)
Return the human readable name of a component.
Definition: h2on2o2.hh:217
static constexpr int liquidPhaseIdx
index of the liquid phase
Definition: h2on2o2.hh:99
static Scalar fugacityCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Returns the fugacity coefficient of a component in a phase.
Definition: h2on2o2.hh:614
static constexpr int gasPhaseMainCompIdx
Definition: h2on2o2.hh:114
static constexpr bool isCompressible(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: h2on2o2.hh:178
static Scalar kelvinVaporPressure(const FluidState &fluidState, const int phaseIdx, const int compIdx, const Scalar radius)
Kelvin equation in .
Definition: h2on2o2.hh:320
static Scalar density(const FluidState &fluidState, int phaseIdx)
Given a phase's composition, temperature, pressure, and the partial pressures of all components,...
Definition: h2on2o2.hh:446
static Scalar acentricFactor(int compIdx)
The acentric factor of a component .
Definition: h2on2o2.hh:296
static constexpr int phase0Idx
index of the first phase
Definition: h2on2o2.hh:101
static void init()
Initialize the fluid system's static parameters generically.
Definition: h2on2o2.hh:394
static constexpr int numPhases
Number of phases in the fluid system.
Definition: h2on2o2.hh:95
static constexpr int liquidPhaseMainCompIdx
Definition: h2on2o2.hh:113
static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
The molar density of a fluid phase in .
Definition: h2on2o2.hh:504
static Scalar kelvinVaporPressure(const FluidState &fluidState, const int phaseIdx, const int compIdx)
Vapor pressure including the Kelvin equation in .
Definition: h2on2o2.hh:349
static bool isIdealMixture(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: h2on2o2.hh:160
static constexpr int gasPhaseIdx
index of the gas phase
Definition: h2on2o2.hh:100
static Scalar heatCapacity(const FluidState &fluidState, int phaseIdx)
Specific isobaric heat capacity of a fluid phase. .
Definition: h2on2o2.hh:873
static Scalar thermalConductivity(const FluidState &fluidState, const int phaseIdx)
Thermal conductivity of a fluid phase .
Definition: h2on2o2.hh:833
static Scalar criticalPressure(int compIdx)
Critical pressure of a component .
Definition: h2on2o2.hh:268
static Scalar diffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Calculate the molecular diffusion coefficient for a component in a fluid phase .
Definition: h2on2o2.hh:664
static Scalar viscosity(const FluidState &fluidState, int phaseIdx)
Calculate the dynamic viscosity of a fluid phase .
Definition: h2on2o2.hh:542
static Scalar enthalpy(const FluidState &fluidState, int phaseIdx)
Given a phase's composition, temperature, pressure and density, calculate its specific enthalpy .
Definition: h2on2o2.hh:755
static constexpr int O2Idx
Definition: h2on2o2.hh:106
static constexpr int comp0Idx
Definition: h2on2o2.hh:108
static Scalar surfaceTension(const FluidState &fluidState)
Calculate the surface tension between water and air in , according to IAPWS Release on Surface Tensio...
Definition: h2on2o2.hh:372
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: h2on2o2.hh:415
static constexpr int comp2Idx
Definition: h2on2o2.hh:110
static constexpr int numComponents
Number of components in the fluid system.
Definition: h2on2o2.hh:96
static Scalar criticalMolarVolume(int compIdx)
Molar volume of a component at the critical point .
Definition: h2on2o2.hh:285
Relations valid for an ideal gas.
Definition: idealgas.hh:37
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.