3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
h2oheavyoil.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_HEAVYOIL_FLUID_SYSTEM_HH
25#define DUMUX_H2O_HEAVYOIL_FLUID_SYSTEM_HH
26
31
33
35
36#include <dumux/io/name.hh>
37
38namespace Dumux {
39namespace FluidSystems {
40
46template <class Scalar,
49 : public Base<Scalar, H2OHeavyOil<Scalar, H2OType> >
50{
53
54public:
56 using H2O = H2OType;
57
58
59 static const int numPhases = 3;
60 static const int numComponents = 2;
61
62 static const int wPhaseIdx = 0; // index of the water phase
63 static const int nPhaseIdx = 1; // index of the NAPL phase
64 static const int gPhaseIdx = 2; // index of the gas phase
65
66 static const int H2OIdx = 0;
67 static const int NAPLIdx = 1;
68
69 // export component indices to indicate the main component
70 // of the corresponding phase at atmospheric pressure 1 bar
71 // and room temperature 20°C:
72 static const int wCompIdx = H2OIdx;
73 static const int nCompIdx = NAPLIdx;
74
81 static void init()
82 {
83 init(/*tempMin=*/273.15,
84 /*tempMax=*/623.15,
85 /*numTemp=*/100,
86 /*pMin=*/0.0,
87 /*pMax=*/20e6,
88 /*numP=*/200);
89 }
90
102 static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp,
103 Scalar pressMin, Scalar pressMax, unsigned nPress)
104 {
105 if (H2O::isTabulated)
106 {
107 H2O::init(tempMin, tempMax, nTemp,
108 pressMin, pressMax, nPress);
109 }
110 }
111
116 static constexpr int getMainComponent(int phaseIdx)
117 {
118 // For the gas phase, choosing a main component appears to be
119 // rather arbitrary. Motivated by the fact that the thermal conductivity
120 // of the gas phase is set to the thermal conductivity of pure water,
121 // water is chosen for now.
122 if (phaseIdx == nPhaseIdx)
123 return nCompIdx;
124 else
125 return wCompIdx;
126 }
127
131 static constexpr bool isMiscible()
132 { return true; }
133
139 static constexpr bool isGas(int phaseIdx)
140 {
141 assert(0 <= phaseIdx && phaseIdx < numPhases);
142 return phaseIdx == gPhaseIdx;
143 }
144
151 static bool isIdealGas(int phaseIdx)
152 { return phaseIdx == gPhaseIdx && H2O::gasIsIdeal() && HeavyOil::gasIsIdeal(); }
153
168 static bool isIdealMixture(int phaseIdx)
169 {
170 assert(0 <= phaseIdx && phaseIdx < numPhases);
171 // we assume Henry's and Raoult's laws for the water phase and
172 // and no interaction between gas molecules of different
173 // components, so all phases are ideal mixtures!
174 return true;
175 }
176
186 static constexpr bool isCompressible(int phaseIdx)
187 {
188 assert(0 <= phaseIdx && phaseIdx < numPhases);
189 // gases are always compressible
190 if (phaseIdx == gPhaseIdx)
191 return true;
192 else if (phaseIdx == wPhaseIdx)
193 // the water component decides for the water phase...
194 return H2O::liquidIsCompressible();
195
196 // the NAPL component decides for the napl phase...
198 }
199
203 static std::string phaseName(int phaseIdx)
204 {
205 assert(0 <= phaseIdx && phaseIdx < numPhases);
206 switch (phaseIdx)
207 {
208 case wPhaseIdx: return IOName::aqueousPhase();
209 case nPhaseIdx: return IOName::naplPhase();
210 case gPhaseIdx: return IOName::gaseousPhase();
211 }
212 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
213 }
214
218 static std::string componentName(int compIdx)
219 {
220 switch (compIdx) {
221 case H2OIdx: return H2O::name();
222 case NAPLIdx: return HeavyOil::name();
223 };
224 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx);
225 }
226
230 static Scalar molarMass(int compIdx)
231 {
232 switch (compIdx) {
233 case H2OIdx: return H2O::molarMass();
234 case NAPLIdx: return HeavyOil::molarMass();
235 };
236 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx);
237 }
238
243 using Base::density;
244 template <class FluidState>
245 static Scalar density(const FluidState &fluidState, int phaseIdx)
246 {
247 if (phaseIdx == wPhaseIdx) {
248 // See: doctoral thesis of Steffen Ochs 2007
249 // Steam injection into saturated porous media : process analysis including experimental and numerical investigations
250 // http://elib.uni-stuttgart.de/bitstream/11682/271/1/Diss_Ochs_OPUS.pdf
251
252 // This assumes each gas molecule displaces exactly one
253 // molecule in the liquid.
254
255 return H2O::liquidMolarDensity(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx))
256 * (H2O::molarMass()*fluidState.moleFraction(wPhaseIdx, H2OIdx)
257 + HeavyOil::molarMass()*fluidState.moleFraction(wPhaseIdx, NAPLIdx));
258 }
259 else if (phaseIdx == nPhaseIdx) {
260 // assume pure NAPL for the NAPL phase
261 Scalar pressure = HeavyOil::liquidIsCompressible()?fluidState.pressure(phaseIdx):1e100;
262 return HeavyOil::liquidDensity(fluidState.temperature(phaseIdx), pressure);
263 }
264
265 assert (phaseIdx == gPhaseIdx);
266 Scalar pH2O =
267 fluidState.moleFraction(gPhaseIdx, H2OIdx) *
268 fluidState.pressure(gPhaseIdx);
269 Scalar pNAPL =
270 fluidState.moleFraction(gPhaseIdx, NAPLIdx) *
271 fluidState.pressure(gPhaseIdx);
272 return H2O::gasDensity(fluidState.temperature(phaseIdx), pH2O)
273 + HeavyOil::gasDensity(fluidState.temperature(phaseIdx), pNAPL);
274 }
275
276 using Base::molarDensity;
286 template <class FluidState>
287 static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
288 {
289 Scalar temperature = fluidState.temperature(phaseIdx);
290 Scalar pressure = fluidState.pressure(phaseIdx);
291 if (phaseIdx == nPhaseIdx)
292 {
294 }
295 else if (phaseIdx == wPhaseIdx)
296 { // This assumes each gas molecule displaces exactly one
297 // molecule in the liquid.
298 return H2O::liquidMolarDensity(temperature, pressure);
299 }
300 else
301 {
302 return H2O::gasMolarDensity(temperature, fluidState.partialPressure(gPhaseIdx, H2OIdx))
303 + HeavyOil::gasMolarDensity(temperature, fluidState.partialPressure(gPhaseIdx, NAPLIdx));
304 }
305 }
306
310 using Base::viscosity;
311 template <class FluidState>
312 static Scalar viscosity(const FluidState &fluidState,
313 int phaseIdx)
314 {
315 if (phaseIdx == wPhaseIdx) {
316 // assume pure water viscosity
317 return H2O::liquidViscosity(fluidState.temperature(phaseIdx),
318 fluidState.pressure(phaseIdx));
319 }
320 else if (phaseIdx == nPhaseIdx) {
321 // assume pure NAPL viscosity
322 return HeavyOil::liquidViscosity(fluidState.temperature(phaseIdx),
323 fluidState.pressure(phaseIdx));
324 }
325
326 assert (phaseIdx == gPhaseIdx);
327
328 /* Wilke method. See:
329 *
330 * See: R. Reid, et al.: The Properties of Gases and Liquids,
331 * 4th edition, McGraw-Hill, 1987, 407-410
332 * 5th edition, McGraw-Hill, 20001, p. 9.21/22
333 *
334 * in this case, we use a simplified version in order to avoid
335 * computationally costly evaluation of sqrt and pow functions and
336 * divisions
337 * -- compare e.g. with Promo Class p. 32/33
338 */
339 const Scalar mu[numComponents] = {
340 h2oGasViscosityInMixture(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx)),
341 HeavyOil::gasViscosity(fluidState.temperature(phaseIdx), HeavyOil::vaporPressure(fluidState.temperature(phaseIdx)))
342 };
343
344 return mu[H2OIdx]*fluidState.moleFraction(gPhaseIdx, H2OIdx)
345 + mu[NAPLIdx]*fluidState.moleFraction(gPhaseIdx, NAPLIdx);
346 }
347
363 template <class FluidState>
364 static Scalar binaryDiffusionCoefficient(const FluidState &fluidState,
365 int phaseIdx,
366 int compIIdx,
367 int compJIdx)
368
369 {
370 const Scalar T = fluidState.temperature(phaseIdx);
371 const Scalar p = fluidState.pressure(phaseIdx);
372
373 // liquid phase
374 if (phaseIdx == wPhaseIdx)
376
377 // gas phase
378 else if (phaseIdx == gPhaseIdx)
380 else
381 DUNE_THROW(Dune::InvalidStateException,
382 "Non-existent binary diffusion coefficient for phase index "
383 << phaseIdx);
384
385 }
386
408 template <class FluidState>
409 static Scalar diffusionCoefficient(const FluidState &fluidState, int phaseIdx)
410 {
411 // liquid phase
412 if (phaseIdx == wPhaseIdx)
413 return binaryDiffusionCoefficient(fluidState, phaseIdx, H2OIdx, NAPLIdx);
414 // gas phase
415 else if (phaseIdx == gPhaseIdx)
416 return binaryDiffusionCoefficient(fluidState, phaseIdx, NAPLIdx, H2OIdx);
417 else
418 DUNE_THROW(Dune::InvalidStateException,
419 "Non-existent diffusion coefficient for phase index "<< phaseIdx);
420 }
421
425 template <class FluidState>
426 static Scalar henryCoefficient(const FluidState &fluidState,
427 int phaseIdx,
428 int compIdx)
429 {
430 assert(0 <= phaseIdx && phaseIdx < numPhases);
431 assert(0 <= compIdx && compIdx < numComponents);
432
433 const Scalar T = fluidState.temperature(phaseIdx);
434
435 if (compIdx == NAPLIdx && phaseIdx == wPhaseIdx)
437
438 else if (phaseIdx == nPhaseIdx && compIdx == H2OIdx)
440
441 else
442 DUNE_THROW(Dune::InvalidStateException, "non-existent henry coefficient for phase index " << phaseIdx
443 << " and component index " << compIdx);
444 }
445
449 template <class FluidState>
450 static Scalar partialPressureGas(const FluidState &fluidState, int phaseIdx,
451 int compIdx)
452 {
453 assert(0 <= compIdx && compIdx < numComponents);
454
455 const Scalar T = fluidState.temperature(phaseIdx);
456 if (compIdx == NAPLIdx)
457 return HeavyOil::vaporPressure(T);
458 else if (compIdx == H2OIdx)
459 return H2O::vaporPressure(T);
460 else
461 DUNE_THROW(Dune::InvalidStateException, "non-existent component index " << compIdx);
462 }
463
467 template <class FluidState>
468 static Scalar inverseVaporPressureCurve(const FluidState &fluidState,
469 int phaseIdx,
470 int compIdx)
471 {
472 assert(0 <= compIdx && compIdx < numComponents);
473
474 const Scalar pressure = fluidState.pressure(phaseIdx);
475 if (compIdx == NAPLIdx)
477 else if (compIdx == H2OIdx)
478 return H2O::vaporTemperature(pressure);
479 else
480 DUNE_THROW(Dune::InvalidStateException, "non-existent component index " << compIdx);
481 }
482
483
484
492 using Base::enthalpy;
493 template <class FluidState>
494 static Scalar enthalpy(const FluidState &fluidState,
495 int phaseIdx)
496 {
497 if (phaseIdx == wPhaseIdx) {
498 return H2O::liquidEnthalpy(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
499 }
500 else if (phaseIdx == nPhaseIdx) {
501 return HeavyOil::liquidEnthalpy(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
502 }
503 else if (phaseIdx == gPhaseIdx) { // gas phase enthalpy depends strongly on composition
504 Scalar hgc = HeavyOil::gasEnthalpy(fluidState.temperature(phaseIdx),
505 fluidState.pressure(phaseIdx));
506 Scalar hgw = H2O::gasEnthalpy(fluidState.temperature(phaseIdx),
507 fluidState.pressure(phaseIdx));
508
509 Scalar result = 0;
510 result += hgw * fluidState.massFraction(gPhaseIdx, H2OIdx);
511 result += hgc * fluidState.massFraction(gPhaseIdx, NAPLIdx);
512
513 return result;
514 }
515 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
516 }
517
524 template <class FluidState>
525 static Scalar componentEnthalpy(const FluidState& fluidState, int phaseIdx, int componentIdx)
526 {
527 const Scalar T = fluidState.temperature(phaseIdx);
528 const Scalar p = fluidState.pressure(phaseIdx);
529
530 if (phaseIdx == wPhaseIdx)
531 {
532 if (componentIdx == H2OIdx)
533 return H2O::liquidEnthalpy(T, p);
534 else if (componentIdx == NAPLIdx)
535 DUNE_THROW(Dune::NotImplemented, "The component enthalpy for NAPL in water is not implemented.");
536 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx);
537 }
538 else if (phaseIdx == nPhaseIdx)
539 {
540 if (componentIdx == H2OIdx)
541 DUNE_THROW(Dune::NotImplemented, "The component enthalpy for water in NAPL is not implemented.");
542 else if (componentIdx == NAPLIdx)
543 return HeavyOil::liquidEnthalpy(T, p);
544 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx);
545 }
546 else if (phaseIdx == gPhaseIdx)
547 {
548 if (componentIdx == H2OIdx)
549 return H2O::gasEnthalpy(T, p);
550 else if (componentIdx == NAPLIdx)
551 return HeavyOil::gasEnthalpy(T, p);
552 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx);
553 }
554 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
555 }
556
557 using Base::heatCapacity;
558 template <class FluidState>
559 static Scalar heatCapacity(const FluidState &fluidState,
560 int phaseIdx)
561 {
562 DUNE_THROW(Dune::NotImplemented, "FluidSystems::H2ONAPL::heatCapacity()");
563 }
564
574 template <class FluidState>
575 static Scalar thermalConductivity(const FluidState &fluidState,
576 int phaseIdx)
577 {
578 const Scalar temperature = fluidState.temperature(phaseIdx) ;
579 const Scalar pressure = fluidState.pressure(phaseIdx);
580 if (phaseIdx == wPhaseIdx)
581 {
582 return H2O::liquidThermalConductivity(temperature, pressure);
583 }
584 else if (phaseIdx == nPhaseIdx)
585 {
587 }
588 else if (phaseIdx == gPhaseIdx)
589 {
590 return H2O::gasThermalConductivity(temperature, pressure);
591 }
592 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
593 }
594
595};
596} // end namespace FluidSystems
597} // end namespace Dumux
598
599#endif
A collection of input/output field names for common physical quantities.
Binary coefficients for water and heavy oil.
Material properties of pure water .
Properties of the component heavyoil.
Tabulates all thermodynamic properties of a given untabulated chemical species.
Relations valid for an ideal gas.
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
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 [m^2/s] for heavy oil in liquid water.
Definition: h2o_heavyoil.hh:90
static Scalar gasDiffCoeff(Scalar temperature, Scalar pressure)
Binary diffusion coefficient [m^2/s] for molecular water and heavy oil.
Definition: h2o_heavyoil.hh:77
static Scalar henryOilInWater(Scalar temperature)
Henry coefficient for heavy oil in liquid water.
Definition: h2o_heavyoil.hh:47
static Scalar henryWaterInOil(Scalar temperature)
Henry coefficient for water in liquid heavy oil.
Definition: h2o_heavyoil.hh:62
Properties of the component heavyoil.
Definition: heavyoil.hh:50
static Scalar liquidViscosity(Scalar temperature, Scalar pressure)
The dynamic viscosity of pure heavyoil.
Definition: heavyoil.hh:434
static constexpr Scalar molarMass()
The molar mass in of heavyoil.
Definition: heavyoil.hh:63
static Scalar vaporPressure(Scalar temperature)
The saturation vapor pressure in of.
Definition: heavyoil.hh:231
static constexpr bool gasIsIdeal()
Returns true if the gas phase is assumed to be ideal.
Definition: heavyoil.hh:390
static Scalar liquidMolarDensity(Scalar temperature, Scalar pressure)
The molar density of pure heavyoil in at a given pressure and temperature.
Definition: heavyoil.hh:378
static Scalar gasMolarDensity(Scalar temperature, Scalar pressure)
The molar density of pure heavyoil in , depending on pressure and temperature.
Definition: heavyoil.hh:351
static Scalar gasViscosity(Scalar temperature, Scalar pressure)
The dynamic viscosity of heavyoil vapor.
Definition: heavyoil.hh:405
static Scalar vaporTemperature(Scalar pressure)
Definition: heavyoil.hh:243
static Scalar liquidThermalConductivity(Scalar temperature, Scalar pressure)
Thermal conductivity of heavy oil.
Definition: heavyoil.hh:475
static Scalar gasDensity(Scalar temperature, Scalar pressure)
The (ideal) gas density of heavyoil vapor at a given temperature and pressure .
Definition: heavyoil.hh:340
static constexpr bool liquidIsCompressible()
Returns true if the liquid phase is assumed to be compressible.
Definition: heavyoil.hh:396
static Scalar liquidEnthalpy(const Scalar temperature, const Scalar pressure)
Specific enthalpy of liquid heavyoil .
Definition: heavyoil.hh:261
static std::string name()
A human readable name for heavyoil.
Definition: heavyoil.hh:57
static Scalar gasEnthalpy(Scalar temperature, Scalar pressure)
Specific enthalpy of heavyoil vapor .
Definition: heavyoil.hh:329
static Scalar liquidDensity(Scalar temperature, Scalar pressure)
The density of pure heavyoil at a given pressure and temperature .
Definition: heavyoil.hh:360
Tabulates all thermodynamic properties of a given untabulated chemical species.
Definition: tabulatedcomponent.hh:82
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 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 compositional fluid system with water and heavy oil components in both the liquid and the gas phase...
Definition: h2oheavyoil.hh:50
static Scalar henryCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Henry coefficients of a component in a phase.
Definition: h2oheavyoil.hh:426
static constexpr bool isCompressible(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: h2oheavyoil.hh:186
static const int gPhaseIdx
Definition: h2oheavyoil.hh:64
static std::string componentName(int compIdx)
Return the human readable name of a component (used in indices)
Definition: h2oheavyoil.hh:218
static const int wCompIdx
Definition: h2oheavyoil.hh:72
static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
The molar density of a fluid phase in .
Definition: h2oheavyoil.hh:287
static bool isIdealMixture(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: h2oheavyoil.hh:168
static const int numPhases
Definition: h2oheavyoil.hh:59
static constexpr bool isMiscible()
Returns whether the fluids are miscible.
Definition: h2oheavyoil.hh:131
static Scalar viscosity(const FluidState &fluidState, int phaseIdx)
Definition: h2oheavyoil.hh:312
static Scalar enthalpy(const FluidState &fluidState, int phaseIdx)
Definition: h2oheavyoil.hh:494
static const int H2OIdx
Definition: h2oheavyoil.hh:66
static Scalar partialPressureGas(const FluidState &fluidState, int phaseIdx, int compIdx)
Partial pressures in the gas phase, taken from saturation vapor pressures.
Definition: h2oheavyoil.hh:450
static Scalar componentEnthalpy(const FluidState &fluidState, int phaseIdx, int componentIdx)
Returns the specific enthalpy of a component in a specific phase.
Definition: h2oheavyoil.hh:525
static const int nPhaseIdx
Definition: h2oheavyoil.hh:63
static Scalar thermalConductivity(const FluidState &fluidState, int phaseIdx)
Definition: h2oheavyoil.hh:575
H2OType H2O
Definition: h2oheavyoil.hh:56
static std::string phaseName(int phaseIdx)
Return the human readable name of a phase (used in indices)
Definition: h2oheavyoil.hh:203
static const int numComponents
Definition: h2oheavyoil.hh:60
static constexpr int getMainComponent(int phaseIdx)
Get the main component of a given phase.
Definition: h2oheavyoil.hh:116
static constexpr bool isGas(int phaseIdx)
Return whether a phase is gaseous.
Definition: h2oheavyoil.hh:139
static const int NAPLIdx
Definition: h2oheavyoil.hh:67
static Scalar molarMass(int compIdx)
Return the molar mass of a component in [kg/mol].
Definition: h2oheavyoil.hh:230
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: h2oheavyoil.hh:102
static Scalar density(const FluidState &fluidState, int phaseIdx)
Definition: h2oheavyoil.hh:245
static Scalar inverseVaporPressureCurve(const FluidState &fluidState, int phaseIdx, int compIdx)
Inverse vapor pressures, taken from inverse saturation vapor pressures.
Definition: h2oheavyoil.hh:468
static const int nCompIdx
Definition: h2oheavyoil.hh:73
static bool isIdealGas(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: h2oheavyoil.hh:151
static Scalar diffusionCoefficient(const FluidState &fluidState, int phaseIdx)
Calculate the binary molecular diffusion coefficient for a component in a fluid phase .
Definition: h2oheavyoil.hh:409
static const int wPhaseIdx
Definition: h2oheavyoil.hh:62
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: h2oheavyoil.hh:364
static void init()
Initialize the fluid system's static parameters generically.
Definition: h2oheavyoil.hh:81
static Scalar heatCapacity(const FluidState &fluidState, int phaseIdx)
Definition: h2oheavyoil.hh:559
Fluid system base class.