Cantera  3.0.0
Loading...
Searching...
No Matches
PengRobinson Class Reference

Implementation of a multi-species Peng-Robinson equation of state. More...

#include <PengRobinson.h>

Inheritance diagram for PengRobinson:
[legend]

Detailed Description

Implementation of a multi-species Peng-Robinson equation of state.

Definition at line 19 of file PengRobinson.h.

Public Member Functions

 PengRobinson (const string &infile="", const string &id="")
 Construct and initialize a PengRobinson object directly from an input file.
 
string type () const override
 String indicating the thermodynamic model implemented.
 
double standardConcentration (size_t k=0) const override
 Returns the standard concentration \( C^0_k \), which is used to normalize the generalized concentration.
 
void getActivityCoefficients (double *ac) const override
 Get the array of non-dimensional activity coefficients at the current solution temperature, pressure, and solution concentration.
 
double speciesCritTemperature (double a, double b) const
 Calculate species-specific critical temperature.
 
double isothermalCompressibility () const override
 Returns the isothermal compressibility. Units: 1/Pa.
 
double thermalExpansionCoeff () const override
 Return the volumetric thermal expansion coefficient. Units: 1/K.
 
double soundSpeed () const override
 Return the speed of sound. Units: m/s.
 
void calculatePressureDerivatives () const
 Calculate \( dp/dV \) and \( dp/dT \) at the current conditions.
 
void updateMixingExpressions () override
 Update the \( a \), \( b \), and \( \alpha \) parameters.
 
void calculateAB (double &aCalc, double &bCalc, double &aAlpha) const
 Calculate the \( a \), \( b \), and \( \alpha \) parameters given the temperature.
 
void calcCriticalConditions (double &pc, double &tc, double &vc) const override
 
int solveCubic (double T, double pres, double a, double b, double aAlpha, double Vroot[3]) const
 Prepare variables and call the function to solve the cubic equation of state.
 
Molar Thermodynamic properties
double cp_mole () const override
 Molar heat capacity at constant pressure. Units: J/kmol/K.
 
double cv_mole () const override
 Molar heat capacity at constant volume. Units: J/kmol/K.
 
Mechanical Properties
double pressure () const override
 Return the thermodynamic pressure (Pa).
 
Partial Molar Properties of the Solution
void getChemPotentials (double *mu) const override
 Get the species chemical potentials. Units: J/kmol.
 
void getPartialMolarEnthalpies (double *hbar) const override
 Returns an array of partial molar enthalpies for the species in the mixture.
 
void getPartialMolarEntropies (double *sbar) const override
 Returns an array of partial molar entropies of the species in the solution.
 
void getPartialMolarIntEnergies (double *ubar) const override
 Return an array of partial molar internal energies for the species in the mixture.
 
void getPartialMolarCp (double *cpbar) const override
 Calculate species-specific molar specific heats.
 
void getPartialMolarVolumes (double *vbar) const override
 Return an array of partial molar volumes for the species in the mixture.
 
Initialization Methods - For Internal use

The following methods are used in the process of constructing the phase and setting its parameters from a specification in an input file.

They are not normally used in application programs. To see how they are used, see importPhase().

bool addSpecies (shared_ptr< Species > spec) override
 Add a Species to this Phase.
 
void initThermo () override
 Initialize the ThermoPhase object after all species have been set up.
 
void getSpeciesParameters (const string &name, AnyMap &speciesNode) const override
 Get phase-specific parameters of a Species object such that an identical one could be reconstructed and added to this phase.
 
void setSpeciesCoeffs (const string &species, double a, double b, double w)
 Set the pure fluid interaction parameters for a species.
 
void setBinaryCoeffs (const string &species_i, const string &species_j, double a)
 Set values for the interaction parameter between two species.
 
- Public Member Functions inherited from MixtureFugacityTP
void setTemperature (const double temp) override
 Set the temperature of the phase.
 
void setPressure (double p) override
 Set the internally stored pressure (Pa) at constant temperature and composition.
 
 MixtureFugacityTP ()=default
 Constructor.
 
int standardStateConvention () const override
 This method returns the convention used in specification of the standard state, of which there are currently two, temperature based, and variable pressure based.
 
void setForcedSolutionBranch (int solnBranch)
 Set the solution branch to force the ThermoPhase to exist on one branch or another.
 
int forcedSolutionBranch () const
 Report the solution branch which the solution is restricted to.
 
int reportSolnBranchActual () const
 Report the solution branch which the solution is actually on.
 
double enthalpy_mole () const override
 Molar enthalpy. Units: J/kmol.
 
double entropy_mole () const override
 Molar entropy. Units: J/kmol/K.
 
void getChemPotentials_RT (double *mu) const override
 Get the array of non-dimensional species chemical potentials These are partial molar Gibbs free energies.
 
void getStandardChemPotentials (double *mu) const override
 Get the array of chemical potentials at unit activity.
 
void getEnthalpy_RT (double *hrt) const override
 Get the nondimensional Enthalpy functions for the species at their standard states at the current T and P of the solution.
 
void getEntropy_R (double *sr) const override
 Get the array of nondimensional Enthalpy functions for the standard state species at the current T and P of the solution.
 
void getGibbs_RT (double *grt) const override
 Get the nondimensional Gibbs functions for the species at their standard states of solution at the current T and P of the solution.
 
void getPureGibbs (double *gpure) const override
 Get the pure Gibbs free energies of each species.
 
void getIntEnergy_RT (double *urt) const override
 Returns the vector of nondimensional internal Energies of the standard state at the current temperature and pressure of the solution for each species.
 
void getCp_R (double *cpr) const override
 Get the nondimensional Heat Capacities at constant pressure for the standard state of the species at the current T and P.
 
void getStandardVolumes (double *vol) const override
 Get the molar volumes of each species in their standard states at the current T and P of the solution.
 
void getEnthalpy_RT_ref (double *hrt) const override
 Returns the vector of nondimensional enthalpies of the reference state at the current temperature of the solution and the reference pressure for the species.
 
void getGibbs_RT_ref (double *grt) const override
 Returns the vector of nondimensional Gibbs Free Energies of the reference state at the current temperature of the solution and the reference pressure for the species.
 
void getGibbs_ref (double *g) const override
 Returns the vector of the Gibbs function of the reference state at the current temperature of the solution and the reference pressure for the species.
 
void getEntropy_R_ref (double *er) const override
 Returns the vector of nondimensional entropies of the reference state at the current temperature of the solution and the reference pressure for each species.
 
void getCp_R_ref (double *cprt) const override
 Returns the vector of nondimensional constant pressure heat capacities of the reference state at the current temperature of the solution and reference pressure for each species.
 
void getStandardVolumes_ref (double *vol) const override
 Get the molar volumes of the species reference states at the current T and P_ref of the solution.
 
int phaseState (bool checkState=false) const
 Returns the Phase State flag for the current state of the object.
 
double calculatePsat (double TKelvin, double &molarVolGas, double &molarVolLiquid)
 Calculate the saturation pressure at the current mixture content for the given temperature.
 
double satPressure (double TKelvin) override
 Calculate the saturation pressure at the current mixture content for the given temperature.
 
void getActivityConcentrations (double *c) const override
 This method returns an array of generalized concentrations.
 
- Public Member Functions inherited from ThermoPhase
 ThermoPhase ()=default
 Constructor.
 
double RT () const
 Return the Gas Constant multiplied by the current temperature.
 
double equivalenceRatio () const
 Compute the equivalence ratio for the current mixture from available oxygen and required oxygen.
 
string type () const override
 String indicating the thermodynamic model implemented.
 
virtual bool isIdeal () const
 Boolean indicating whether phase is ideal.
 
virtual string phaseOfMatter () const
 String indicating the mechanical phase of the matter in this Phase.
 
virtual double refPressure () const
 Returns the reference pressure in Pa.
 
virtual double minTemp (size_t k=npos) const
 Minimum temperature for which the thermodynamic data for the species or phase are valid.
 
double Hf298SS (const size_t k) const
 Report the 298 K Heat of Formation of the standard state of one species (J kmol-1)
 
virtual void modifyOneHf298SS (const size_t k, const double Hf298New)
 Modify the value of the 298 K Heat of Formation of one species in the phase (J kmol-1)
 
virtual void resetHf298 (const size_t k=npos)
 Restore the original heat of formation of one or more species.
 
virtual double maxTemp (size_t k=npos) const
 Maximum temperature for which the thermodynamic data for the species are valid.
 
bool chargeNeutralityNecessary () const
 Returns the chargeNeutralityNecessity boolean.
 
virtual double intEnergy_mole () const
 Molar internal energy. Units: J/kmol.
 
virtual double gibbs_mole () const
 Molar Gibbs function. Units: J/kmol.
 
void setElectricPotential (double v)
 Set the electric potential of this phase (V).
 
double electricPotential () const
 Returns the electric potential of this phase (V).
 
virtual int activityConvention () const
 This method returns the convention used in specification of the activities, of which there are currently two, molar- and molality-based conventions.
 
virtual Units standardConcentrationUnits () const
 Returns the units of the "standard concentration" for this phase.
 
virtual double logStandardConc (size_t k=0) const
 Natural logarithm of the standard concentration of the kth species.
 
virtual void getActivities (double *a) const
 Get the array of non-dimensional activities at the current solution temperature, pressure, and solution concentration.
 
virtual void getLnActivityCoefficients (double *lnac) const
 Get the array of non-dimensional molar-based ln activity coefficients at the current solution temperature, pressure, and solution concentration.
 
void getElectrochemPotentials (double *mu) const
 Get the species electrochemical potentials.
 
virtual void getIntEnergy_RT_ref (double *urt) const
 Returns the vector of nondimensional internal Energies of the reference state at the current temperature of the solution and the reference pressure for each species.
 
double enthalpy_mass () const
 Specific enthalpy. Units: J/kg.
 
double intEnergy_mass () const
 Specific internal energy. Units: J/kg.
 
double entropy_mass () const
 Specific entropy. Units: J/kg/K.
 
double gibbs_mass () const
 Specific Gibbs function. Units: J/kg.
 
double cp_mass () const
 Specific heat at constant pressure. Units: J/kg/K.
 
double cv_mass () const
 Specific heat at constant volume. Units: J/kg/K.
 
virtual void setState_TPX (double t, double p, const double *x)
 Set the temperature (K), pressure (Pa), and mole fractions.
 
virtual void setState_TPX (double t, double p, const Composition &x)
 Set the temperature (K), pressure (Pa), and mole fractions.
 
virtual void setState_TPX (double t, double p, const string &x)
 Set the temperature (K), pressure (Pa), and mole fractions.
 
virtual void setState_TPY (double t, double p, const double *y)
 Set the internally stored temperature (K), pressure (Pa), and mass fractions of the phase.
 
virtual void setState_TPY (double t, double p, const Composition &y)
 Set the internally stored temperature (K), pressure (Pa), and mass fractions of the phase.
 
virtual void setState_TPY (double t, double p, const string &y)
 Set the internally stored temperature (K), pressure (Pa), and mass fractions of the phase.
 
virtual void setState_TP (double t, double p)
 Set the temperature (K) and pressure (Pa)
 
virtual void setState_PX (double p, double *x)
 Set the pressure (Pa) and mole fractions.
 
virtual void setState_PY (double p, double *y)
 Set the internally stored pressure (Pa) and mass fractions.
 
virtual void setState_HP (double h, double p, double tol=1e-9)
 Set the internally stored specific enthalpy (J/kg) and pressure (Pa) of the phase.
 
virtual void setState_UV (double u, double v, double tol=1e-9)
 Set the specific internal energy (J/kg) and specific volume (m^3/kg).
 
virtual void setState_SP (double s, double p, double tol=1e-9)
 Set the specific entropy (J/kg/K) and pressure (Pa).
 
virtual void setState_SV (double s, double v, double tol=1e-9)
 Set the specific entropy (J/kg/K) and specific volume (m^3/kg).
 
virtual void setState_ST (double s, double t, double tol=1e-9)
 Set the specific entropy (J/kg/K) and temperature (K).
 
virtual void setState_TV (double t, double v, double tol=1e-9)
 Set the temperature (K) and specific volume (m^3/kg).
 
virtual void setState_PV (double p, double v, double tol=1e-9)
 Set the pressure (Pa) and specific volume (m^3/kg).
 
virtual void setState_UP (double u, double p, double tol=1e-9)
 Set the specific internal energy (J/kg) and pressure (Pa).
 
virtual void setState_VH (double v, double h, double tol=1e-9)
 Set the specific volume (m^3/kg) and the specific enthalpy (J/kg)
 
virtual void setState_TH (double t, double h, double tol=1e-9)
 Set the temperature (K) and the specific enthalpy (J/kg)
 
virtual void setState_SH (double s, double h, double tol=1e-9)
 Set the specific entropy (J/kg/K) and the specific enthalpy (J/kg)
 
void setState_RP (double rho, double p)
 Set the density (kg/m**3) and pressure (Pa) at constant composition.
 
virtual void setState_DP (double rho, double p)
 Set the density (kg/m**3) and pressure (Pa) at constant composition.
 
virtual void setState_RPX (double rho, double p, const double *x)
 Set the density (kg/m**3), pressure (Pa) and mole fractions.
 
virtual void setState_RPX (double rho, double p, const Composition &x)
 Set the density (kg/m**3), pressure (Pa) and mole fractions.
 
virtual void setState_RPX (double rho, double p, const string &x)
 Set the density (kg/m**3), pressure (Pa) and mole fractions.
 
virtual void setState_RPY (double rho, double p, const double *y)
 Set the density (kg/m**3), pressure (Pa) and mass fractions.
 
virtual void setState_RPY (double rho, double p, const Composition &y)
 Set the density (kg/m**3), pressure (Pa) and mass fractions.
 
virtual void setState_RPY (double rho, double p, const string &y)
 Set the density (kg/m**3), pressure (Pa) and mass fractions.
 
virtual void setState (const AnyMap &state)
 Set the state using an AnyMap containing any combination of properties supported by the thermodynamic model.
 
void setMixtureFraction (double mixFrac, const double *fuelComp, const double *oxComp, ThermoBasis basis=ThermoBasis::molar)
 Set the mixture composition according to the mixture fraction = kg fuel / (kg oxidizer + kg fuel)
 
void setMixtureFraction (double mixFrac, const string &fuelComp, const string &oxComp, ThermoBasis basis=ThermoBasis::molar)
 Set the mixture composition according to the mixture fraction = kg fuel / (kg oxidizer + kg fuel)
 
void setMixtureFraction (double mixFrac, const Composition &fuelComp, const Composition &oxComp, ThermoBasis basis=ThermoBasis::molar)
 Set the mixture composition according to the mixture fraction = kg fuel / (kg oxidizer + kg fuel)
 
double mixtureFraction (const double *fuelComp, const double *oxComp, ThermoBasis basis=ThermoBasis::molar, const string &element="Bilger") const
 Compute the mixture fraction = kg fuel / (kg oxidizer + kg fuel) for the current mixture given fuel and oxidizer compositions.
 
double mixtureFraction (const string &fuelComp, const string &oxComp, ThermoBasis basis=ThermoBasis::molar, const string &element="Bilger") const
 Compute the mixture fraction = kg fuel / (kg oxidizer + kg fuel) for the current mixture given fuel and oxidizer compositions.
 
double mixtureFraction (const Composition &fuelComp, const Composition &oxComp, ThermoBasis basis=ThermoBasis::molar, const string &element="Bilger") const
 Compute the mixture fraction = kg fuel / (kg oxidizer + kg fuel) for the current mixture given fuel and oxidizer compositions.
 
void setEquivalenceRatio (double phi, const double *fuelComp, const double *oxComp, ThermoBasis basis=ThermoBasis::molar)
 Set the mixture composition according to the equivalence ratio.
 
void setEquivalenceRatio (double phi, const string &fuelComp, const string &oxComp, ThermoBasis basis=ThermoBasis::molar)
 Set the mixture composition according to the equivalence ratio.
 
void setEquivalenceRatio (double phi, const Composition &fuelComp, const Composition &oxComp, ThermoBasis basis=ThermoBasis::molar)
 Set the mixture composition according to the equivalence ratio.
 
double equivalenceRatio (const double *fuelComp, const double *oxComp, ThermoBasis basis=ThermoBasis::molar) const
 Compute the equivalence ratio for the current mixture given the compositions of fuel and oxidizer.
 
double equivalenceRatio (const string &fuelComp, const string &oxComp, ThermoBasis basis=ThermoBasis::molar) const
 Compute the equivalence ratio for the current mixture given the compositions of fuel and oxidizer.
 
double equivalenceRatio (const Composition &fuelComp, const Composition &oxComp, ThermoBasis basis=ThermoBasis::molar) const
 Compute the equivalence ratio for the current mixture given the compositions of fuel and oxidizer.
 
double stoichAirFuelRatio (const double *fuelComp, const double *oxComp, ThermoBasis basis=ThermoBasis::molar) const
 Compute the stoichiometric air to fuel ratio (kg oxidizer / kg fuel) given fuel and oxidizer compositions.
 
double stoichAirFuelRatio (const string &fuelComp, const string &oxComp, ThermoBasis basis=ThermoBasis::molar) const
 Compute the stoichiometric air to fuel ratio (kg oxidizer / kg fuel) given fuel and oxidizer compositions.
 
double stoichAirFuelRatio (const Composition &fuelComp, const Composition &oxComp, ThermoBasis basis=ThermoBasis::molar) const
 Compute the stoichiometric air to fuel ratio (kg oxidizer / kg fuel) given fuel and oxidizer compositions.
 
void equilibrate (const string &XY, const string &solver="auto", double rtol=1e-9, int max_steps=50000, int max_iter=100, int estimate_equil=0, int log_level=0)
 Equilibrate a ThermoPhase object.
 
virtual void setToEquilState (const double *mu_RT)
 This method is used by the ChemEquil equilibrium solver.
 
virtual bool compatibleWithMultiPhase () const
 Indicates whether this phase type can be used with class MultiPhase for equilibrium calculations.
 
virtual double satTemperature (double p) const
 Return the saturation temperature given the pressure.
 
virtual double vaporFraction () const
 Return the fraction of vapor at the current conditions.
 
virtual void setState_Tsat (double t, double x)
 Set the state to a saturated system at a particular temperature.
 
virtual void setState_Psat (double p, double x)
 Set the state to a saturated system at a particular pressure.
 
void setState_TPQ (double T, double P, double Q)
 Set the temperature, pressure, and vapor fraction (quality).
 
bool addSpecies (shared_ptr< Species > spec) override
 Add a Species to this Phase.
 
void modifySpecies (size_t k, shared_ptr< Species > spec) override
 Modify the thermodynamic data associated with a species.
 
virtual MultiSpeciesThermospeciesThermo (int k=-1)
 Return a changeable reference to the calculation manager for species reference-state thermodynamic properties.
 
virtual const MultiSpeciesThermospeciesThermo (int k=-1) const
 
void initThermoFile (const string &inputFile, const string &id)
 Initialize a ThermoPhase object using an input file.
 
virtual void setParameters (const AnyMap &phaseNode, const AnyMap &rootNode=AnyMap())
 Set equation of state parameters from an AnyMap phase description.
 
AnyMap parameters (bool withInput=true) const
 Returns the parameters of a ThermoPhase object such that an identical one could be reconstructed using the newThermo(AnyMap&) function.
 
const AnyMapinput () const
 Access input data associated with the phase description.
 
AnyMapinput ()
 
void invalidateCache () override
 Invalidate any cached values which are normally updated only when a change in state is detected.
 
virtual void getdlnActCoeffds (const double dTds, const double *const dXds, double *dlnActCoeffds) const
 Get the change in activity coefficients wrt changes in state (temp, mole fraction, etc) along a line in parameter space or along a line in physical space.
 
virtual void getdlnActCoeffdlnX_diag (double *dlnActCoeffdlnX_diag) const
 Get the array of ln mole fraction derivatives of the log activity coefficients - diagonal component only.
 
virtual void getdlnActCoeffdlnN_diag (double *dlnActCoeffdlnN_diag) const
 Get the array of log species mole number derivatives of the log activity coefficients.
 
virtual void getdlnActCoeffdlnN (const size_t ld, double *const dlnActCoeffdlnN)
 Get the array of derivatives of the log activity coefficients with respect to the log of the species mole numbers.
 
virtual void getdlnActCoeffdlnN_numderiv (const size_t ld, double *const dlnActCoeffdlnN)
 
virtual string report (bool show_thermo=true, double threshold=-1e-14) const
 returns a summary of the state of the phase as a string
 
virtual void reportCSV (std::ofstream &csvFile) const
 returns a summary of the state of the phase to a comma separated file.
 
- Public Member Functions inherited from Phase
 Phase ()=default
 Default constructor.
 
 Phase (const Phase &)=delete
 
Phaseoperator= (const Phase &)=delete
 
virtual bool isPure () const
 Return whether phase represents a pure (single species) substance.
 
virtual bool hasPhaseTransition () const
 Return whether phase represents a substance with phase transitions.
 
virtual bool isCompressible () const
 Return whether phase represents a compressible substance.
 
virtual map< string, size_t > nativeState () const
 Return a map of properties defining the native state of a substance.
 
string nativeMode () const
 Return string acronym representing the native state of a Phase.
 
virtual vector< string > fullStates () const
 Return a vector containing full states defining a phase.
 
virtual vector< string > partialStates () const
 Return a vector of settable partial property sets within a phase.
 
virtual size_t stateSize () const
 Return size of vector defining internal state of the phase.
 
void saveState (vector< double > &state) const
 Save the current internal state of the phase.
 
virtual void saveState (size_t lenstate, double *state) const
 Write to array 'state' the current internal state.
 
void restoreState (const vector< double > &state)
 Restore a state saved on a previous call to saveState.
 
virtual void restoreState (size_t lenstate, const double *state)
 Restore the state of the phase from a previously saved state vector.
 
double molecularWeight (size_t k) const
 Molecular weight of species k.
 
void getMolecularWeights (vector< double > &weights) const
 Copy the vector of molecular weights into vector weights.
 
void getMolecularWeights (double *weights) const
 Copy the vector of molecular weights into array weights.
 
const vector< double > & molecularWeights () const
 Return a const reference to the internal vector of molecular weights.
 
const vector< double > & inverseMolecularWeights () const
 Return a const reference to the internal vector of molecular weights.
 
void getCharges (double *charges) const
 Copy the vector of species charges into array charges.
 
virtual void setMolesNoTruncate (const double *const N)
 Set the state of the object with moles in [kmol].
 
double elementalMassFraction (const size_t m) const
 Elemental mass fraction of element m.
 
double elementalMoleFraction (const size_t m) const
 Elemental mole fraction of element m.
 
const double * moleFractdivMMW () const
 Returns a const pointer to the start of the moleFraction/MW array.
 
double charge (size_t k) const
 Dimensionless electrical charge of a single molecule of species k The charge is normalized by the the magnitude of the electron charge.
 
double chargeDensity () const
 Charge density [C/m^3].
 
size_t nDim () const
 Returns the number of spatial dimensions (1, 2, or 3)
 
void setNDim (size_t ndim)
 Set the number of spatial dimensions (1, 2, or 3).
 
virtual bool ready () const
 Returns a bool indicating whether the object is ready for use.
 
int stateMFNumber () const
 Return the State Mole Fraction Number.
 
virtual void invalidateCache ()
 Invalidate any cached values which are normally updated only when a change in state is detected.
 
bool caseSensitiveSpecies () const
 Returns true if case sensitive species names are enforced.
 
void setCaseSensitiveSpecies (bool cflag=true)
 Set flag that determines whether case sensitive species are enforced in look-up operations, for example speciesIndex.
 
vector< double > getCompositionFromMap (const Composition &comp) const
 Converts a Composition to a vector with entries for each species Species that are not specified are set to zero in the vector.
 
void massFractionsToMoleFractions (const double *Y, double *X) const
 Converts a mixture composition from mole fractions to mass fractions.
 
void moleFractionsToMassFractions (const double *X, double *Y) const
 Converts a mixture composition from mass fractions to mole fractions.
 
string name () const
 Return the name of the phase.
 
void setName (const string &nm)
 Sets the string name for the phase.
 
string elementName (size_t m) const
 Name of the element with index m.
 
size_t elementIndex (const string &name) const
 Return the index of element named 'name'.
 
const vector< string > & elementNames () const
 Return a read-only reference to the vector of element names.
 
double atomicWeight (size_t m) const
 Atomic weight of element m.
 
double entropyElement298 (size_t m) const
 Entropy of the element in its standard state at 298 K and 1 bar.
 
int atomicNumber (size_t m) const
 Atomic number of element m.
 
int elementType (size_t m) const
 Return the element constraint type Possible types include:
 
int changeElementType (int m, int elem_type)
 Change the element type of the mth constraint Reassigns an element type.
 
const vector< double > & atomicWeights () const
 Return a read-only reference to the vector of atomic weights.
 
size_t nElements () const
 Number of elements.
 
void checkElementIndex (size_t m) const
 Check that the specified element index is in range.
 
void checkElementArraySize (size_t mm) const
 Check that an array size is at least nElements().
 
double nAtoms (size_t k, size_t m) const
 Number of atoms of element m in species k.
 
void getAtoms (size_t k, double *atomArray) const
 Get a vector containing the atomic composition of species k.
 
size_t speciesIndex (const string &name) const
 Returns the index of a species named 'name' within the Phase object.
 
string speciesName (size_t k) const
 Name of the species with index k.
 
string speciesSPName (int k) const
 Returns the expanded species name of a species, including the phase name This is guaranteed to be unique within a Cantera problem.
 
const vector< string > & speciesNames () const
 Return a const reference to the vector of species names.
 
size_t nSpecies () const
 Returns the number of species in the phase.
 
void checkSpeciesIndex (size_t k) const
 Check that the specified species index is in range.
 
void checkSpeciesArraySize (size_t kk) const
 Check that an array size is at least nSpecies().
 
void setMoleFractionsByName (const Composition &xMap)
 Set the species mole fractions by name.
 
void setMoleFractionsByName (const string &x)
 Set the mole fractions of a group of species by name.
 
void setMassFractionsByName (const Composition &yMap)
 Set the species mass fractions by name.
 
void setMassFractionsByName (const string &x)
 Set the species mass fractions by name.
 
void setState_TRX (double t, double dens, const double *x)
 Set the internally stored temperature (K), density, and mole fractions.
 
void setState_TRX (double t, double dens, const Composition &x)
 Set the internally stored temperature (K), density, and mole fractions.
 
void setState_TRY (double t, double dens, const double *y)
 Set the internally stored temperature (K), density, and mass fractions.
 
void setState_TRY (double t, double dens, const Composition &y)
 Set the internally stored temperature (K), density, and mass fractions.
 
void setState_TNX (double t, double n, const double *x)
 Set the internally stored temperature (K), molar density (kmol/m^3), and mole fractions.
 
void setState_TR (double t, double rho)
 Set the internally stored temperature (K) and density (kg/m^3)
 
void setState_TD (double t, double rho)
 Set the internally stored temperature (K) and density (kg/m^3)
 
void setState_TX (double t, double *x)
 Set the internally stored temperature (K) and mole fractions.
 
void setState_TY (double t, double *y)
 Set the internally stored temperature (K) and mass fractions.
 
void setState_RX (double rho, double *x)
 Set the density (kg/m^3) and mole fractions.
 
void setState_RY (double rho, double *y)
 Set the density (kg/m^3) and mass fractions.
 
Composition getMoleFractionsByName (double threshold=0.0) const
 Get the mole fractions by name.
 
double moleFraction (size_t k) const
 Return the mole fraction of a single species.
 
double moleFraction (const string &name) const
 Return the mole fraction of a single species.
 
Composition getMassFractionsByName (double threshold=0.0) const
 Get the mass fractions by name.
 
double massFraction (size_t k) const
 Return the mass fraction of a single species.
 
double massFraction (const string &name) const
 Return the mass fraction of a single species.
 
void getMoleFractions (double *const x) const
 Get the species mole fraction vector.
 
virtual void setMoleFractions (const double *const x)
 Set the mole fractions to the specified values.
 
virtual void setMoleFractions_NoNorm (const double *const x)
 Set the mole fractions to the specified values without normalizing.
 
void getMassFractions (double *const y) const
 Get the species mass fractions.
 
const double * massFractions () const
 Return a const pointer to the mass fraction array.
 
virtual void setMassFractions (const double *const y)
 Set the mass fractions to the specified values and normalize them.
 
virtual void setMassFractions_NoNorm (const double *const y)
 Set the mass fractions to the specified values without normalizing.
 
virtual void getConcentrations (double *const c) const
 Get the species concentrations (kmol/m^3).
 
virtual double concentration (const size_t k) const
 Concentration of species k.
 
virtual void setConcentrations (const double *const conc)
 Set the concentrations to the specified values within the phase.
 
virtual void setConcentrationsNoNorm (const double *const conc)
 Set the concentrations without ignoring negative concentrations.
 
double temperature () const
 Temperature (K).
 
virtual double electronTemperature () const
 Electron Temperature (K)
 
virtual double density () const
 Density (kg/m^3).
 
virtual double molarDensity () const
 Molar density (kmol/m^3).
 
virtual double molarVolume () const
 Molar volume (m^3/kmol).
 
virtual void setDensity (const double density_)
 Set the internally stored density (kg/m^3) of the phase.
 
virtual void setMolarDensity (const double molarDensity)
 Set the internally stored molar density (kmol/m^3) of the phase.
 
virtual void setElectronTemperature (double etemp)
 Set the internally stored electron temperature of the phase (K).
 
double mean_X (const double *const Q) const
 Evaluate the mole-fraction-weighted mean of an array Q.
 
double mean_X (const vector< double > &Q) const
 Evaluate the mole-fraction-weighted mean of an array Q.
 
double meanMolecularWeight () const
 The mean molecular weight. Units: (kg/kmol)
 
double sum_xlogx () const
 Evaluate \( \sum_k X_k \ln X_k \).
 
size_t addElement (const string &symbol, double weight=-12345.0, int atomicNumber=0, double entropy298=ENTROPY298_UNKNOWN, int elem_type=CT_ELEM_TYPE_ABSPOS)
 Add an element.
 
void addSpeciesAlias (const string &name, const string &alias)
 Add a species alias (that is, a user-defined alternative species name).
 
virtual vector< string > findIsomers (const Composition &compMap) const
 Return a vector with isomers names matching a given composition map.
 
virtual vector< string > findIsomers (const string &comp) const
 Return a vector with isomers names matching a given composition string.
 
shared_ptr< Speciesspecies (const string &name) const
 Return the Species object for the named species.
 
shared_ptr< Speciesspecies (size_t k) const
 Return the Species object for species whose index is k.
 
void ignoreUndefinedElements ()
 Set behavior when adding a species containing undefined elements to just skip the species.
 
void addUndefinedElements ()
 Set behavior when adding a species containing undefined elements to add those elements to the phase.
 
void throwUndefinedElements ()
 Set the behavior when adding a species containing undefined elements to throw an exception.
 

Protected Types

enum class  CoeffSource { EoS , CritProps , Database }
 

Protected Member Functions

double sresid () const override
 Calculate the deviation terms for the total entropy of the mixture from the ideal gas mixture.
 
double hresid () const override
 Calculate the deviation terms for the total enthalpy of the mixture from the ideal gas mixture.
 
double liquidVolEst (double T, double &pres) const override
 Estimate for the molar volume of the liquid.
 
double densityCalc (double T, double pressure, int phase, double rhoguess) override
 Calculates the density given the temperature and the pressure and a guess at the density.
 
double densSpinodalLiquid () const override
 Return the value of the density at the liquid spinodal point (on the liquid side) for the current temperature.
 
double densSpinodalGas () const override
 Return the value of the density at the gas spinodal point (on the gas side) for the current temperature.
 
double dpdVCalc (double T, double molarVol, double &presCalc) const override
 Calculate the pressure and the pressure derivative given the temperature and the molar volume.
 
double daAlpha_dT () const
 Calculate temperature derivative \( d(a \alpha)/dT \).
 
double d2aAlpha_dT2 () const
 Calculate second derivative \( d^2(a \alpha)/dT^2 \).
 
- Protected Member Functions inherited from MixtureFugacityTP
void compositionChanged () override
 Apply changes to the state which are needed after the composition changes.
 
virtual void _updateReferenceStateThermo () const
 Updates the reference state thermodynamic functions at the current T of the solution.
 
double critTemperature () const override
 Critical temperature (K).
 
double critPressure () const override
 Critical pressure (Pa).
 
double critVolume () const override
 Critical volume (m3/kmol).
 
double critCompressibility () const override
 Critical compressibility (unitless).
 
double critDensity () const override
 Critical density (kg/m3).
 
int solveCubic (double T, double pres, double a, double b, double aAlpha, double Vroot[3], double an, double bn, double cn, double dn, double tc, double vc) const
 Solve the cubic equation of state.
 
const vector< double > & gibbs_RT_ref () const
 Returns the vector of nondimensional Gibbs free energies of the reference state at the current temperature of the solution and the reference pressure for the species.
 
double z () const
 Calculate the value of z.
 
virtual double psatEst (double TKelvin) const
 Estimate for the saturation pressure.
 
int corr0 (double TKelvin, double pres, double &densLiq, double &densGas, double &liqGRT, double &gasGRT)
 Utility routine in the calculation of the saturation pressure.
 
- Protected Member Functions inherited from ThermoPhase
virtual void getParameters (AnyMap &phaseNode) const
 Store the parameters of a ThermoPhase object such that an identical one could be reconstructed using the newThermo(AnyMap&) function.
 
virtual void getCsvReportData (vector< string > &names, vector< vector< double > > &data) const
 Fills names and data with the column names and species thermo properties to be included in the output of the reportCSV method.
 
- Protected Member Functions inherited from Phase
void assertCompressible (const string &setter) const
 Ensure that phase is compressible.
 
void assignDensity (const double density_)
 Set the internally stored constant density (kg/m^3) of the phase.
 
void setMolecularWeight (const int k, const double mw)
 Set the molecular weight of a single species to a given value.
 
virtual void compositionChanged ()
 Apply changes to the state which are needed after the composition changes.
 

Protected Attributes

double m_b = 0.0
 Value of \( b \) in the equation of state.
 
double m_a = 0.0
 Value of \( a \) in the equation of state.
 
double m_aAlpha_mix = 0.0
 Value of \( a \alpha \) in the equation of state.
 
vector< double > m_b_coeffs
 
vector< double > m_kappa
 
vector< double > m_acentric
 acentric factor for each species, length m_kk
 
vector< double > m_dalphadT
 
vector< double > m_d2alphadT2
 
vector< double > m_alpha
 
Array2D m_a_coeffs
 
Array2D m_aAlpha_binary
 
map< string, map< string, double > > m_binaryParameters
 Explicitly-specified binary interaction parameters, to enable serialization.
 
int m_NSolns = 0
 
double m_Vroot [3] = {0.0, 0.0, 0.0}
 
vector< double > m_pp
 Temporary storage - length = m_kk.
 
vector< double > m_partialMolarVolumes
 
double m_dpdV = 0.0
 The derivative of the pressure with respect to the volume.
 
double m_dpdT = 0.0
 The derivative of the pressure with respect to the temperature.
 
vector< double > m_dpdni
 Vector of derivatives of pressure with respect to mole number.
 
vector< CoeffSource > m_coeffSource
 For each species, specifies the source of the a, b, and omega coefficients.
 
- Protected Attributes inherited from MixtureFugacityTP
vector< double > m_tmpV
 Temporary storage - length = m_kk.
 
vector< double > moleFractions_
 Storage for the current values of the mole fractions of the species.
 
int iState_ = FLUID_GAS
 Current state of the fluid.
 
int forcedState_ = FLUID_UNDEFINED
 Force the system to be on a particular side of the spinodal curve.
 
vector< double > m_h0_RT
 Temporary storage for dimensionless reference state enthalpies.
 
vector< double > m_cp0_R
 Temporary storage for dimensionless reference state heat capacities.
 
vector< double > m_g0_RT
 Temporary storage for dimensionless reference state Gibbs energies.
 
vector< double > m_s0_R
 Temporary storage for dimensionless reference state entropies.
 
- Protected Attributes inherited from ThermoPhase
MultiSpeciesThermo m_spthermo
 Pointer to the calculation manager for species reference-state thermodynamic properties.
 
AnyMap m_input
 Data supplied via setParameters.
 
double m_phi = 0.0
 Stored value of the electric potential for this phase. Units are Volts.
 
bool m_chargeNeutralityNecessary = false
 Boolean indicating whether a charge neutrality condition is a necessity.
 
int m_ssConvention = cSS_CONVENTION_TEMPERATURE
 Contains the standard state convention.
 
double m_tlast = 0.0
 last value of the temperature processed by reference state
 
- Protected Attributes inherited from Phase
ValueCache m_cache
 Cached for saved calculations within each ThermoPhase.
 
size_t m_kk = 0
 Number of species in the phase.
 
size_t m_ndim = 3
 Dimensionality of the phase.
 
vector< double > m_speciesComp
 Atomic composition of the species.
 
vector< double > m_speciesCharge
 Vector of species charges. length m_kk.
 
map< string, shared_ptr< Species > > m_species
 
UndefElement::behavior m_undefinedElementBehavior = UndefElement::add
 Flag determining behavior when adding species with an undefined element.
 
bool m_caseSensitiveSpecies = false
 Flag determining whether case sensitive species names are enforced.
 

Static Private Attributes

static const double omega_a = 4.5723552892138218E-01
 Omega constant: a0 (= omega_a) used in Peng-Robinson equation of state.
 
static const double omega_b = 7.77960739038885E-02
 Omega constant: b0 (= omega_b) used in Peng-Robinson equation of state.
 
static const double omega_vc = 3.07401308698703833E-01
 Omega constant for the critical molar volume.
 

Member Enumeration Documentation

◆ CoeffSource

enum class CoeffSource
strongprotected

Definition at line 305 of file PengRobinson.h.

Constructor & Destructor Documentation

◆ PengRobinson()

PengRobinson ( const string &  infile = "",
const string &  id = "" 
)
explicit

Construct and initialize a PengRobinson object directly from an input file.

Parameters
infileName of the input file containing the phase YAML data. If blank, an empty phase will be created.
idID of the phase in the input file. If empty, the first phase definition in the input file will be used.

Definition at line 24 of file PengRobinson.cpp.

Member Function Documentation

◆ type()

string type ( ) const
inlineoverridevirtual

String indicating the thermodynamic model implemented.

Usually corresponds to the name of the derived class, less any suffixes such as "Phase", TP", "VPSS", etc.

Since
Starting in Cantera 3.0, the name returned by this method corresponds to the canonical name used in the YAML input format.

Reimplemented from MixtureFugacityTP.

Definition at line 33 of file PengRobinson.h.

◆ cp_mole()

double cp_mole ( ) const
overridevirtual

Molar heat capacity at constant pressure. Units: J/kmol/K.

Reimplemented from ThermoPhase.

Definition at line 98 of file PengRobinson.cpp.

◆ cv_mole()

double cv_mole ( ) const
overridevirtual

Molar heat capacity at constant volume. Units: J/kmol/K.

Reimplemented from ThermoPhase.

Definition at line 112 of file PengRobinson.cpp.

◆ pressure()

double pressure ( ) const
overridevirtual

Return the thermodynamic pressure (Pa).

Since the mass density, temperature, and mass fractions are stored, this method uses these values to implement the mechanical equation of state \( P(T, \rho, Y_1, \dots, Y_K) \).

\[ P = \frac{RT}{v-b_{mix}} - \frac{\left(\alpha a\right)_{mix}}{v^2 + 2b_{mix}v - b_{mix}^2} \]

where:

\[ \alpha = \left[ 1 + \kappa \left(1-T_r^{0.5}\right)\right]^2 \]

and

\[ \kappa = \left(0.37464 + 1.54226\omega - 0.26992\omega^2\right), \qquad \qquad \text{For } \omega <= 0.491 \\ \kappa = \left(0.379642 + 1.487503\omega - 0.164423\omega^2 + 0.016667\omega^3 \right), \qquad \text{For } \omega > 0.491 \]

Coefficients \( a_{mix}, b_{mix} \) and \( (a \alpha)_{mix} \) are calculated as

\[ a_{mix} = \sum_i \sum_j X_i X_j a_{i, j} = \sum_i \sum_j X_i X_j \sqrt{a_i a_j} \]

\[ b_{mix} = \sum_i X_i b_i \]

\[ {a \alpha}_{mix} = \sum_i \sum_j X_i X_j {a \alpha}_{i, j} = \sum_i \sum_j X_i X_j \sqrt{a_i a_j} \sqrt{\alpha_i \alpha_j} \]

Reimplemented from Phase.

Definition at line 120 of file PengRobinson.cpp.

◆ standardConcentration()

double standardConcentration ( size_t  k = 0) const
overridevirtual

Returns the standard concentration \( C^0_k \), which is used to normalize the generalized concentration.

This is defined as the concentration by which the generalized concentration is normalized to produce the activity. The ideal gas mixture is considered as the standard or reference state here. Since the activity for an ideal gas mixture is simply the mole fraction, for an ideal gas, \( C^0_k = P/\hat R T \).

Parameters
kOptional parameter indicating the species. The default is to assume this refers to species 0.
Returns
Returns the standard Concentration in units of m^3 / kmol.

Reimplemented from ThermoPhase.

Definition at line 130 of file PengRobinson.cpp.

◆ getActivityCoefficients()

void getActivityCoefficients ( double *  ac) const
overridevirtual

Get the array of non-dimensional activity coefficients at the current solution temperature, pressure, and solution concentration.

For all objects with the Mixture Fugacity approximation, we define the standard state as an ideal gas at the current temperature and pressure of the solution. The activities are based on this standard state.

Parameters
acOutput vector of activity coefficients. Length: m_kk.

Reimplemented from ThermoPhase.

Definition at line 136 of file PengRobinson.cpp.

◆ getChemPotentials()

void getChemPotentials ( double *  mu) const
overridevirtual

Get the species chemical potentials. Units: J/kmol.

This function returns a vector of chemical potentials of the species in solution at the current temperature, pressure and mole fraction of the solution.

Parameters
muOutput vector of species chemical potentials. Length: m_kk. Units: J/kmol

Reimplemented from ThermoPhase.

Definition at line 169 of file PengRobinson.cpp.

◆ getPartialMolarEnthalpies()

void getPartialMolarEnthalpies ( double *  hbar) const
overridevirtual

Returns an array of partial molar enthalpies for the species in the mixture.

Units (J/kmol)

Parameters
hbarOutput vector of species partial molar enthalpies. Length: m_kk. units are J/kmol.

Reimplemented from ThermoPhase.

Definition at line 204 of file PengRobinson.cpp.

◆ getPartialMolarEntropies()

void getPartialMolarEntropies ( double *  sbar) const
overridevirtual

Returns an array of partial molar entropies of the species in the solution.

Units: J/kmol/K.

Parameters
sbarOutput vector of species partial molar entropies. Length = m_kk. units are J/kmol/K.

Reimplemented from ThermoPhase.

Definition at line 254 of file PengRobinson.cpp.

◆ getPartialMolarIntEnergies()

void getPartialMolarIntEnergies ( double *  ubar) const
overridevirtual

Return an array of partial molar internal energies for the species in the mixture.

Units: J/kmol.

Parameters
ubarOutput vector of species partial molar internal energies. Length = m_kk. units are J/kmol.

Reimplemented from ThermoPhase.

Definition at line 265 of file PengRobinson.cpp.

◆ getPartialMolarCp()

void getPartialMolarCp ( double *  cpbar) const
overridevirtual

Calculate species-specific molar specific heats.

This function is currently not implemented for Peng-Robinson phase.

Reimplemented from ThermoPhase.

Definition at line 276 of file PengRobinson.cpp.

◆ getPartialMolarVolumes()

void getPartialMolarVolumes ( double *  vbar) const
overridevirtual

Return an array of partial molar volumes for the species in the mixture.

Units: m^3/kmol.

Parameters
vbarOutput vector of species partial molar volumes. Length = m_kk. units are m^3/kmol.

Reimplemented from ThermoPhase.

Definition at line 281 of file PengRobinson.cpp.

◆ speciesCritTemperature()

double speciesCritTemperature ( double  a,
double  b 
) const

Calculate species-specific critical temperature.

The temperature dependent parameter in P-R EoS is calculated as

\[ T_{crit} = (0.0778 a)/(0.4572 b R) \]

Units: Kelvin

Parameters
aspecies-specific coefficients used in P-R EoS
bspecies-specific coefficients used in P-R EoS

Definition at line 307 of file PengRobinson.cpp.

◆ addSpecies()

bool addSpecies ( shared_ptr< Species spec)
overridevirtual

Add a Species to this Phase.

Returns true if the species was successfully added, or false if the species was ignored.

Derived classes which need to size arrays according to the number of species should overload this method. The derived class implementation should call the base class method, and, if this returns true (indicating that the species has been added), adjust their array sizes accordingly.

See also
ignoreUndefinedElements addUndefinedElements throwUndefinedElements

Reimplemented from MixtureFugacityTP.

Definition at line 318 of file PengRobinson.cpp.

◆ initThermo()

void initThermo ( )
overridevirtual

Initialize the ThermoPhase object after all species have been set up.

This method is provided to allow subclasses to perform any initialization required after all species have been added. For example, it might be used to resize internal work arrays that must have an entry for each species. The base class implementation does nothing, and subclasses that do not require initialization do not need to overload this method. Derived classes which do override this function should call their parent class's implementation of this function as their last action.

When importing from an AnyMap phase description (or from a YAML file), setupPhase() adds all the species, stores the input data in m_input, and then calls this method to set model parameters from the data stored in m_input.

Reimplemented from ThermoPhase.

Definition at line 338 of file PengRobinson.cpp.

◆ getSpeciesParameters()

void getSpeciesParameters ( const string &  name,
AnyMap speciesNode 
) const
overridevirtual

Get phase-specific parameters of a Species object such that an identical one could be reconstructed and added to this phase.

Parameters
nameName of the species
speciesNodeMapping to be populated with parameters

Reimplemented from ThermoPhase.

Definition at line 423 of file PengRobinson.cpp.

◆ setSpeciesCoeffs()

void setSpeciesCoeffs ( const string &  species,
double  a,
double  b,
double  w 
)

Set the pure fluid interaction parameters for a species.

Parameters
speciesName of the species
a\( a \) parameter in the Peng-Robinson model [Pa-m^6/kmol^2]
b\( a \) parameter in the Peng-Robinson model [m^3/kmol]
wacentric factor

Definition at line 29 of file PengRobinson.cpp.

◆ setBinaryCoeffs()

void setBinaryCoeffs ( const string &  species_i,
const string &  species_j,
double  a 
)

Set values for the interaction parameter between two species.

Parameters
species_iName of one species
species_jName of the other species
a\( a \) parameter in the Peng-Robinson model [Pa-m^6/kmol^2]

Definition at line 74 of file PengRobinson.cpp.

◆ sresid()

double sresid ( ) const
overrideprotectedvirtual

Calculate the deviation terms for the total entropy of the mixture from the ideal gas mixture.

Here we use the current state conditions

Returns
the change in entropy in units of J kmol-1 K-1.

Reimplemented from MixtureFugacityTP.

Definition at line 459 of file PengRobinson.cpp.

◆ hresid()

double hresid ( ) const
overrideprotectedvirtual

Calculate the deviation terms for the total enthalpy of the mixture from the ideal gas mixture.

Here we use the current state conditions

Returns
the change in entropy in units of J kmol-1.

Reimplemented from MixtureFugacityTP.

Definition at line 472 of file PengRobinson.cpp.

◆ liquidVolEst()

double liquidVolEst ( double  TKelvin,
double &  pres 
) const
overrideprotectedvirtual

Estimate for the molar volume of the liquid.

Note: this is only used as a starting guess for later routines that actually calculate an accurate value for the liquid molar volume. This routine doesn't change the state of the system.

Parameters
TKelvintemperature in kelvin
presPressure in Pa. This is used as an initial guess. If the routine needs to change the pressure to find a stable liquid state, the new pressure is returned in this variable.
Returns
the estimate of the liquid volume. If the liquid can't be found, this routine returns -1.

Reimplemented from MixtureFugacityTP.

Definition at line 485 of file PengRobinson.cpp.

◆ densityCalc()

double densityCalc ( double  TKelvin,
double  pressure,
int  phaseRequested,
double  rhoguess 
)
overrideprotectedvirtual

Calculates the density given the temperature and the pressure and a guess at the density.

Note, below T_c, this is a multivalued function. We do not cross the vapor dome in this. This is protected because it is called during setState_TP() routines. Infinite loops would result if it were not protected.

Parameters
TKelvinTemperature in Kelvin
pressurePressure in Pascals (Newton/m**2)
phaseRequestedint representing the phase whose density we are requesting. If we put a gas or liquid phase here, we will attempt to find a volume in that part of the volume space, only, in this routine. A value of FLUID_UNDEFINED means that we will accept anything.
rhoguessGuessed density of the fluid. A value of -1.0 indicates that there is no guessed density
Returns
We return the density of the fluid at the requested phase. If we have not found any acceptable density we return a -1. If we have found an acceptable density at a different phase, we return a -2.

Reimplemented from MixtureFugacityTP.

Definition at line 516 of file PengRobinson.cpp.

◆ densSpinodalLiquid()

double densSpinodalLiquid ( ) const
overrideprotectedvirtual

Return the value of the density at the liquid spinodal point (on the liquid side) for the current temperature.

Returns
the density with units of kg m-3

Reimplemented from MixtureFugacityTP.

Definition at line 574 of file PengRobinson.cpp.

◆ densSpinodalGas()

double densSpinodalGas ( ) const
overrideprotectedvirtual

Return the value of the density at the gas spinodal point (on the gas side) for the current temperature.

Returns
the density with units of kg m-3

Reimplemented from MixtureFugacityTP.

Definition at line 596 of file PengRobinson.cpp.

◆ dpdVCalc()

double dpdVCalc ( double  TKelvin,
double  molarVol,
double &  presCalc 
) const
overrideprotectedvirtual

Calculate the pressure and the pressure derivative given the temperature and the molar volume.

Temperature and mole number are held constant

Parameters
TKelvintemperature in kelvin
molarVolmolar volume ( m3/kmol)
presCalcReturns the pressure.
Returns
the derivative of the pressure wrt the molar volume

Reimplemented from MixtureFugacityTP.

Definition at line 618 of file PengRobinson.cpp.

◆ daAlpha_dT()

double daAlpha_dT ( ) const
protected

Calculate temperature derivative \( d(a \alpha)/dT \).

These are stored internally.

Definition at line 690 of file PengRobinson.cpp.

◆ d2aAlpha_dT2()

double d2aAlpha_dT2 ( ) const
protected

Calculate second derivative \( d^2(a \alpha)/dT^2 \).

These are stored internally.

Definition at line 713 of file PengRobinson.cpp.

◆ isothermalCompressibility()

double isothermalCompressibility ( ) const
overridevirtual

Returns the isothermal compressibility. Units: 1/Pa.

The isothermal compressibility is defined as

\[ \kappa_T = -\frac{1}{v}\left(\frac{\partial v}{\partial P}\right)_T \]

or

\[ \kappa_T = \frac{1}{\rho}\left(\frac{\partial \rho}{\partial P}\right)_T \]

Reimplemented from ThermoPhase.

Definition at line 626 of file PengRobinson.cpp.

◆ thermalExpansionCoeff()

double thermalExpansionCoeff ( ) const
overridevirtual

Return the volumetric thermal expansion coefficient. Units: 1/K.

The thermal expansion coefficient is defined as

\[ \beta = \frac{1}{v}\left(\frac{\partial v}{\partial T}\right)_P \]

Reimplemented from ThermoPhase.

Definition at line 632 of file PengRobinson.cpp.

◆ soundSpeed()

double soundSpeed ( ) const
overridevirtual

Return the speed of sound. Units: m/s.

The speed of sound is defined as

\[ c = \sqrt{\left(\frac{\partial P}{\partial\rho}\right)_s} \]

Reimplemented from ThermoPhase.

Definition at line 638 of file PengRobinson.cpp.

◆ calculatePressureDerivatives()

void calculatePressureDerivatives ( ) const

Calculate \( dp/dV \) and \( dp/dT \) at the current conditions.

These are stored internally.

Definition at line 644 of file PengRobinson.cpp.

◆ updateMixingExpressions()

void updateMixingExpressions ( )
overridevirtual

Update the \( a \), \( b \), and \( \alpha \) parameters.

The \( a \) and the \( b \) parameters depend on the mole fraction and the parameter \( \alpha \) depends on the temperature. This function updates the internal numbers based on the state of the object.

Reimplemented from MixtureFugacityTP.

Definition at line 656 of file PengRobinson.cpp.

◆ calculateAB()

void calculateAB ( double &  aCalc,
double &  bCalc,
double &  aAlpha 
) const

Calculate the \( a \), \( b \), and \( \alpha \) parameters given the temperature.

This function doesn't change the internal state of the object, so it is a const function. It does use the stored mole fractions in the object.

Parameters
aCalc(output) Returns the a value
bCalc(output) Returns the b value.
aAlpha(output) Returns the (a*alpha) value.

Definition at line 676 of file PengRobinson.cpp.

◆ calcCriticalConditions()

void calcCriticalConditions ( double &  pc,
double &  tc,
double &  vc 
) const
overridevirtual

Reimplemented from MixtureFugacityTP.

Definition at line 744 of file PengRobinson.cpp.

◆ solveCubic()

int solveCubic ( double  T,
double  pres,
double  a,
double  b,
double  aAlpha,
double  Vroot[3] 
) const

Prepare variables and call the function to solve the cubic equation of state.

Definition at line 763 of file PengRobinson.cpp.

Member Data Documentation

◆ m_b

double m_b = 0.0
protected

Value of \( b \) in the equation of state.

m_b is a function of the mole fractions and species-specific b values.

Definition at line 243 of file PengRobinson.h.

◆ m_a

double m_a = 0.0
protected

Value of \( a \) in the equation of state.

m_a depends only on the mole fractions.

Definition at line 249 of file PengRobinson.h.

◆ m_aAlpha_mix

double m_aAlpha_mix = 0.0
protected

Value of \( a \alpha \) in the equation of state.

m_aAlpha_mix is a function of the temperature and the mole fractions.

Definition at line 255 of file PengRobinson.h.

◆ m_b_coeffs

vector<double> m_b_coeffs
protected

Definition at line 259 of file PengRobinson.h.

◆ m_kappa

vector<double> m_kappa
protected

Definition at line 260 of file PengRobinson.h.

◆ m_acentric

vector<double> m_acentric
protected

acentric factor for each species, length m_kk

Definition at line 261 of file PengRobinson.h.

◆ m_dalphadT

vector<double> m_dalphadT
mutableprotected

Definition at line 262 of file PengRobinson.h.

◆ m_d2alphadT2

vector<double> m_d2alphadT2
mutableprotected

Definition at line 263 of file PengRobinson.h.

◆ m_alpha

vector<double> m_alpha
protected

Definition at line 264 of file PengRobinson.h.

◆ m_a_coeffs

Array2D m_a_coeffs
protected

Definition at line 268 of file PengRobinson.h.

◆ m_aAlpha_binary

Array2D m_aAlpha_binary
protected

Definition at line 269 of file PengRobinson.h.

◆ m_binaryParameters

map<string, map<string, double> > m_binaryParameters
protected

Explicitly-specified binary interaction parameters, to enable serialization.

Definition at line 272 of file PengRobinson.h.

◆ m_NSolns

int m_NSolns = 0
protected

Definition at line 274 of file PengRobinson.h.

◆ m_Vroot

double m_Vroot[3] = {0.0, 0.0, 0.0}
protected

Definition at line 276 of file PengRobinson.h.

◆ m_pp

vector<double> m_pp
mutableprotected

Temporary storage - length = m_kk.

Definition at line 279 of file PengRobinson.h.

◆ m_partialMolarVolumes

vector<double> m_partialMolarVolumes
mutableprotected

Definition at line 282 of file PengRobinson.h.

◆ m_dpdV

double m_dpdV = 0.0
mutableprotected

The derivative of the pressure with respect to the volume.

Calculated at the current conditions. temperature and mole number kept constant

Definition at line 289 of file PengRobinson.h.

◆ m_dpdT

double m_dpdT = 0.0
mutableprotected

The derivative of the pressure with respect to the temperature.

Calculated at the current conditions. Total volume and mole number kept constant

Definition at line 296 of file PengRobinson.h.

◆ m_dpdni

vector<double> m_dpdni
mutableprotected

Vector of derivatives of pressure with respect to mole number.

Calculated at the current conditions. Total volume, temperature and other mole number kept constant

Definition at line 303 of file PengRobinson.h.

◆ m_coeffSource

vector<CoeffSource> m_coeffSource
protected

For each species, specifies the source of the a, b, and omega coefficients.

Definition at line 307 of file PengRobinson.h.

◆ omega_a

const double omega_a = 4.5723552892138218E-01
staticprivate

Omega constant: a0 (= omega_a) used in Peng-Robinson equation of state.

This value is calculated by solving P-R cubic equation at the critical point.

Definition at line 314 of file PengRobinson.h.

◆ omega_b

const double omega_b = 7.77960739038885E-02
staticprivate

Omega constant: b0 (= omega_b) used in Peng-Robinson equation of state.

This value is calculated by solving P-R cubic equation at the critical point.

Definition at line 320 of file PengRobinson.h.

◆ omega_vc

const double omega_vc = 3.07401308698703833E-01
staticprivate

Omega constant for the critical molar volume.

Definition at line 323 of file PengRobinson.h.


The documentation for this class was generated from the following files: