12#ifndef DUMUX_SPE5_FLUID_SYSTEM_HH
13#define DUMUX_SPE5_FLUID_SYSTEM_HH
22namespace FluidSystems {
43template <
class Scalar>
77 assert(0 <= phaseIdx && phaseIdx <
numPhases);
84 DUNE_THROW(Dune::InvalidStateException,
"Invalid phase index " << phaseIdx);
97 static constexpr bool isGas(
int phaseIdx)
99 assert(0 <= phaseIdx && phaseIdx <
numPhases);
111 assert(0 <= phaseIdx && phaseIdx <
numPhases);
144 assert(0 <= phaseIdx && phaseIdx <
numPhases);
169 static std::string name[] = {
180 return name[compIdx];
189 static const Scalar M[] = {
209 static const Scalar Tcrit[] = {
220 return Tcrit[compIdx];
229 static const Scalar pcrit[] = {
240 return pcrit[compIdx];
249 static const Scalar R = 8.314472;
250 static const Scalar vcrit[] = {
261 return vcrit[compIdx];
270 static const Scalar accFac[] = {
281 return accFac[compIdx];
295 int i = min(comp1Idx, comp2Idx);
296 int j = max(comp1Idx, comp2Idx);
320 Scalar minT = 273.15, maxT = 373.15;
321 Scalar minP = 1e4, maxP = 100e6;
323 Scalar minA = 1e100, maxA = -1e100;
324 Scalar minB = 1e100, maxB = -1e100;
330 minA = min(prParams.
pureParams(compIdx).
a(), minA);
331 maxA = max(prParams.
pureParams(compIdx).
a(), maxA);
332 minB = min(prParams.
pureParams(compIdx).
b(), minB);
333 maxB = max(prParams.
pureParams(compIdx).
b(), maxB);
338 minA = min(prParams.
pureParams(compIdx).
a(), minA);
339 maxA = max(prParams.
pureParams(compIdx).
a(), maxA);
340 minB = min(prParams.
pureParams(compIdx).
b(), minB);
341 maxB = max(prParams.
pureParams(compIdx).
b(), maxB);
346 minA = min(prParams.
pureParams(compIdx).
a(), minA);
347 maxA = max(prParams.
pureParams(compIdx).
a(), maxA);
348 minB = min(prParams.
pureParams(compIdx).
b(), minB);
349 maxB = max(prParams.
pureParams(compIdx).
b(), maxB);
354 minA = min(prParams.
pureParams(compIdx).
a(), minA);
355 maxA = max(prParams.
pureParams(compIdx).
a(), maxA);
356 minB = min(prParams.
pureParams(compIdx).
b(), minB);
357 maxB = max(prParams.
pureParams(compIdx).
b(), maxB);
372 template <
class Flu
idState>
373 static Scalar
density(
const FluidState &fluidState,
377 assert(0 <= phaseIdx && phaseIdx <
numPhases);
379 return fluidState.averageMolarMass(phaseIdx)/paramCache.
molarVolume(phaseIdx);
391 template <
class Flu
idState>
403 template <
class Flu
idState>
408 assert(0 <= phaseIdx && phaseIdx <=
numPhases);
413 return 0.0170e-2 * 0.1;
422 return 0.208e-2 * 0.1;
441 template <
class Flu
idState>
447 assert(0 <= phaseIdx && phaseIdx <=
numPhases);
457 return henryCoeffWater_(compIdx, fs.temperature(
wPhaseIdx))
486 template <
class Flu
idState>
491 { DUNE_THROW(Dune::NotImplemented,
"Diffusion coefficients"); }
504 template <
class Flu
idState>
510 { DUNE_THROW(Dune::NotImplemented,
"Binary diffusion coefficients"); }
519 template <
class Flu
idState>
523 { DUNE_THROW(Dune::NotImplemented,
"Enthalpies"); }
532 template <
class Flu
idState>
536 { DUNE_THROW(Dune::NotImplemented,
"Thermal conductivities"); }
545 template <
class Flu
idState>
549 { DUNE_THROW(Dune::NotImplemented,
"Heat capacities"); }
553 static Scalar henryCoeffWater_(
int compIdx, Scalar
temperature)
562 case C1Idx:
return 5.57601e+09;
563 case C3Idx:
return 1e10;
564 case C6Idx:
return 1e10;
568 default: DUNE_THROW(Dune::InvalidStateException,
"Unknown component index " << compIdx);
Material properties of pure water .
Definition: h2o.hh:49
static Scalar vaporPressure(Scalar T)
The vapor pressure in of pure water at a given temperature.
Definition: h2o.hh:119
static constexpr Scalar criticalTemperature()
Returns the critical temperature of water.
Definition: h2o.hh:81
static std::string name()
A human readable name for the water.
Definition: h2o.hh:63
static constexpr Scalar molarMass()
The molar mass in of water.
Definition: h2o.hh:69
static constexpr Scalar criticalPressure()
Returns the critical pressure of water.
Definition: h2o.hh:87
static constexpr Scalar criticalMolarVolume()
Returns the molar volume of water at the critical point.
Definition: h2o.hh:93
The fluid system for the SPE-5 benchmark problem.
Definition: spe5.hh:45
static Scalar molarDensity(const FluidState &fluidState, const ParameterCache ¶mCache, int phaseIdx)
The molar density of a fluid phase in .
Definition: spe5.hh:392
static bool isIdealMixture(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: spe5.hh:128
static Scalar criticalTemperature(int compIdx)
Critical temperature of a component .
Definition: spe5.hh:207
static const int C15Idx
Definition: spe5.hh:160
static Scalar binaryDiffusionCoefficient(const FluidState &fluidState, const ParameterCache ¶mCache, int phaseIdx, int compIIdx, int compJIdx)
Given a phase's composition, temperature and pressure, return the binary diffusion coefficient for c...
Definition: spe5.hh:505
static Scalar heatCapacity(const FluidState &fluidState, const ParameterCache ¶mCache, int phaseIdx)
Given a phase's composition, temperature and pressure, calculate its heat capacity .
Definition: spe5.hh:546
static const int C3Idx
Definition: spe5.hh:157
static Scalar fugacityCoefficient(const FluidState &fs, const ParameterCache ¶mCache, int phaseIdx, int compIdx)
Calculate the fugacity coefficient of an individual component in a fluid phase.
Definition: spe5.hh:442
static Scalar thermalConductivity(const FluidState &fluidState, const ParameterCache ¶mCache, int phaseIdx)
Given a phase's composition, temperature and pressure, calculate its thermal conductivity .
Definition: spe5.hh:533
static const int H2OIdx
Definition: spe5.hh:155
static std::string componentName(int compIdx)
Return the human readable name of a component.
Definition: spe5.hh:167
static Scalar diffusionCoefficient(const FluidState &fs, const ParameterCache ¶mCache, int phaseIdx, int compIdx)
Calculate the binary molecular diffusion coefficient for a component in a fluid phase .
Definition: spe5.hh:487
static const int C6Idx
Definition: spe5.hh:158
static const int C20Idx
Definition: spe5.hh:161
static Scalar molarMass(int compIdx)
Return the molar mass of a component in .
Definition: spe5.hh:187
static constexpr bool isMiscible()
Returns whether the fluids are miscible.
Definition: spe5.hh:90
static bool isIdealGas(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: spe5.hh:142
static const int C1Idx
Definition: spe5.hh:156
static const int numComponents
Number of components in the fluid system.
Definition: spe5.hh:153
static const int oPhaseIdx
Index of the oil phase.
Definition: spe5.hh:66
static Scalar density(const FluidState &fluidState, const ParameterCache ¶mCache, int phaseIdx)
Calculate the density of a fluid phase.
Definition: spe5.hh:373
static const int C10Idx
Definition: spe5.hh:159
static Scalar acentricFactor(int compIdx)
The acentric factor of a component .
Definition: spe5.hh:268
static const int gPhaseIdx
Index of the gas phase.
Definition: spe5.hh:62
static Scalar viscosity(const FluidState &fs, const ParameterCache ¶mCache, int phaseIdx)
Calculate the dynamic viscosity of a fluid phase .
Definition: spe5.hh:404
static void init()
Initialize the fluid system's static parameters.
Definition: spe5.hh:311
static Scalar interactionCoefficient(int comp1Idx, int comp2Idx)
Returns the interaction coefficient for two components.
Definition: spe5.hh:291
static const int numPhases
Number of phases in the fluid system.
Definition: spe5.hh:59
static std::string phaseName(int phaseIdx)
Return the human readable name of a fluid phase.
Definition: spe5.hh:75
static const int wPhaseIdx
Index of the water phase.
Definition: spe5.hh:64
static constexpr bool isGas(int phaseIdx)
Return whether a phase is gaseous.
Definition: spe5.hh:97
static Scalar criticalPressure(int compIdx)
Critical pressure of a component .
Definition: spe5.hh:227
static Scalar criticalMolarVolume(int compIdx)
Molar volume of a component at the critical point .
Definition: spe5.hh:247
static Scalar enthalpy(const FluidState &fs, const ParameterCache ¶mCache, int phaseIdx)
Given a phase's composition, temperature and pressure, calculate its specific enthalpy .
Definition: spe5.hh:520
static constexpr bool isCompressible(int phaseIdx)
Return whether a phase is compressible.
Definition: spe5.hh:109
Implements the Peng-Robinson equation of state for liquids and gases.
Definition: pengrobinson.hh:48
static void init(Scalar aMin, Scalar aMax, int na, Scalar bMin, Scalar bMax, int nb)
Definition: pengrobinson.hh:53
Implements the Peng-Robinson equation of state for a mixture.
Definition: pengrobinsonmixture.hh:29
static Scalar computeFugacityCoefficient(const FluidState &fs, const Params ¶ms, int phaseIdx, int compIdx)
Returns the fugacity coefficient of an individual component in the phase.
Definition: pengrobinsonmixture.hh:61
Scalar a() const
Returns the attractive parameter 'a' of the Peng-Robinson fluid.
Definition: pengrobinsonparams.hh:38
Scalar b() const
Returns the repulsive parameter 'b' of the Peng-Robinson fluid.
Definition: pengrobinsonparams.hh:46
The mixing rule for the oil and the gas phases of the SPE5 problem.
Definition: pengrobinsonparamsmixture.hh:51
const PureParams & pureParams(int compIdx) const
Return the Peng-Robinson parameters of a pure substance,.
Definition: pengrobinsonparamsmixture.hh:180
void updatePure(const FluidState &fluidState)
Update Peng-Robinson parameters for the pure components.
Definition: pengrobinsonparamsmixture.hh:64
Specifies the parameters required by the SPE5 problem which are despondent on the thermodynamic state...
Definition: spe5parametercache.hh:33
Scalar molarVolume(int phaseIdx) const
Returns the molar volume of a phase .
Definition: spe5parametercache.hh:186
A collection of input/output field names for common physical quantities.
std::string temperature() noexcept
I/O name of temperature for equilibrium models.
Definition: name.hh:39
std::string gaseousPhase() noexcept
I/O name of gaseous phase.
Definition: name.hh:111
std::string naplPhase() noexcept
I/O name of napl phase.
Definition: name.hh:119
std::string aqueousPhase() noexcept
I/O name of aqueous phase.
Definition: name.hh:115
Implements the Peng-Robinson equation of state for a mixture.
Specifies the parameters required by the SPE5 problem which are despondent on the thermodynamic state...
Provides 3rd order polynomial splines.