24#ifndef DUMUX_SPE5_FLUID_SYSTEM_HH
25#define DUMUX_SPE5_FLUID_SYSTEM_HH
34namespace FluidSystems {
55template <
class Scalar>
89 assert(0 <= phaseIdx && phaseIdx <
numPhases);
96 DUNE_THROW(Dune::InvalidStateException,
"Invalid phase index " << phaseIdx);
109 static constexpr bool isGas(
int phaseIdx)
111 assert(0 <= phaseIdx && phaseIdx <
numPhases);
123 assert(0 <= phaseIdx && phaseIdx <
numPhases);
156 assert(0 <= phaseIdx && phaseIdx <
numPhases);
181 static std::string name[] = {
192 return name[compIdx];
201 static const Scalar M[] = {
221 static const Scalar Tcrit[] = {
232 return Tcrit[compIdx];
241 static const Scalar pcrit[] = {
252 return pcrit[compIdx];
261 static const Scalar R = 8.314472;
262 static const Scalar vcrit[] = {
273 return vcrit[compIdx];
282 static const Scalar accFac[] = {
293 return accFac[compIdx];
307 int i = min(comp1Idx, comp2Idx);
308 int j = max(comp1Idx, comp2Idx);
332 Scalar minT = 273.15, maxT = 373.15;
333 Scalar minP = 1e4, maxP = 100e6;
335 Scalar minA = 1e100, maxA = -1e100;
336 Scalar minB = 1e100, maxB = -1e100;
342 minA = min(prParams.
pureParams(compIdx).
a(), minA);
343 maxA = max(prParams.
pureParams(compIdx).
a(), maxA);
344 minB = min(prParams.
pureParams(compIdx).
b(), minB);
345 maxB = max(prParams.
pureParams(compIdx).
b(), maxB);
350 minA = min(prParams.
pureParams(compIdx).
a(), minA);
351 maxA = max(prParams.
pureParams(compIdx).
a(), maxA);
352 minB = min(prParams.
pureParams(compIdx).
b(), minB);
353 maxB = max(prParams.
pureParams(compIdx).
b(), maxB);
358 minA = min(prParams.
pureParams(compIdx).
a(), minA);
359 maxA = max(prParams.
pureParams(compIdx).
a(), maxA);
360 minB = min(prParams.
pureParams(compIdx).
b(), minB);
361 maxB = max(prParams.
pureParams(compIdx).
b(), maxB);
366 minA = min(prParams.
pureParams(compIdx).
a(), minA);
367 maxA = max(prParams.
pureParams(compIdx).
a(), maxA);
368 minB = min(prParams.
pureParams(compIdx).
b(), minB);
369 maxB = max(prParams.
pureParams(compIdx).
b(), maxB);
384 template <
class Flu
idState>
385 static Scalar
density(
const FluidState &fluidState,
389 assert(0 <= phaseIdx && phaseIdx <
numPhases);
391 return fluidState.averageMolarMass(phaseIdx)/paramCache.
molarVolume(phaseIdx);
403 template <
class Flu
idState>
415 template <
class Flu
idState>
420 assert(0 <= phaseIdx && phaseIdx <=
numPhases);
425 return 0.0170e-2 * 0.1;
434 return 0.208e-2 * 0.1;
453 template <
class Flu
idState>
459 assert(0 <= phaseIdx && phaseIdx <=
numPhases);
469 return henryCoeffWater_(compIdx, fs.temperature(
wPhaseIdx))
498 template <
class Flu
idState>
503 { DUNE_THROW(Dune::NotImplemented,
"Diffusion coefficients"); }
516 template <
class Flu
idState>
522 { DUNE_THROW(Dune::NotImplemented,
"Binary diffusion coefficients"); }
531 template <
class Flu
idState>
535 { DUNE_THROW(Dune::NotImplemented,
"Enthalpies"); }
544 template <
class Flu
idState>
548 { DUNE_THROW(Dune::NotImplemented,
"Thermal conductivities"); }
557 template <
class Flu
idState>
561 { DUNE_THROW(Dune::NotImplemented,
"Heat capacities"); }
565 static Scalar henryCoeffWater_(
int compIdx, Scalar
temperature)
574 case C1Idx:
return 5.57601e+09;
575 case C3Idx:
return 1e10;
576 case C6Idx:
return 1e10;
580 default: DUNE_THROW(Dune::InvalidStateException,
"Unknown component index " << compIdx);
Provides 3rd order polynomial splines.
A collection of input/output field names for common physical quantities.
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...
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
std::string temperature() noexcept
I/O name of temperature for equilibrium models.
Definition: name.hh:51
std::string gaseousPhase() noexcept
I/O name of gaseous phase.
Definition: name.hh:123
std::string naplPhase() noexcept
I/O name of napl phase.
Definition: name.hh:131
std::string aqueousPhase() noexcept
I/O name of aqueous phase.
Definition: name.hh:127
Material properties of pure water .
Definition: h2o.hh:61
static Scalar vaporPressure(Scalar T)
The vapor pressure in of pure water at a given temperature.
Definition: h2o.hh:131
static constexpr Scalar criticalTemperature()
Returns the critical temperature of water.
Definition: h2o.hh:93
static std::string name()
A human readable name for the water.
Definition: h2o.hh:75
static constexpr Scalar molarMass()
The molar mass in of water.
Definition: h2o.hh:81
static constexpr Scalar criticalPressure()
Returns the critical pressure of water.
Definition: h2o.hh:99
static constexpr Scalar criticalMolarVolume()
Returns the molar volume of water at the critical point.
Definition: h2o.hh:105
Implements the Peng-Robinson equation of state for liquids and gases.
Definition: pengrobinson.hh:57
static void init(Scalar aMin, Scalar aMax, int na, Scalar bMin, Scalar bMax, int nb)
Definition: pengrobinson.hh:62
Implements the Peng-Robinson equation of state for a mixture.
Definition: pengrobinsonmixture.hh:41
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:90
Scalar a() const
Returns the attractive parameter 'a' of the Peng-Robinson fluid.
Definition: pengrobinsonparams.hh:52
Scalar b() const
Returns the repulsive parameter 'b' of the Peng-Robinson fluid.
Definition: pengrobinsonparams.hh:60
The mixing rule for the oil and the gas phases of the SPE5 problem.
Definition: pengrobinsonparamsmixture.hh:63
const PureParams & pureParams(int compIdx) const
Return the Peng-Robinson parameters of a pure substance,.
Definition: pengrobinsonparamsmixture.hh:204
void updatePure(const FluidState &fluidState)
Update Peng-Robinson parameters for the pure components.
Definition: pengrobinsonparamsmixture.hh:76
The fluid system for the SPE-5 benchmark problem.
Definition: spe5.hh:57
static Scalar molarDensity(const FluidState &fluidState, const ParameterCache ¶mCache, int phaseIdx)
The molar density of a fluid phase in .
Definition: spe5.hh:404
static bool isIdealMixture(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: spe5.hh:140
static Scalar criticalTemperature(int compIdx)
Critical temperature of a component .
Definition: spe5.hh:219
static const int C15Idx
Definition: spe5.hh:172
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:517
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:558
static const int C3Idx
Definition: spe5.hh:169
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:454
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:545
static const int H2OIdx
Definition: spe5.hh:167
static std::string componentName(int compIdx)
Return the human readable name of a component.
Definition: spe5.hh:179
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:499
static const int C6Idx
Definition: spe5.hh:170
static const int C20Idx
Definition: spe5.hh:173
static Scalar molarMass(int compIdx)
Return the molar mass of a component in .
Definition: spe5.hh:199
static constexpr bool isMiscible()
Returns whether the fluids are miscible.
Definition: spe5.hh:102
static bool isIdealGas(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: spe5.hh:154
static const int C1Idx
Definition: spe5.hh:168
static const int numComponents
Number of components in the fluid system.
Definition: spe5.hh:165
static const int oPhaseIdx
Index of the oil phase.
Definition: spe5.hh:78
static Scalar density(const FluidState &fluidState, const ParameterCache ¶mCache, int phaseIdx)
Calculate the density of a fluid phase.
Definition: spe5.hh:385
static const int C10Idx
Definition: spe5.hh:171
static Scalar acentricFactor(int compIdx)
The acentric factor of a component .
Definition: spe5.hh:280
static const int gPhaseIdx
Index of the gas phase.
Definition: spe5.hh:74
static Scalar viscosity(const FluidState &fs, const ParameterCache ¶mCache, int phaseIdx)
Calculate the dynamic viscosity of a fluid phase .
Definition: spe5.hh:416
static void init()
Initialize the fluid system's static parameters.
Definition: spe5.hh:323
static Scalar interactionCoefficient(int comp1Idx, int comp2Idx)
Returns the interaction coefficient for two components.
Definition: spe5.hh:303
static const int numPhases
Number of phases in the fluid system.
Definition: spe5.hh:71
static std::string phaseName(int phaseIdx)
Return the human readable name of a fluid phase.
Definition: spe5.hh:87
static const int wPhaseIdx
Index of the water phase.
Definition: spe5.hh:76
static constexpr bool isGas(int phaseIdx)
Return whether a phase is gaseous.
Definition: spe5.hh:109
static Scalar criticalPressure(int compIdx)
Critical pressure of a component .
Definition: spe5.hh:239
static Scalar criticalMolarVolume(int compIdx)
Molar volume of a component at the critical point .
Definition: spe5.hh:259
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:532
static constexpr bool isCompressible(int phaseIdx)
Return whether a phase is compressible.
Definition: spe5.hh:121
Specifies the parameters required by the SPE5 problem which are despondent on the thermodynamic state...
Definition: spe5parametercache.hh:45
Scalar molarVolume(int phaseIdx) const
Returns the molar volume of a phase .
Definition: spe5parametercache.hh:200