Cantera 2.6.0
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
MixtureFugacityTP Class Reference

This is a filter class for ThermoPhase that implements some preparatory steps for efficiently handling mixture of gases that whose standard states are defined as ideal gases, but which describe also non-ideal solutions. More...

#include <MixtureFugacityTP.h>

Inheritance diagram for MixtureFugacityTP:
[legend]
Collaboration diagram for MixtureFugacityTP:
[legend]

Public Member Functions

virtual void setTemperature (const doublereal temp)
 Set the temperature of the phase. More...
 
virtual void setPressure (doublereal p)
 Set the internally stored pressure (Pa) at constant temperature and composition. More...
 
Constructors and Duplicators for MixtureFugacityTP
 MixtureFugacityTP ()
 Constructor. More...
 
Utilities
virtual std::string type () const
 String indicating the thermodynamic model implemented. More...
 
virtual int standardStateConvention () const
 This method returns the convention used in specification of the standard state, of which there are currently two, temperature based, and variable pressure based. More...
 
virtual void setForcedSolutionBranch (int solnBranch)
 Set the solution branch to force the ThermoPhase to exist on one branch or another. More...
 
virtual int forcedSolutionBranch () const
 Report the solution branch which the solution is restricted to. More...
 
virtual int reportSolnBranchActual () const
 Report the solution branch which the solution is actually on. More...
 
virtual void getdlnActCoeffdlnN_diag (doublereal *dlnActCoeffdlnN_diag) const
 Get the array of log species mole number derivatives of the log activity coefficients. More...
 
Molar Thermodynamic properties
virtual double enthalpy_mole () const
 Molar enthalpy. Units: J/kmol. More...
 
virtual double entropy_mole () const
 Molar entropy. Units: J/kmol/K. More...
 
Partial Molar Properties of the Solution
virtual void getChemPotentials_RT (doublereal *mu) const
 Get the array of non-dimensional species chemical potentials These are partial molar Gibbs free energies. More...
 
Properties of the Standard State of the Species in the Solution

Within MixtureFugacityTP, these properties are calculated via a common routine, _updateStandardStateThermo(), which must be overloaded in inherited objects. The values are cached within this object, and are not recalculated unless the temperature or pressure changes.

virtual void getStandardChemPotentials (doublereal *mu) const
 Get the array of chemical potentials at unit activity. More...
 
virtual void getEnthalpy_RT (doublereal *hrt) const
 Get the nondimensional Enthalpy functions for the species at their standard states at the current T and P of the solution. More...
 
virtual void getEntropy_R (doublereal *sr) const
 Get the array of nondimensional Enthalpy functions for the standard state species at the current T and P of the solution. More...
 
virtual void getGibbs_RT (doublereal *grt) const
 Get the nondimensional Gibbs functions for the species at their standard states of solution at the current T and P of the solution. More...
 
virtual void getPureGibbs (doublereal *gpure) const
 Get the pure Gibbs free energies of each species. More...
 
virtual void getIntEnergy_RT (doublereal *urt) const
 Returns the vector of nondimensional internal Energies of the standard state at the current temperature and pressure of the solution for each species. More...
 
virtual void getCp_R (doublereal *cpr) const
 Get the nondimensional Heat Capacities at constant pressure for the standard state of the species at the current T and P. More...
 
virtual void getStandardVolumes (doublereal *vol) const
 Get the molar volumes of each species in their standard states at the current T and P of the solution. More...
 
Initialization Methods - For Internal use
virtual bool addSpecies (shared_ptr< Species > spec)
 
virtual void setStateFromXML (const XML_Node &state)
 Set the initial state of the phase to the conditions specified in the state XML element. More...
 
- Public Member Functions inherited from ThermoPhase
 ThermoPhase ()
 Constructor. More...
 
doublereal RT () const
 Return the Gas Constant multiplied by the current temperature. More...
 
double equivalenceRatio () const
 Compute the equivalence ratio for the current mixture from available oxygen and required oxygen. More...
 
virtual bool isIdeal () const
 Boolean indicating whether phase is ideal. More...
 
virtual std::string phaseOfMatter () const
 String indicating the mechanical phase of the matter in this Phase. More...
 
virtual doublereal refPressure () const
 Returns the reference pressure in Pa. More...
 
virtual doublereal minTemp (size_t k=npos) const
 Minimum temperature for which the thermodynamic data for the species or phase are valid. More...
 
doublereal Hf298SS (const size_t k) const
 Report the 298 K Heat of Formation of the standard state of one species (J kmol-1) More...
 
virtual void modifyOneHf298SS (const size_t k, const doublereal Hf298New)
 Modify the value of the 298 K Heat of Formation of one species in the phase (J kmol-1) More...
 
virtual void resetHf298 (const size_t k=npos)
 Restore the original heat of formation of one or more species. More...
 
virtual doublereal maxTemp (size_t k=npos) const
 Maximum temperature for which the thermodynamic data for the species are valid. More...
 
bool chargeNeutralityNecessary () const
 Returns the chargeNeutralityNecessity boolean. More...
 
virtual doublereal intEnergy_mole () const
 Molar internal energy. Units: J/kmol. More...
 
virtual doublereal gibbs_mole () const
 Molar Gibbs function. Units: J/kmol. More...
 
virtual doublereal cp_mole () const
 Molar heat capacity at constant pressure. Units: J/kmol/K. More...
 
virtual doublereal cv_mole () const
 Molar heat capacity at constant volume. Units: J/kmol/K. More...
 
virtual doublereal isothermalCompressibility () const
 Returns the isothermal compressibility. Units: 1/Pa. More...
 
virtual doublereal thermalExpansionCoeff () const
 Return the volumetric thermal expansion coefficient. Units: 1/K. More...
 
void setElectricPotential (doublereal v)
 Set the electric potential of this phase (V). More...
 
doublereal electricPotential () const
 Returns the electric potential of this phase (V). More...
 
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. More...
 
virtual Units standardConcentrationUnits () const
 Returns the units of the "standard concentration" for this phase. More...
 
virtual doublereal standardConcentration (size_t k=0) const
 Return the standard concentration for the kth species. More...
 
virtual doublereal logStandardConc (size_t k=0) const
 Natural logarithm of the standard concentration of the kth species. More...
 
virtual void getActivities (doublereal *a) const
 Get the array of non-dimensional activities at the current solution temperature, pressure, and solution concentration. More...
 
virtual void getActivityCoefficients (doublereal *ac) const
 Get the array of non-dimensional molar-based activity coefficients at the current solution temperature, pressure, and solution concentration. More...
 
virtual void getLnActivityCoefficients (doublereal *lnac) const
 Get the array of non-dimensional molar-based ln activity coefficients at the current solution temperature, pressure, and solution concentration. More...
 
virtual void getChemPotentials (doublereal *mu) const
 Get the species chemical potentials. Units: J/kmol. More...
 
void getElectrochemPotentials (doublereal *mu) const
 Get the species electrochemical potentials. More...
 
virtual void getPartialMolarEnthalpies (doublereal *hbar) const
 Returns an array of partial molar enthalpies for the species in the mixture. More...
 
virtual void getPartialMolarEntropies (doublereal *sbar) const
 Returns an array of partial molar entropies of the species in the solution. More...
 
virtual void getPartialMolarIntEnergies (doublereal *ubar) const
 Return an array of partial molar internal energies for the species in the mixture. More...
 
virtual void getPartialMolarCp (doublereal *cpbar) const
 Return an array of partial molar heat capacities for the species in the mixture. More...
 
virtual void getPartialMolarVolumes (doublereal *vbar) const
 Return an array of partial molar volumes for the species in the mixture. More...
 
virtual void getIntEnergy_RT_ref (doublereal *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. More...
 
doublereal enthalpy_mass () const
 Specific enthalpy. Units: J/kg. More...
 
doublereal intEnergy_mass () const
 Specific internal energy. Units: J/kg. More...
 
doublereal entropy_mass () const
 Specific entropy. Units: J/kg/K. More...
 
doublereal gibbs_mass () const
 Specific Gibbs function. Units: J/kg. More...
 
doublereal cp_mass () const
 Specific heat at constant pressure. Units: J/kg/K. More...
 
doublereal cv_mass () const
 Specific heat at constant volume. Units: J/kg/K. More...
 
virtual void setState_TPX (doublereal t, doublereal p, const doublereal *x)
 Set the temperature (K), pressure (Pa), and mole fractions. More...
 
virtual void setState_TPX (doublereal t, doublereal p, const compositionMap &x)
 Set the temperature (K), pressure (Pa), and mole fractions. More...
 
virtual void setState_TPX (doublereal t, doublereal p, const std::string &x)
 Set the temperature (K), pressure (Pa), and mole fractions. More...
 
virtual void setState_TPY (doublereal t, doublereal p, const doublereal *y)
 Set the internally stored temperature (K), pressure (Pa), and mass fractions of the phase. More...
 
virtual void setState_TPY (doublereal t, doublereal p, const compositionMap &y)
 Set the internally stored temperature (K), pressure (Pa), and mass fractions of the phase. More...
 
virtual void setState_TPY (doublereal t, doublereal p, const std::string &y)
 Set the internally stored temperature (K), pressure (Pa), and mass fractions of the phase. More...
 
virtual void setState_TP (doublereal t, doublereal p)
 Set the temperature (K) and pressure (Pa) More...
 
virtual void setState_PX (doublereal p, doublereal *x)
 Set the pressure (Pa) and mole fractions. More...
 
virtual void setState_PY (doublereal p, doublereal *y)
 Set the internally stored pressure (Pa) and mass fractions. More...
 
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. More...
 
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). More...
 
virtual void setState_SP (double s, double p, double tol=1e-9)
 Set the specific entropy (J/kg/K) and pressure (Pa). More...
 
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). More...
 
virtual void setState_ST (double s, double t, double tol=1e-9)
 Set the specific entropy (J/kg/K) and temperature (K). More...
 
virtual void setState_TV (double t, double v, double tol=1e-9)
 Set the temperature (K) and specific volume (m^3/kg). More...
 
virtual void setState_PV (double p, double v, double tol=1e-9)
 Set the pressure (Pa) and specific volume (m^3/kg). More...
 
virtual void setState_UP (double u, double p, double tol=1e-9)
 Set the specific internal energy (J/kg) and pressure (Pa). More...
 
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) More...
 
virtual void setState_TH (double t, double h, double tol=1e-9)
 Set the temperature (K) and the specific enthalpy (J/kg) More...
 
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) More...
 
virtual void setState_RP (doublereal rho, doublereal p)
 Set the density (kg/m**3) and pressure (Pa) at constant composition. More...
 
virtual void setState_RPX (doublereal rho, doublereal p, const doublereal *x)
 Set the density (kg/m**3), pressure (Pa) and mole fractions. More...
 
virtual void setState_RPX (doublereal rho, doublereal p, const compositionMap &x)
 Set the density (kg/m**3), pressure (Pa) and mole fractions. More...
 
virtual void setState_RPX (doublereal rho, doublereal p, const std::string &x)
 Set the density (kg/m**3), pressure (Pa) and mole fractions. More...
 
virtual void setState_RPY (doublereal rho, doublereal p, const doublereal *y)
 Set the density (kg/m**3), pressure (Pa) and mass fractions. More...
 
virtual void setState_RPY (doublereal rho, doublereal p, const compositionMap &y)
 Set the density (kg/m**3), pressure (Pa) and mass fractions. More...
 
virtual void setState_RPY (doublereal rho, doublereal p, const std::string &y)
 Set the density (kg/m**3), pressure (Pa) and mass fractions. More...
 
virtual void setState (const AnyMap &state)
 Set the state using an AnyMap containing any combination of properties supported by the thermodynamic model. More...
 
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) More...
 
void setMixtureFraction (double mixFrac, const std::string &fuelComp, const std::string &oxComp, ThermoBasis basis=ThermoBasis::molar)
 Set the mixture composition according to the mixture fraction = kg fuel / (kg oxidizer + kg fuel) More...
 
void setMixtureFraction (double mixFrac, const compositionMap &fuelComp, const compositionMap &oxComp, ThermoBasis basis=ThermoBasis::molar)
 Set the mixture composition according to the mixture fraction = kg fuel / (kg oxidizer + kg fuel) More...
 
double mixtureFraction (const double *fuelComp, const double *oxComp, ThermoBasis basis=ThermoBasis::molar, const std::string &element="Bilger") const
 Compute the mixture fraction = kg fuel / (kg oxidizer + kg fuel) for the current mixture given fuel and oxidizer compositions. More...
 
double mixtureFraction (const std::string &fuelComp, const std::string &oxComp, ThermoBasis basis=ThermoBasis::molar, const std::string &element="Bilger") const
 Compute the mixture fraction = kg fuel / (kg oxidizer + kg fuel) for the current mixture given fuel and oxidizer compositions. More...
 
double mixtureFraction (const compositionMap &fuelComp, const compositionMap &oxComp, ThermoBasis basis=ThermoBasis::molar, const std::string &element="Bilger") const
 Compute the mixture fraction = kg fuel / (kg oxidizer + kg fuel) for the current mixture given fuel and oxidizer compositions. More...
 
void setEquivalenceRatio (double phi, const double *fuelComp, const double *oxComp, ThermoBasis basis=ThermoBasis::molar)
 Set the mixture composition according to the equivalence ratio. More...
 
void setEquivalenceRatio (double phi, const std::string &fuelComp, const std::string &oxComp, ThermoBasis basis=ThermoBasis::molar)
 Set the mixture composition according to the equivalence ratio. More...
 
void setEquivalenceRatio (double phi, const compositionMap &fuelComp, const compositionMap &oxComp, ThermoBasis basis=ThermoBasis::molar)
 Set the mixture composition according to the equivalence ratio. More...
 
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. More...
 
double equivalenceRatio (const std::string &fuelComp, const std::string &oxComp, ThermoBasis basis=ThermoBasis::molar) const
 Compute the equivalence ratio for the current mixture given the compositions of fuel and oxidizer. More...
 
double equivalenceRatio (const compositionMap &fuelComp, const compositionMap &oxComp, ThermoBasis basis=ThermoBasis::molar) const
 Compute the equivalence ratio for the current mixture given the compositions of fuel and oxidizer. More...
 
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. More...
 
double stoichAirFuelRatio (const std::string &fuelComp, const std::string &oxComp, ThermoBasis basis=ThermoBasis::molar) const
 Compute the stoichiometric air to fuel ratio (kg oxidizer / kg fuel) given fuel and oxidizer compositions. More...
 
double stoichAirFuelRatio (const compositionMap &fuelComp, const compositionMap &oxComp, ThermoBasis basis=ThermoBasis::molar) const
 Compute the stoichiometric air to fuel ratio (kg oxidizer / kg fuel) given fuel and oxidizer compositions. More...
 
void equilibrate (const std::string &XY, const std::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. More...
 
virtual void setToEquilState (const doublereal *mu_RT)
 This method is used by the ChemEquil equilibrium solver. More...
 
virtual bool compatibleWithMultiPhase () const
 Indicates whether this phase type can be used with class MultiPhase for equilibrium calculations. More...
 
virtual doublereal satTemperature (doublereal p) const
 Return the saturation temperature given the pressure. More...
 
virtual doublereal vaporFraction () const
 Return the fraction of vapor at the current conditions. More...
 
virtual void setState_Tsat (doublereal t, doublereal x)
 Set the state to a saturated system at a particular temperature. More...
 
virtual void setState_Psat (doublereal p, doublereal x)
 Set the state to a saturated system at a particular pressure. More...
 
void setState_TPQ (double T, double P, double Q)
 Set the temperature, pressure, and vapor fraction (quality). More...
 
virtual void modifySpecies (size_t k, shared_ptr< Species > spec)
 Modify the thermodynamic data associated with a species. More...
 
void saveSpeciesData (const size_t k, const XML_Node *const data)
 Store a reference pointer to the XML tree containing the species data for this phase. More...
 
const std::vector< const XML_Node * > & speciesData () const
 Return a pointer to the vector of XML nodes containing the species data for this phase. More...
 
virtual MultiSpeciesThermospeciesThermo (int k=-1)
 Return a changeable reference to the calculation manager for species reference-state thermodynamic properties. More...
 
virtual const MultiSpeciesThermospeciesThermo (int k=-1) const
 
void initThermoFile (const std::string &inputFile, const std::string &id)
 
virtual void initThermoXML (XML_Node &phaseNode, const std::string &id)
 Import and initialize a ThermoPhase object using an XML tree. More...
 
virtual void initThermo ()
 Initialize the ThermoPhase object after all species have been set up. More...
 
virtual void setParameters (int n, doublereal *const c)
 Set the equation of state parameters. More...
 
virtual void getParameters (int &n, doublereal *const c) const
 Get the equation of state parameters in a vector. More...
 
virtual void setParameters (const AnyMap &phaseNode, const AnyMap &rootNode=AnyMap())
 Set equation of state parameters from an AnyMap phase description. More...
 
AnyMap parameters (bool withInput=true) const
 Returns the parameters of a ThermoPhase object such that an identical one could be reconstructed using the newPhase(AnyMap&) function. More...
 
virtual void getSpeciesParameters (const std::string &name, AnyMap &speciesNode) const
 Get phase-specific parameters of a Species object such that an identical one could be reconstructed and added to this phase. More...
 
const AnyMapinput () const
 Access input data associated with the phase description. More...
 
AnyMapinput ()
 
virtual void setParametersFromXML (const XML_Node &eosdata)
 Set equation of state parameter values from XML entries. More...
 
virtual void invalidateCache ()
 Invalidate any cached values which are normally updated only when a change in state is detected. More...
 
virtual void getdlnActCoeffds (const doublereal dTds, const doublereal *const dXds, doublereal *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. More...
 
virtual void getdlnActCoeffdlnX_diag (doublereal *dlnActCoeffdlnX_diag) const
 Get the array of ln mole fraction derivatives of the log activity coefficients - diagonal component only. More...
 
virtual void getdlnActCoeffdlnN (const size_t ld, doublereal *const dlnActCoeffdlnN)
 Get the array of derivatives of the log activity coefficients with respect to the log of the species mole numbers. More...
 
virtual void getdlnActCoeffdlnN_numderiv (const size_t ld, doublereal *const dlnActCoeffdlnN)
 
virtual std::string report (bool show_thermo=true, doublereal threshold=-1e-14) const
 returns a summary of the state of the phase as a string More...
 
virtual void reportCSV (std::ofstream &csvFile) const
 returns a summary of the state of the phase to a comma separated file. More...
 
- Public Member Functions inherited from Phase
 Phase ()
 Default constructor. More...
 
 Phase (const Phase &)=delete
 
Phaseoperator= (const Phase &)=delete
 
XML_Nodexml () const
 Returns a const reference to the XML_Node that describes the phase. More...
 
void setXMLdata (XML_Node &xmlPhase)
 Stores the XML tree information for the current phase. More...
 
virtual bool isPure () const
 Return whether phase represents a pure (single species) substance. More...
 
virtual bool hasPhaseTransition () const
 Return whether phase represents a substance with phase transitions. More...
 
virtual bool isCompressible () const
 Return whether phase represents a compressible substance. More...
 
virtual std::map< std::string, size_t > nativeState () const
 Return a map of properties defining the native state of a substance. More...
 
virtual std::vector< std::string > fullStates () const
 Return a vector containing full states defining a phase. More...
 
virtual std::vector< std::string > partialStates () const
 Return a vector of settable partial property sets within a phase. More...
 
virtual size_t stateSize () const
 Return size of vector defining internal state of the phase. More...
 
void saveState (vector_fp &state) const
 Save the current internal state of the phase. More...
 
virtual void saveState (size_t lenstate, doublereal *state) const
 Write to array 'state' the current internal state. More...
 
void restoreState (const vector_fp &state)
 Restore a state saved on a previous call to saveState. More...
 
virtual void restoreState (size_t lenstate, const doublereal *state)
 Restore the state of the phase from a previously saved state vector. More...
 
doublereal molecularWeight (size_t k) const
 Molecular weight of species k. More...
 
void getMolecularWeights (vector_fp &weights) const
 Copy the vector of molecular weights into vector weights. More...
 
void getMolecularWeights (doublereal *weights) const
 Copy the vector of molecular weights into array weights. More...
 
const vector_fpmolecularWeights () const
 Return a const reference to the internal vector of molecular weights. More...
 
void getCharges (double *charges) const
 Copy the vector of species charges into array charges. More...
 
doublereal elementalMassFraction (const size_t m) const
 Elemental mass fraction of element m. More...
 
doublereal elementalMoleFraction (const size_t m) const
 Elemental mole fraction of element m. More...
 
const double * moleFractdivMMW () const
 Returns a const pointer to the start of the moleFraction/MW array. More...
 
doublereal 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. More...
 
doublereal chargeDensity () const
 Charge density [C/m^3]. More...
 
size_t nDim () const
 Returns the number of spatial dimensions (1, 2, or 3) More...
 
void setNDim (size_t ndim)
 Set the number of spatial dimensions (1, 2, or 3). More...
 
virtual bool ready () const
 Returns a bool indicating whether the object is ready for use. More...
 
int stateMFNumber () const
 Return the State Mole Fraction Number. More...
 
bool caseSensitiveSpecies () const
 Returns true if case sensitive species names are enforced. More...
 
void setCaseSensitiveSpecies (bool cflag=true)
 Set flag that determines whether case sensitive species are enforced in look-up operations, for example speciesIndex. More...
 
virtual void setRoot (std::shared_ptr< Solution > root)
 Set root Solution holding all phase information. More...
 
vector_fp getCompositionFromMap (const compositionMap &comp) const
 Converts a compositionMap to a vector with entries for each species Species that are not specified are set to zero in the vector. More...
 
void massFractionsToMoleFractions (const double *Y, double *X) const
 Converts a mixture composition from mole fractions to mass fractions. More...
 
void moleFractionsToMassFractions (const double *X, double *Y) const
 Converts a mixture composition from mass fractions to mole fractions. More...
 
std::string name () const
 Return the name of the phase. More...
 
void setName (const std::string &nm)
 Sets the string name for the phase. More...
 
std::string elementName (size_t m) const
 Name of the element with index m. More...
 
size_t elementIndex (const std::string &name) const
 Return the index of element named 'name'. More...
 
const std::vector< std::string > & elementNames () const
 Return a read-only reference to the vector of element names. More...
 
doublereal atomicWeight (size_t m) const
 Atomic weight of element m. More...
 
doublereal entropyElement298 (size_t m) const
 Entropy of the element in its standard state at 298 K and 1 bar. More...
 
int atomicNumber (size_t m) const
 Atomic number of element m. More...
 
int elementType (size_t m) const
 Return the element constraint type Possible types include: More...
 
int changeElementType (int m, int elem_type)
 Change the element type of the mth constraint Reassigns an element type. More...
 
const vector_fpatomicWeights () const
 Return a read-only reference to the vector of atomic weights. More...
 
size_t nElements () const
 Number of elements. More...
 
void checkElementIndex (size_t m) const
 Check that the specified element index is in range. More...
 
void checkElementArraySize (size_t mm) const
 Check that an array size is at least nElements(). More...
 
doublereal nAtoms (size_t k, size_t m) const
 Number of atoms of element m in species k. More...
 
void getAtoms (size_t k, double *atomArray) const
 Get a vector containing the atomic composition of species k. More...
 
size_t speciesIndex (const std::string &name) const
 Returns the index of a species named 'name' within the Phase object. More...
 
std::string speciesName (size_t k) const
 Name of the species with index k. More...
 
std::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. More...
 
const std::vector< std::string > & speciesNames () const
 Return a const reference to the vector of species names. More...
 
size_t nSpecies () const
 Returns the number of species in the phase. More...
 
void checkSpeciesIndex (size_t k) const
 Check that the specified species index is in range. More...
 
void checkSpeciesArraySize (size_t kk) const
 Check that an array size is at least nSpecies(). More...
 
void setMoleFractionsByName (const compositionMap &xMap)
 Set the species mole fractions by name. More...
 
void setMoleFractionsByName (const std::string &x)
 Set the mole fractions of a group of species by name. More...
 
void setMassFractionsByName (const compositionMap &yMap)
 Set the species mass fractions by name. More...
 
void setMassFractionsByName (const std::string &x)
 Set the species mass fractions by name. More...
 
void setState_TRX (doublereal t, doublereal dens, const doublereal *x)
 Set the internally stored temperature (K), density, and mole fractions. More...
 
void setState_TRX (doublereal t, doublereal dens, const compositionMap &x)
 Set the internally stored temperature (K), density, and mole fractions. More...
 
void setState_TRY (doublereal t, doublereal dens, const doublereal *y)
 Set the internally stored temperature (K), density, and mass fractions. More...
 
void setState_TRY (doublereal t, doublereal dens, const compositionMap &y)
 Set the internally stored temperature (K), density, and mass fractions. More...
 
void setState_TNX (doublereal t, doublereal n, const doublereal *x)
 Set the internally stored temperature (K), molar density (kmol/m^3), and mole fractions. More...
 
void setState_TR (doublereal t, doublereal rho)
 Set the internally stored temperature (K) and density (kg/m^3) More...
 
void setState_TX (doublereal t, doublereal *x)
 Set the internally stored temperature (K) and mole fractions. More...
 
void setState_TY (doublereal t, doublereal *y)
 Set the internally stored temperature (K) and mass fractions. More...
 
void setState_RX (doublereal rho, doublereal *x)
 Set the density (kg/m^3) and mole fractions. More...
 
void setState_RY (doublereal rho, doublereal *y)
 Set the density (kg/m^3) and mass fractions. More...
 
compositionMap getMoleFractionsByName (double threshold=0.0) const
 Get the mole fractions by name. More...
 
double moleFraction (size_t k) const
 Return the mole fraction of a single species. More...
 
double moleFraction (const std::string &name) const
 Return the mole fraction of a single species. More...
 
compositionMap getMassFractionsByName (double threshold=0.0) const
 Get the mass fractions by name. More...
 
double massFraction (size_t k) const
 Return the mass fraction of a single species. More...
 
double massFraction (const std::string &name) const
 Return the mass fraction of a single species. More...
 
void getMoleFractions (double *const x) const
 Get the species mole fraction vector. More...
 
virtual void setMoleFractions (const double *const x)
 Set the mole fractions to the specified values. More...
 
virtual void setMoleFractions_NoNorm (const double *const x)
 Set the mole fractions to the specified values without normalizing. More...
 
void getMassFractions (double *const y) const
 Get the species mass fractions. More...
 
const double * massFractions () const
 Return a const pointer to the mass fraction array. More...
 
virtual void setMassFractions (const double *const y)
 Set the mass fractions to the specified values and normalize them. More...
 
virtual void setMassFractions_NoNorm (const double *const y)
 Set the mass fractions to the specified values without normalizing. More...
 
void getConcentrations (double *const c) const
 Get the species concentrations (kmol/m^3). More...
 
double concentration (const size_t k) const
 Concentration of species k. More...
 
virtual void setConcentrations (const double *const conc)
 Set the concentrations to the specified values within the phase. More...
 
virtual void setConcentrationsNoNorm (const double *const conc)
 Set the concentrations without ignoring negative concentrations. More...
 
doublereal temperature () const
 Temperature (K). More...
 
virtual double electronTemperature () const
 Electron Temperature (K) More...
 
virtual double pressure () const
 Return the thermodynamic pressure (Pa). More...
 
virtual double density () const
 Density (kg/m^3). More...
 
double molarDensity () const
 Molar density (kmol/m^3). More...
 
double molarVolume () const
 Molar volume (m^3/kmol). More...
 
virtual void setDensity (const double density_)
 Set the internally stored density (kg/m^3) of the phase. More...
 
virtual void setMolarDensity (const double molarDensity)
 Set the internally stored molar density (kmol/m^3) of the phase. More...
 
virtual void setElectronTemperature (double etemp)
 Set the internally stored electron temperature of the phase (K). More...
 
doublereal mean_X (const doublereal *const Q) const
 Evaluate the mole-fraction-weighted mean of an array Q. More...
 
doublereal mean_X (const vector_fp &Q) const
 Evaluate the mole-fraction-weighted mean of an array Q. More...
 
doublereal meanMolecularWeight () const
 The mean molecular weight. Units: (kg/kmol) More...
 
doublereal sum_xlogx () const
 Evaluate \( \sum_k X_k \log X_k \). More...
 
size_t addElement (const std::string &symbol, doublereal weight=-12345.0, int atomicNumber=0, doublereal entropy298=ENTROPY298_UNKNOWN, int elem_type=CT_ELEM_TYPE_ABSPOS)
 Add an element. More...
 
void addSpeciesAlias (const std::string &name, const std::string &alias)
 Add a species alias (that is, a user-defined alternative species name). More...
 
virtual std::vector< std::string > findIsomers (const compositionMap &compMap) const
 Return a vector with isomers names matching a given composition map. More...
 
virtual std::vector< std::string > findIsomers (const std::string &comp) const
 Return a vector with isomers names matching a given composition string. More...
 
shared_ptr< Speciesspecies (const std::string &name) const
 Return the Species object for the named species. More...
 
shared_ptr< Speciesspecies (size_t k) const
 Return the Species object for species whose index is k. More...
 
void ignoreUndefinedElements ()
 Set behavior when adding a species containing undefined elements to just skip the species. More...
 
void addUndefinedElements ()
 Set behavior when adding a species containing undefined elements to add those elements to the phase. More...
 
void throwUndefinedElements ()
 Set the behavior when adding a species containing undefined elements to throw an exception. More...
 

Protected Member Functions

virtual void compositionChanged ()
 Apply changes to the state which are needed after the composition changes. More...
 
virtual void _updateReferenceStateThermo () const
 Updates the reference state thermodynamic functions at the current T of the solution. More...
 
Critical State Properties.
virtual double critTemperature () const
 Critical temperature (K). More...
 
virtual double critPressure () const
 Critical pressure (Pa). More...
 
virtual double critVolume () const
 Critical volume (m3/kmol). More...
 
virtual double critCompressibility () const
 Critical compressibility (unitless). More...
 
virtual double critDensity () const
 Critical density (kg/m3). More...
 
virtual void calcCriticalConditions (double &pc, double &tc, double &vc) const
 
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. More...
 
- 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 newPhase(AnyMap&) function. More...
 
virtual void getCsvReportData (std::vector< std::string > &names, std::vector< vector_fp > &data) const
 Fills names and data with the column names and species thermo properties to be included in the output of the reportCSV method. More...
 
- Protected Member Functions inherited from Phase
void assertCompressible (const std::string &setter) const
 Ensure that phase is compressible. More...
 
void assignDensity (const double density_)
 Set the internally stored constant density (kg/m^3) of the phase. More...
 
void setMolecularWeight (const int k, const double mw)
 Set the molecular weight of a single species to a given value. More...
 

Protected Attributes

vector_fp m_tmpV
 Temporary storage - length = m_kk. More...
 
vector_fp moleFractions_
 Storage for the current values of the mole fractions of the species. More...
 
int iState_
 Current state of the fluid. More...
 
int forcedState_
 Force the system to be on a particular side of the spinodal curve. More...
 
vector_fp m_h0_RT
 Temporary storage for dimensionless reference state enthalpies. More...
 
vector_fp m_cp0_R
 Temporary storage for dimensionless reference state heat capacities. More...
 
vector_fp m_g0_RT
 Temporary storage for dimensionless reference state Gibbs energies. More...
 
vector_fp m_s0_R
 Temporary storage for dimensionless reference state entropies. More...
 
- Protected Attributes inherited from ThermoPhase
MultiSpeciesThermo m_spthermo
 Pointer to the calculation manager for species reference-state thermodynamic properties. More...
 
AnyMap m_input
 Data supplied via setParameters. More...
 
std::vector< const XML_Node * > m_speciesData
 Vector of pointers to the species databases. More...
 
doublereal m_phi
 Stored value of the electric potential for this phase. Units are Volts. More...
 
bool m_chargeNeutralityNecessary
 Boolean indicating whether a charge neutrality condition is a necessity. More...
 
int m_ssConvention
 Contains the standard state convention. More...
 
doublereal m_tlast
 last value of the temperature processed by reference state More...
 
- Protected Attributes inherited from Phase
ValueCache m_cache
 Cached for saved calculations within each ThermoPhase. More...
 
size_t m_kk
 Number of species in the phase. More...
 
size_t m_ndim
 Dimensionality of the phase. More...
 
vector_fp m_speciesComp
 Atomic composition of the species. More...
 
vector_fp m_speciesCharge
 Vector of species charges. length m_kk. More...
 
std::map< std::string, shared_ptr< Species > > m_species
 
UndefElement::behavior m_undefinedElementBehavior
 Flag determining behavior when adding species with an undefined element. More...
 
bool m_caseSensitiveSpecies
 Flag determining whether case sensitive species names are enforced. More...
 

Thermodynamic Values for the Species Reference States

virtual void getEnthalpy_RT_ref (doublereal *hrt) const
 
virtual void getGibbs_RT_ref (doublereal *grt) 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. More...
 
virtual void getGibbs_ref (doublereal *g) const
 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. More...
 
virtual void getEntropy_R_ref (doublereal *er) const
 Returns the vector of nondimensional entropies of the reference state at the current temperature of the solution and the reference pressure for each species. More...
 
virtual void getCp_R_ref (doublereal *cprt) const
 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. More...
 
virtual void getStandardVolumes_ref (doublereal *vol) const
 Get the molar volumes of the species reference states at the current T and P_ref of the solution. More...
 
const vector_fpgibbs_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. More...
 

Special Functions for fugacity classes

virtual doublereal liquidVolEst (doublereal TKelvin, doublereal &pres) const
 Estimate for the molar volume of the liquid. More...
 
virtual doublereal densityCalc (doublereal TKelvin, doublereal pressure, int phaseRequested, doublereal rhoguess)
 Calculates the density given the temperature and the pressure and a guess at the density. More...
 
int phaseState (bool checkState=false) const
 Returns the Phase State flag for the current state of the object. More...
 
virtual doublereal densSpinodalLiquid () const
 Return the value of the density at the liquid spinodal point (on the liquid side) for the current temperature. More...
 
virtual doublereal densSpinodalGas () const
 Return the value of the density at the gas spinodal point (on the gas side) for the current temperature. More...
 
doublereal calculatePsat (doublereal TKelvin, doublereal &molarVolGas, doublereal &molarVolLiquid)
 Calculate the saturation pressure at the current mixture content for the given temperature. More...
 
virtual doublereal satPressure (doublereal TKelvin)
 Calculate the saturation pressure at the current mixture content for the given temperature. More...
 
virtual void getActivityConcentrations (double *c) const
 This method returns an array of generalized concentrations. More...
 
doublereal z () const
 Calculate the value of z. More...
 
virtual doublereal sresid () const
 Calculate the deviation terms for the total entropy of the mixture from the ideal gas mixture. More...
 
virtual doublereal hresid () const
 Calculate the deviation terms for the total enthalpy of the mixture from the ideal gas mixture. More...
 
virtual doublereal psatEst (doublereal TKelvin) const
 Estimate for the saturation pressure. More...
 
int corr0 (doublereal TKelvin, doublereal pres, doublereal &densLiq, doublereal &densGas, doublereal &liqGRT, doublereal &gasGRT)
 Utility routine in the calculation of the saturation pressure. More...
 
virtual doublereal dpdVCalc (doublereal TKelvin, doublereal molarVol, doublereal &presCalc) const
 Calculate the pressure and the pressure derivative given the temperature and the molar volume. More...
 
virtual void updateMixingExpressions ()
 

Detailed Description

This is a filter class for ThermoPhase that implements some preparatory steps for efficiently handling mixture of gases that whose standard states are defined as ideal gases, but which describe also non-ideal solutions.

In addition a multicomponent liquid phase below the critical temperature of the mixture is also allowed. The main subclass is currently a mixture Redlich- Kwong class.

Attention
This class currently does not have any test cases or examples. Its implementation may be incomplete, and future changes to Cantera may unexpectedly cause this class to stop working. If you use this class, please consider contributing examples or test cases. In the absence of new tests or examples, this class may be deprecated and removed in a future version of Cantera. See https://github.com/Cantera/cantera/issues/267 for additional information.

Several concepts are introduced. The first concept is there are temporary variables for holding the species standard state values of Cp, H, S, G, and V at the last temperature and pressure called. These functions are not recalculated if a new call is made using the previous temperature and pressure.

The other concept is that the current state of the mixture is tracked. The state variable is either GAS, LIQUID, or SUPERCRIT fluid. Additionally, the variable LiquidContent is used and may vary between 0 and 1.

Typically, only one liquid phase is allowed to be formed within these classes. Additionally, there is an inherent contradiction between three phase models and the ThermoPhase class. The ThermoPhase class is really only meant to represent a single instantiation of a phase. The three phase models may be in equilibrium with multiple phases of the fluid in equilibrium with each other. This has yet to be resolved.

This class is usually used for non-ideal gases.

Definition at line 74 of file MixtureFugacityTP.h.

Constructor & Destructor Documentation

◆ MixtureFugacityTP()

Constructor.

Definition at line 20 of file MixtureFugacityTP.cpp.

Member Function Documentation

◆ type()

virtual std::string type ( ) const
inlinevirtual

String indicating the thermodynamic model implemented.

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

Reimplemented from ThermoPhase.

Reimplemented in PengRobinson, and RedlichKwongMFTP.

Definition at line 87 of file MixtureFugacityTP.h.

◆ standardStateConvention()

int standardStateConvention ( ) const
virtual

This method returns the convention used in specification of the standard state, of which there are currently two, temperature based, and variable pressure based.

Currently, there are two standard state conventions:

  • Temperature-based activities cSS_CONVENTION_TEMPERATURE 0
    • default
  • Variable Pressure and Temperature -based activities cSS_CONVENTION_VPSS 1
  • Thermodynamics is set via slave ThermoPhase objects with nothing being carried out at this ThermoPhase object level cSS_CONVENTION_SLAVE 2

Reimplemented from ThermoPhase.

Definition at line 26 of file MixtureFugacityTP.cpp.

References Cantera::cSS_CONVENTION_TEMPERATURE.

◆ setForcedSolutionBranch()

void setForcedSolutionBranch ( int  solnBranch)
virtual

Set the solution branch to force the ThermoPhase to exist on one branch or another.

Parameters
solnBranchBranch that the solution is restricted to. the value -1 means gas. The value -2 means unrestricted. Values of zero or greater refer to species dominated condensed phases.

Definition at line 31 of file MixtureFugacityTP.cpp.

References MixtureFugacityTP::forcedState_.

◆ forcedSolutionBranch()

int forcedSolutionBranch ( ) const
virtual

Report the solution branch which the solution is restricted to.

Returns
Branch that the solution is restricted to. the value -1 means gas. The value -2 means unrestricted. Values of zero or greater refer to species dominated condensed phases.

Definition at line 36 of file MixtureFugacityTP.cpp.

References MixtureFugacityTP::forcedState_.

◆ reportSolnBranchActual()

int reportSolnBranchActual ( ) const
virtual

Report the solution branch which the solution is actually on.

Returns
Branch that the solution is restricted to. the value -1 means gas. The value -2 means superfluid.. Values of zero or greater refer to species dominated condensed phases.

Definition at line 41 of file MixtureFugacityTP.cpp.

References MixtureFugacityTP::iState_.

◆ getdlnActCoeffdlnN_diag()

virtual void getdlnActCoeffdlnN_diag ( doublereal *  dlnActCoeffdlnN_diag) const
inlinevirtual

Get the array of log species mole number derivatives of the log activity coefficients.

For ideal mixtures (unity activity coefficients), this can return zero. Implementations should take the derivative of the logarithm of the activity coefficient with respect to the logarithm of the concentration- like variable (for example, moles) that represents the standard state. This quantity is to be used in conjunction with derivatives of that species mole number variable when the derivative of the chemical potential is taken.

units = dimensionless

Parameters
dlnActCoeffdlnN_diagOutput vector of derivatives of the log Activity Coefficients. length = m_kk

Reimplemented from ThermoPhase.

Definition at line 118 of file MixtureFugacityTP.h.

◆ enthalpy_mole()

double enthalpy_mole ( ) const
virtual

Molar enthalpy. Units: J/kmol.

Reimplemented from ThermoPhase.

Definition at line 47 of file MixtureFugacityTP.cpp.

References MixtureFugacityTP::hresid(), MixtureFugacityTP::m_h0_RT, Phase::mean_X(), and ThermoPhase::RT().

◆ entropy_mole()

double entropy_mole ( ) const
virtual

◆ getChemPotentials_RT()

void getChemPotentials_RT ( doublereal *  mu) const
virtual

Get the array of non-dimensional species chemical potentials These are partial molar Gibbs free energies.

\( \mu_k / \hat R T \). Units: unitless

We close the loop on this function, here, calling getChemPotentials() and then dividing by RT. No need for child classes to handle.

Parameters
muOutput vector of non-dimensional species chemical potentials Length: m_kk.

Reimplemented from ThermoPhase.

Reimplemented in RedlichKwongMFTP.

Definition at line 65 of file MixtureFugacityTP.cpp.

References ThermoPhase::getChemPotentials(), Phase::m_kk, and ThermoPhase::RT().

◆ getStandardChemPotentials()

void getStandardChemPotentials ( doublereal *  mu) const
virtual

Get the array of chemical potentials at unit activity.

These are the standard state chemical potentials \( \mu^0_k(T,P) \). The values are evaluated at the current temperature and pressure.

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.

Parameters
muOutput vector of standard state chemical potentials. length = m_kk. units are J / kmol.

Reimplemented from ThermoPhase.

Definition at line 75 of file MixtureFugacityTP.cpp.

References MixtureFugacityTP::m_g0_RT, Phase::m_kk, Phase::pressure(), ThermoPhase::refPressure(), and ThermoPhase::RT().

◆ getEnthalpy_RT()

void getEnthalpy_RT ( doublereal *  hrt) const
virtual

Get the nondimensional Enthalpy functions for the species at their standard states at the current T and P of the solution.

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.

Parameters
hrtOutput vector of standard state enthalpies. length = m_kk. units are unitless.

Reimplemented from ThermoPhase.

Definition at line 84 of file MixtureFugacityTP.cpp.

References MixtureFugacityTP::getEnthalpy_RT_ref().

◆ getEntropy_R()

void getEntropy_R ( doublereal *  sr) const
virtual

Get the array of nondimensional Enthalpy functions for the standard state species at the current T and P of the solution.

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.

Parameters
srOutput vector of nondimensional standard state entropies. length = m_kk.

Reimplemented from ThermoPhase.

Definition at line 89 of file MixtureFugacityTP.cpp.

References Phase::m_kk, MixtureFugacityTP::m_s0_R, Phase::pressure(), and ThermoPhase::refPressure().

◆ getGibbs_RT()

void getGibbs_RT ( doublereal *  grt) const
virtual

Get the nondimensional Gibbs functions for the species at their standard states of solution at the current T and P of the solution.

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.

Parameters
grtOutput vector of nondimensional standard state Gibbs free energies. length = m_kk.

Reimplemented from ThermoPhase.

Definition at line 98 of file MixtureFugacityTP.cpp.

References MixtureFugacityTP::m_g0_RT, Phase::m_kk, Phase::pressure(), and ThermoPhase::refPressure().

◆ getPureGibbs()

void getPureGibbs ( doublereal *  gpure) const
virtual

Get the pure Gibbs free energies of each species.

Species are assumed to be in their standard states.

This is the same as getStandardChemPotentials().

Parameters
[out]gpureArray of standard state Gibbs free energies. length = m_kk. units are J/kmol.

Reimplemented from ThermoPhase.

Definition at line 107 of file MixtureFugacityTP.cpp.

References MixtureFugacityTP::m_g0_RT, Phase::m_kk, Phase::pressure(), ThermoPhase::refPressure(), ThermoPhase::RT(), and Cantera::scale().

◆ getIntEnergy_RT()

void getIntEnergy_RT ( doublereal *  urt) const
virtual

Returns the vector of nondimensional internal Energies of the standard state at the current temperature and pressure of the solution for each species.

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.

\[ u^{ss}_k(T,P) = h^{ss}_k(T) - P * V^{ss}_k \]

Parameters
urtOutput vector of nondimensional standard state internal energies. length = m_kk.

Reimplemented from ThermoPhase.

Definition at line 116 of file MixtureFugacityTP.cpp.

References MixtureFugacityTP::m_h0_RT, and Phase::m_kk.

◆ getCp_R()

void getCp_R ( doublereal *  cpr) const
virtual

Get the nondimensional Heat Capacities at constant pressure for the standard state of the species at the current T and P.

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.

Parameters
cprOutput vector containing the the nondimensional Heat Capacities at constant pressure for the standard state of the species. Length: m_kk.

Reimplemented from ThermoPhase.

Definition at line 124 of file MixtureFugacityTP.cpp.

References MixtureFugacityTP::m_cp0_R.

◆ getStandardVolumes()

void getStandardVolumes ( doublereal *  vol) const
virtual

Get the molar volumes of each species in their standard states at the current T and P of the solution.

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.

units = m^3 / kmol

Parameters
volOutput vector of species volumes. length = m_kk. units = m^3 / kmol

Reimplemented from ThermoPhase.

Definition at line 129 of file MixtureFugacityTP.cpp.

References Phase::m_kk, Phase::pressure(), and ThermoPhase::RT().

Referenced by PengRobinson::standardConcentration(), and RedlichKwongMFTP::standardConcentration().

◆ setTemperature()

void setTemperature ( const doublereal  temp)
virtual

Set the temperature of the phase.

Currently this passes down to setState_TP(). It does not make sense to calculate the standard state without first setting T and P.

Parameters
tempTemperature (kelvin)

Reimplemented from Phase.

Definition at line 227 of file MixtureFugacityTP.cpp.

References MixtureFugacityTP::_updateReferenceStateThermo(), MixtureFugacityTP::iState_, MixtureFugacityTP::phaseState(), and Phase::setTemperature().

Referenced by MixtureFugacityTP::calculatePsat(), and RedlichKwongMFTP::densityCalc().

◆ setPressure()

void setPressure ( doublereal  p)
virtual

Set the internally stored pressure (Pa) at constant temperature and composition.

Currently this passes down to setState_TP(). It does not make sense to calculate the standard state without first setting T and P.

Parameters
pinput Pressure (Pa)

Reimplemented from Phase.

Definition at line 236 of file MixtureFugacityTP.cpp.

References Phase::density(), MixtureFugacityTP::forcedState_, and Phase::temperature().

◆ compositionChanged()

void compositionChanged ( )
protectedvirtual

Apply changes to the state which are needed after the composition changes.

This function is called after any call to setMassFractions(), setMoleFractions(), or similar. For phases which need to execute a callback after any change to the composition, it should be done by overriding this function rather than overriding all of the composition- setting functions. Derived class implementations of compositionChanged() should call the parent class method as well.

Reimplemented from Phase.

Definition at line 300 of file MixtureFugacityTP.cpp.

References Phase::compositionChanged(), Phase::getMoleFractions(), and MixtureFugacityTP::moleFractions_.

◆ _updateReferenceStateThermo()

void _updateReferenceStateThermo ( ) const
protectedvirtual

Updates the reference state thermodynamic functions at the current T of the solution.

This function must be called for every call to functions in this class. It checks to see whether the temperature has changed and thus the ss thermodynamics functions for all of the species must be recalculated.

This function is responsible for updating the following internal members:

  • m_h0_RT;
  • m_cp0_R;
  • m_g0_RT;
  • m_s0_R;

Definition at line 770 of file MixtureFugacityTP.cpp.

References MixtureFugacityTP::m_cp0_R, MixtureFugacityTP::m_g0_RT, MixtureFugacityTP::m_h0_RT, Phase::m_kk, MixtureFugacityTP::m_s0_R, ThermoPhase::m_spthermo, ThermoPhase::m_tlast, ThermoPhase::refPressure(), Phase::temperature(), and MultiSpeciesThermo::update().

Referenced by PengRobinson::cp_mole(), RedlichKwongMFTP::cp_mole(), PengRobinson::cv_mole(), RedlichKwongMFTP::cv_mole(), PengRobinson::pressure(), RedlichKwongMFTP::pressure(), and MixtureFugacityTP::setTemperature().

◆ getEnthalpy_RT_ref()

void getEnthalpy_RT_ref ( doublereal *  hrt) const
virtual

There are also temporary variables for holding the species reference- state values of Cp, H, S, and V at the last temperature and reference pressure called. These functions are not recalculated if a new call is made using the previous temperature. All calculations are done within the routine _updateRefStateThermo().

Reimplemented from ThermoPhase.

Definition at line 138 of file MixtureFugacityTP.cpp.

References MixtureFugacityTP::m_h0_RT.

Referenced by MixtureFugacityTP::getEnthalpy_RT(), PengRobinson::getPartialMolarEnthalpies(), and RedlichKwongMFTP::getPartialMolarEnthalpies().

◆ getGibbs_RT_ref()

void getGibbs_RT_ref ( doublereal *  grt) const
virtual

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.

Parameters
grtOutput vector containing the nondimensional reference state Gibbs Free energies. Length: m_kk.

Reimplemented from ThermoPhase.

Definition at line 143 of file MixtureFugacityTP.cpp.

References MixtureFugacityTP::m_g0_RT.

◆ gibbs_RT_ref()

const vector_fp & gibbs_RT_ref ( ) const
protected

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.

Returns
Output vector contains the nondimensional Gibbs free energies of the reference state of the species length = m_kk, units = dimensionless.

Definition at line 154 of file MixtureFugacityTP.cpp.

References MixtureFugacityTP::m_g0_RT.

Referenced by MixtureFugacityTP::getGibbs_ref().

◆ getGibbs_ref()

void getGibbs_ref ( doublereal *  g) const
virtual

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.

Parameters
gOutput vector containing the reference state Gibbs Free energies. Length: m_kk. Units: J/kmol.

Reimplemented from ThermoPhase.

Definition at line 148 of file MixtureFugacityTP.cpp.

References MixtureFugacityTP::gibbs_RT_ref(), ThermoPhase::RT(), and Cantera::scale().

Referenced by PengRobinson::getChemPotentials(), and RedlichKwongMFTP::getChemPotentials().

◆ getEntropy_R_ref()

void getEntropy_R_ref ( doublereal *  er) const
virtual

Returns the vector of nondimensional entropies of the reference state at the current temperature of the solution and the reference pressure for each species.

Parameters
erOutput vector containing the nondimensional reference state entropies. Length: m_kk.

Reimplemented from ThermoPhase.

Definition at line 159 of file MixtureFugacityTP.cpp.

References MixtureFugacityTP::m_s0_R.

Referenced by RedlichKwongMFTP::getPartialMolarEntropies().

◆ getCp_R_ref()

void getCp_R_ref ( doublereal *  cprt) const
virtual

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.

Parameters
cprtOutput vector of nondimensional reference state heat capacities at constant pressure for the species. Length: m_kk

Reimplemented from ThermoPhase.

Definition at line 164 of file MixtureFugacityTP.cpp.

References MixtureFugacityTP::m_cp0_R.

◆ getStandardVolumes_ref()

void getStandardVolumes_ref ( doublereal *  vol) const
virtual

Get the molar volumes of the species reference states at the current T and P_ref of the solution.

units = m^3 / kmol

Parameters
volOutput vector containing the standard state volumes. Length: m_kk.

Reimplemented from ThermoPhase.

Definition at line 169 of file MixtureFugacityTP.cpp.

References Phase::m_kk, ThermoPhase::refPressure(), and ThermoPhase::RT().

◆ addSpecies()

bool addSpecies ( shared_ptr< Species spec)
virtual

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().

Reimplemented from ThermoPhase.

Reimplemented in PengRobinson, and RedlichKwongMFTP.

Definition at line 209 of file MixtureFugacityTP.cpp.

References ThermoPhase::addSpecies(), MixtureFugacityTP::m_cp0_R, MixtureFugacityTP::m_g0_RT, MixtureFugacityTP::m_h0_RT, Phase::m_kk, MixtureFugacityTP::m_s0_R, MixtureFugacityTP::m_tmpV, and MixtureFugacityTP::moleFractions_.

Referenced by PengRobinson::addSpecies(), and RedlichKwongMFTP::addSpecies().

◆ setStateFromXML()

void setStateFromXML ( const XML_Node state)
virtual

Set the initial state of the phase to the conditions specified in the state XML element.

This method sets the temperature, pressure, and mole fraction vector to a set default value.

Parameters
stateAN XML_Node object corresponding to the "state" entry for this phase in the input file.
Deprecated:
The XML input format is deprecated and will be removed in Cantera 3.0.

Reimplemented from ThermoPhase.

Definition at line 176 of file MixtureFugacityTP.cpp.

References Phase::density(), Cantera::getChildValue(), Cantera::getFloat(), XML_Node::hasChild(), Phase::setMassFractionsByName(), Phase::setMoleFractionsByName(), ThermoPhase::setState_TP(), Phase::setState_TR(), and Phase::temperature().

◆ z()

doublereal z ( ) const
protected

Calculate the value of z.

\[ z = \frac{P v}{R T} \]

returns the value of z

Definition at line 316 of file MixtureFugacityTP.cpp.

References Phase::density(), Phase::meanMolecularWeight(), Phase::pressure(), and ThermoPhase::RT().

Referenced by PengRobinson::hresid(), RedlichKwongMFTP::hresid(), PengRobinson::sresid(), and RedlichKwongMFTP::sresid().

◆ sresid()

doublereal sresid ( ) const
protectedvirtual

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

Reimplemented in PengRobinson, and RedlichKwongMFTP.

Definition at line 321 of file MixtureFugacityTP.cpp.

Referenced by MixtureFugacityTP::entropy_mole().

◆ hresid()

doublereal hresid ( ) const
protectedvirtual

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

Reimplemented in PengRobinson, and RedlichKwongMFTP.

Definition at line 326 of file MixtureFugacityTP.cpp.

Referenced by MixtureFugacityTP::enthalpy_mole().

◆ psatEst()

doublereal psatEst ( doublereal  TKelvin) const
protectedvirtual

Estimate for the saturation pressure.

Note: this is only used as a starting guess for later routines that actually calculate an accurate value for the saturation pressure.

Parameters
TKelvintemperature in kelvin
Returns
the estimated saturation pressure at the given temperature

Definition at line 331 of file MixtureFugacityTP.cpp.

References MixtureFugacityTP::critPressure(), and MixtureFugacityTP::critTemperature().

Referenced by MixtureFugacityTP::calculatePsat(), and RedlichKwongMFTP::liquidVolEst().

◆ liquidVolEst()

doublereal liquidVolEst ( doublereal  TKelvin,
doublereal &  pres 
) const
virtual

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 in PengRobinson, and RedlichKwongMFTP.

Definition at line 342 of file MixtureFugacityTP.cpp.

Referenced by MixtureFugacityTP::calculatePsat().

◆ densityCalc()

doublereal densityCalc ( doublereal  TKelvin,
doublereal  pressure,
int  phaseRequested,
doublereal  rhoguess 
)
virtual

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 in PengRobinson, and RedlichKwongMFTP.

Definition at line 347 of file MixtureFugacityTP.cpp.

References MixtureFugacityTP::critTemperature(), Cantera::GasConstant, and Phase::meanMolecularWeight().

Referenced by MixtureFugacityTP::calculatePsat(), and MixtureFugacityTP::corr0().

◆ corr0()

int corr0 ( doublereal  TKelvin,
doublereal  pres,
doublereal &  densLiq,
doublereal &  densGas,
doublereal &  liqGRT,
doublereal &  gasGRT 
)
protected

Utility routine in the calculation of the saturation pressure.

Parameters
TKelvintemperature (kelvin)
prespressure (Pascal)
[out]densLiqdensity of liquid
[out]densGasdensity of gas
[out]liqGRTdeltaG/RT of liquid
[out]gasGRTdeltaG/RT of gas

Definition at line 488 of file MixtureFugacityTP.cpp.

References MixtureFugacityTP::densityCalc().

◆ phaseState()

int phaseState ( bool  checkState = false) const

Returns the Phase State flag for the current state of the object.

Parameters
checkStateIf true, this function does a complete check to see where in parameters space we are

There are three values:

  • WATER_GAS below the critical temperature but below the critical density
  • WATER_LIQUID below the critical temperature but above the critical density
  • WATER_SUPERCRIT above the critical temperature

Definition at line 517 of file MixtureFugacityTP.cpp.

References MixtureFugacityTP::critDensity(), MixtureFugacityTP::critTemperature(), MixtureFugacityTP::iState_, and Phase::temperature().

Referenced by MixtureFugacityTP::setTemperature().

◆ densSpinodalLiquid()

doublereal densSpinodalLiquid ( ) const
virtual

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 in PengRobinson, and RedlichKwongMFTP.

Definition at line 558 of file MixtureFugacityTP.cpp.

◆ densSpinodalGas()

doublereal densSpinodalGas ( ) const
virtual

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 in PengRobinson, and RedlichKwongMFTP.

Definition at line 563 of file MixtureFugacityTP.cpp.

◆ calculatePsat()

doublereal calculatePsat ( doublereal  TKelvin,
doublereal &  molarVolGas,
doublereal &  molarVolLiquid 
)

Calculate the saturation pressure at the current mixture content for the given temperature.

Parameters
TKelvin(input) Temperature (Kelvin)
molarVolGas(return) Molar volume of the gas
molarVolLiquid(return) Molar volume of the liquid
Returns
the saturation pressure at the given temperature

Definition at line 575 of file MixtureFugacityTP.cpp.

References MixtureFugacityTP::critTemperature(), Phase::density(), MixtureFugacityTP::densityCalc(), Cantera::GasConstant, MixtureFugacityTP::liquidVolEst(), Phase::meanMolecularWeight(), MixtureFugacityTP::psatEst(), MixtureFugacityTP::setTemperature(), and Phase::temperature().

Referenced by MixtureFugacityTP::satPressure().

◆ satPressure()

doublereal satPressure ( doublereal  TKelvin)
virtual

Calculate the saturation pressure at the current mixture content for the given temperature.

Parameters
TKelvinTemperature (Kelvin)
Returns
The saturation pressure at the given temperature

Reimplemented from ThermoPhase.

Definition at line 568 of file MixtureFugacityTP.cpp.

References MixtureFugacityTP::calculatePsat().

◆ getActivityConcentrations()

void getActivityConcentrations ( double *  c) const
virtual

This method returns an array of generalized concentrations.

\( C^a_k\) are defined such that \( a_k = C^a_k / C^0_k, \) where \( C^0_k \) is a standard concentration defined below and \( a_k \) are activities used in the thermodynamic functions. These activity (or generalized) concentrations are used by kinetics manager classes to compute the forward and reverse rates of elementary reactions. Note that they may or may not have units of concentration — they might be partial pressures, mole fractions, or surface coverages, for example.

Parameters
cOutput array of generalized concentrations. The units depend upon the implementation of the reaction rate expressions within the phase.

Reimplemented from ThermoPhase.

Definition at line 307 of file MixtureFugacityTP.cpp.

References ThermoPhase::getActivityCoefficients(), Phase::m_kk, Phase::moleFraction(), Phase::pressure(), and ThermoPhase::RT().

◆ dpdVCalc()

doublereal dpdVCalc ( doublereal  TKelvin,
doublereal  molarVol,
doublereal &  presCalc 
) const
protectedvirtual

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 in PengRobinson, and RedlichKwongMFTP.

Definition at line 765 of file MixtureFugacityTP.cpp.

◆ updateMixingExpressions()

void updateMixingExpressions ( )
protectedvirtual

Reimplemented in PengRobinson, and RedlichKwongMFTP.

Definition at line 484 of file MixtureFugacityTP.cpp.

◆ critTemperature()

double critTemperature ( ) const
protectedvirtual

◆ critPressure()

double critPressure ( ) const
protectedvirtual

Critical pressure (Pa).

Reimplemented from ThermoPhase.

Definition at line 799 of file MixtureFugacityTP.cpp.

Referenced by RedlichKwongMFTP::liquidVolEst(), and MixtureFugacityTP::psatEst().

◆ critVolume()

double critVolume ( ) const
protectedvirtual

Critical volume (m3/kmol).

Reimplemented from ThermoPhase.

Definition at line 806 of file MixtureFugacityTP.cpp.

◆ critCompressibility()

double critCompressibility ( ) const
protectedvirtual

Critical compressibility (unitless).

Reimplemented from ThermoPhase.

Definition at line 813 of file MixtureFugacityTP.cpp.

References Cantera::GasConstant.

◆ critDensity()

double critDensity ( ) const
protectedvirtual

Critical density (kg/m3).

Reimplemented from ThermoPhase.

Definition at line 820 of file MixtureFugacityTP.cpp.

References Phase::meanMolecularWeight().

Referenced by MixtureFugacityTP::phaseState().

◆ calcCriticalConditions()

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

Definition at line 828 of file MixtureFugacityTP.cpp.

◆ solveCubic()

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
protected

Solve the cubic equation of state.

Returns the number of solutions found. For the gas phase solution, it returns a positive number (1 or 2). If it only finds the liquid branch solution, it will return -1 or -2 instead of 1 or 2. If it returns 0, then there is an error. The cubic equation is solved using Nickall's method (Ref: The Mathematical Gazette(1993), 77(November), 354–359, https://www.jstor.org/stable/3619777)

Parameters
Ttemperature (kelvin)
prespressure (Pa)
a"a" parameter in the non-ideal EoS [Pa-m^6/kmol^2]
b"b" parameter in the non-ideal EoS [m^3/kmol]
aAlphaa*alpha (temperature dependent function for P-R EoS, 1 for R-K EoS)
VrootRoots of the cubic equation for molar volume (m3/kmol)
anconstant used in cubic equation
bnconstant used in cubic equation
cnconstant used in cubic equation
dnconstant used in cubic equation
tcCritical temperature (kelvin)
vcCritical volume
Returns
the number of solutions found

Definition at line 833 of file MixtureFugacityTP.cpp.

Member Data Documentation

◆ m_tmpV

vector_fp m_tmpV
mutableprotected

◆ moleFractions_

vector_fp moleFractions_
protected

◆ iState_

int iState_
protected

Current state of the fluid.

There are three possible states of the fluid:

  • FLUID_GAS
  • FLUID_LIQUID
  • FLUID_SUPERCRIT

Definition at line 564 of file MixtureFugacityTP.h.

Referenced by MixtureFugacityTP::phaseState(), MixtureFugacityTP::reportSolnBranchActual(), and MixtureFugacityTP::setTemperature().

◆ forcedState_

int forcedState_
protected

Force the system to be on a particular side of the spinodal curve.

Definition at line 567 of file MixtureFugacityTP.h.

Referenced by MixtureFugacityTP::forcedSolutionBranch(), MixtureFugacityTP::setForcedSolutionBranch(), and MixtureFugacityTP::setPressure().

◆ m_h0_RT

vector_fp m_h0_RT
mutableprotected

◆ m_cp0_R

vector_fp m_cp0_R
mutableprotected

◆ m_g0_RT

vector_fp m_g0_RT
mutableprotected

◆ m_s0_R

vector_fp m_s0_R
mutableprotected

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