Cantera
2.5.1
|
RedlichKisterVPSSTP is a derived class of GibbsExcessVPSSTP that employs the Redlich-Kister approximation for the excess Gibbs free energy. More...
#include <RedlichKisterVPSSTP.h>
Public Member Functions | |
RedlichKisterVPSSTP () | |
Constructor. More... | |
RedlichKisterVPSSTP (const std::string &inputFile, const std::string &id="") | |
Construct a RedlichKisterVPSSTP object from an input file. More... | |
RedlichKisterVPSSTP (XML_Node &phaseRef, const std::string &id="") | |
Construct and initialize a RedlichKisterVPSSTP ThermoPhase object directly from an XML database. More... | |
virtual std::string | type () const |
String indicating the thermodynamic model implemented. More... | |
Molar Thermodynamic Properties | |
virtual doublereal | enthalpy_mole () const |
Molar enthalpy. Units: J/kmol. More... | |
virtual doublereal | entropy_mole () const |
Molar entropy. Units: J/kmol/K. 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... | |
Activities, Standard States, and Activity Concentrations | |
The activity \(a_k\) of a species in solution is related to the chemical potential by \[ \mu_k = \mu_k^0(T) + \hat R T \log a_k. \] The quantity \(\mu_k^0(T,P)\) is the chemical potential at unit activity, which depends only on temperature and pressure. | |
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... | |
Partial Molar Properties of the Solution | |
virtual void | getChemPotentials (doublereal *mu) const |
Get the species chemical potentials. Units: J/kmol. 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 for the species in the mixture. More... | |
virtual void | getPartialMolarCp (doublereal *cpbar) const |
Returns an array of partial molar entropies 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 | getd2lnActCoeffdT2 (doublereal *d2lnActCoeffdT2) const |
Get the array of temperature second derivatives of the log activity coefficients. More... | |
virtual void | getdlnActCoeffdT (doublereal *dlnActCoeffdT) const |
Get the array of temperature derivatives of the log activity coefficients. More... | |
Initialization | |
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(). | |
virtual void | initThermo () |
virtual void | initThermoXML (XML_Node &phaseNode, const std::string &id) |
Import and initialize a ThermoPhase object using an XML tree. More... | |
void | addBinaryInteraction (const std::string &speciesA, const std::string &speciesB, const double *excess_enthalpy, size_t n_enthalpy, const double *excess_entropy, size_t n_entropy) |
Add a binary species interaction with the specified parameters. More... | |
Public Member Functions inherited from GibbsExcessVPSSTP | |
GibbsExcessVPSSTP () | |
virtual Units | standardConcentrationUnits () const |
Returns the units of the "standard concentration" for this phase. More... | |
virtual void | getActivityConcentrations (doublereal *c) const |
This method returns an array of generalized concentrations. More... | |
virtual doublereal | standardConcentration (size_t k=0) const |
The standard concentration \( C^0_k \) used to normalize the generalized concentration. More... | |
virtual doublereal | logStandardConc (size_t k=0) const |
Natural logarithm of the standard concentration of the kth species. More... | |
virtual void | getActivities (doublereal *ac) const |
Get the array of non-dimensional activities (molality based for this class and classes that derive from it) 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 | getdlnActCoeffdlnX (doublereal *dlnActCoeffdlnX) const |
Get the array of log concentration-like derivatives of the log activity coefficients. More... | |
virtual const vector_fp & | getPartialMolarVolumesVector () const |
virtual bool | addSpecies (shared_ptr< Species > spec) |
Public Member Functions inherited from VPStandardStateTP | |
VPStandardStateTP () | |
Constructor. More... | |
virtual bool | isCompressible () const |
Return whether phase represents a compressible substance. 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 | getChemPotentials_RT (doublereal *mu) const |
Get the array of non-dimensional species chemical potentials. More... | |
void | installPDSS (size_t k, std::unique_ptr< PDSS > &&pdss) |
Install a PDSS object for species k More... | |
PDSS * | providePDSS (size_t k) |
const PDSS * | providePDSS (size_t k) const |
virtual bool | addSpecies (shared_ptr< Species > spec) |
Add a Species to this Phase. More... | |
virtual void | getStandardChemPotentials (doublereal *mu) const |
Get the array of chemical potentials at unit activity for the species at their standard states at the current T and P of the solution. 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 Entropy 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 in their standard states at the current T and P of the solution. More... | |
virtual void | getPureGibbs (doublereal *gpure) const |
Get the Gibbs functions for the standard state of the species at the current T and P of the solution. More... | |
virtual void | getIntEnergy_RT (doublereal *urt) const |
Returns the vector of nondimensional Internal Energies of the standard state species at the current T and P of the solution. More... | |
virtual void | getCp_R (doublereal *cpr) const |
Get the nondimensional Heat Capacities at constant pressure for the species standard states at the current T and P of the solution. More... | |
virtual void | getStandardVolumes (doublereal *vol) const |
Get the molar volumes of the species standard states at the current T and P of the solution. More... | |
virtual const vector_fp & | getStandardVolumes () const |
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... | |
virtual void | setState_TP (doublereal T, doublereal pres) |
Set the temperature and pressure at the same time. More... | |
virtual doublereal | pressure () const |
Returns the current pressure of the phase. More... | |
virtual void | updateStandardStateThermo () const |
Updates the standard state thermodynamic functions at the current T and P of the solution. More... | |
virtual double | minTemp (size_t k=npos) const |
Minimum temperature for which the thermodynamic data for the species or phase are valid. More... | |
virtual double | maxTemp (size_t k=npos) const |
Maximum temperature for which the thermodynamic data for the species are valid. More... | |
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... | |
Public Member Functions inherited from ThermoPhase | |
ThermoPhase () | |
Constructor. 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... | |
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... | |
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 | 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... | |
void | getElectrochemPotentials (doublereal *mu) const |
Get the species electrochemical potentials. More... | |
virtual void | getPartialMolarIntEnergies (doublereal *ubar) const |
Return an array of partial molar internal energies 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... | |
doublereal | RT () const |
Return the Gas Constant multiplied by the current temperature. 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_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 | equivalenceRatio () const |
Compute the equivalence ratio for the current mixture from available oxygen and required oxygen. 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 | critTemperature () const |
Critical temperature (K). More... | |
virtual doublereal | critPressure () const |
Critical pressure (Pa). More... | |
virtual doublereal | critVolume () const |
Critical volume (m3/kmol). More... | |
virtual doublereal | critCompressibility () const |
Critical compressibility (unitless). More... | |
virtual doublereal | critDensity () const |
Critical density (kg/m3). More... | |
virtual doublereal | satTemperature (doublereal p) const |
Return the saturation temperature given the pressure. More... | |
virtual doublereal | satPressure (doublereal t) |
Return the saturation pressure given the temperature. 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 MultiSpeciesThermo & | speciesThermo (int k=-1) |
Return a changeable reference to the calculation manager for species reference-state thermodynamic properties. More... | |
virtual const MultiSpeciesThermo & | speciesThermo (int k=-1) const |
virtual void | initThermoFile (const std::string &inputFile, const std::string &id) |
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... | |
const AnyMap & | input () const |
Access input data associated with the phase description. More... | |
AnyMap & | input () |
virtual void | setParametersFromXML (const XML_Node &eosdata) |
Set equation of state parameter values from XML entries. More... | |
virtual void | setStateFromXML (const XML_Node &state) |
Set the initial state of the phase to the conditions specified in the state XML element. 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... | |
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... | |
Public Member Functions inherited from Phase | |
Phase () | |
Default constructor. More... | |
Phase (const Phase &)=delete | |
Phase & | operator= (const Phase &)=delete |
XML_Node & | xml () 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 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_fp & | molecularWeights () 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... | |
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, e.g. 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 | id () const |
Return the string id for the phase. More... | |
void | setID (const std::string &id) |
Set the string id for the phase. 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_fp & | atomicWeights () 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 | 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... | |
doublereal | temperature () const |
Temperature (K). 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... | |
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 (i.e. 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< Species > | species (const std::string &name) const |
Return the Species object for the named species. More... | |
shared_ptr< Species > | species (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... | |
Derivatives of Thermodynamic Variables needed for Applications | |
size_t | numBinaryInteractions_ |
number of binary interaction expressions More... | |
std::vector< size_t > | m_pSpecies_A_ij |
vector of species indices representing species A in the interaction More... | |
std::vector< size_t > | m_pSpecies_B_ij |
vector of species indices representing species B in the interaction More... | |
std::vector< size_t > | m_N_ij |
Vector of the length of the polynomial for the interaction. More... | |
std::vector< vector_fp > | m_HE_m_ij |
Enthalpy term for the binary mole fraction interaction of the excess Gibbs free energy expression. More... | |
std::vector< vector_fp > | m_SE_m_ij |
Entropy term for the binary mole fraction interaction of the excess Gibbs free energy expression. More... | |
int | formRedlichKister_ |
form of the RedlichKister interaction expression. More... | |
int | formTempModel_ |
form of the temperature dependence of the Redlich-Kister interaction expression. More... | |
Array2D | dlnActCoeff_dX_ |
Two dimensional array of derivatives of activity coefficients wrt mole fractions. 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_diag (doublereal *dlnActCoeffdlnN_diag) const |
Get the array of log species mole number derivatives of the log activity coefficients. 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... | |
void | readXMLBinarySpecies (XML_Node &xmlBinarySpecies) |
Process an XML node called "binaryNeutralSpeciesParameters". More... | |
void | initLengths () |
Initialize lengths of local variables after all species have been identified. More... | |
void | s_update_lnActCoeff () const |
Update the activity coefficients. More... | |
void | s_update_dlnActCoeff_dT () const |
Update the derivative of the log of the activity coefficients wrt T. More... | |
void | s_update_dlnActCoeff_dX_ () const |
Internal routine that calculates the derivative of the activity coefficients wrt the mole fractions. More... | |
void | s_update_dlnActCoeff_dlnX_diag () const |
Internal routine that calculates the total derivative of the activity coefficients with respect to the log of the mole fractions. More... | |
Additional Inherited Members | |
Protected Member Functions inherited from GibbsExcessVPSSTP | |
void | calcDensity () |
Calculate the density of the mixture using the partial molar volumes and mole fractions as input. More... | |
virtual void | compositionChanged () |
Apply changes to the state which are needed after the composition changes. More... | |
double | checkMFSum (const doublereal *const x) const |
utility routine to check mole fraction sum More... | |
Protected Member Functions inherited from VPStandardStateTP | |
virtual void | invalidateCache () |
Invalidate any cached values which are normally updated only when a change in state is detected. More... | |
virtual void | _updateStandardStateThermo () const |
Updates the standard state thermodynamic functions at the current T and P of the solution. More... | |
const vector_fp & | Gibbs_RT_ref () const |
Protected Member Functions inherited from ThermoPhase | |
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 inherited from GibbsExcessVPSSTP | |
vector_fp | moleFractions_ |
Storage for the current values of the mole fractions of the species. More... | |
vector_fp | lnActCoeff_Scaled_ |
Storage for the current values of the activity coefficients of the species. More... | |
vector_fp | dlnActCoeffdT_Scaled_ |
Storage for the current derivative values of the gradients with respect to temperature of the log of the activity coefficients of the species. More... | |
vector_fp | d2lnActCoeffdT2_Scaled_ |
Storage for the current derivative values of the gradients with respect to temperature of the log of the activity coefficients of the species. More... | |
vector_fp | dlnActCoeffdlnN_diag_ |
Storage for the current derivative values of the gradients with respect to logarithm of the mole fraction of the log of the activity coefficients of the species. More... | |
vector_fp | dlnActCoeffdlnX_diag_ |
Storage for the current derivative values of the gradients with respect to logarithm of the mole fraction of the log of the activity coefficients of the species. More... | |
Array2D | dlnActCoeffdlnN_ |
Storage for the current derivative values of the gradients with respect to logarithm of the species mole number of the log of the activity coefficients of the species. More... | |
Protected Attributes inherited from VPStandardStateTP | |
doublereal | m_Pcurrent |
Current value of the pressure - state variable. More... | |
double | m_minTemp |
The minimum temperature at which data for all species is valid. More... | |
double | m_maxTemp |
The maximum temperature at which data for all species is valid. More... | |
doublereal | m_Tlast_ss |
The last temperature at which the standard state thermodynamic properties were calculated at. More... | |
doublereal | m_Plast_ss |
The last pressure at which the Standard State thermodynamic properties were calculated at. More... | |
std::vector< std::unique_ptr< PDSS > > | m_PDSS_storage |
Storage for the PDSS objects for the species. More... | |
vector_fp | m_h0_RT |
Vector containing the species reference enthalpies at T = m_tlast and P = p_ref. More... | |
vector_fp | m_cp0_R |
Vector containing the species reference constant pressure heat capacities at T = m_tlast and P = p_ref. More... | |
vector_fp | m_g0_RT |
Vector containing the species reference Gibbs functions at T = m_tlast and P = p_ref. More... | |
vector_fp | m_s0_R |
Vector containing the species reference entropies at T = m_tlast and P = p_ref. More... | |
vector_fp | m_V0 |
Vector containing the species reference molar volumes. More... | |
vector_fp | m_hss_RT |
Vector containing the species Standard State enthalpies at T = m_tlast and P = m_plast. More... | |
vector_fp | m_cpss_R |
Vector containing the species Standard State constant pressure heat capacities at T = m_tlast and P = m_plast. More... | |
vector_fp | m_gss_RT |
Vector containing the species Standard State Gibbs functions at T = m_tlast and P = m_plast. More... | |
vector_fp | m_sss_R |
Vector containing the species Standard State entropies at T = m_tlast and P = m_plast. More... | |
vector_fp | m_Vss |
Vector containing the species standard state volumes at T = m_tlast and P = m_plast. 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... | |
RedlichKisterVPSSTP is a derived class of GibbsExcessVPSSTP that employs the Redlich-Kister approximation for the excess Gibbs free energy.
RedlichKisterVPSSTP derives from class GibbsExcessVPSSTP which is derived from VPStandardStateTP, and overloads the virtual methods defined there with ones that use expressions appropriate for the Redlich Kister Excess Gibbs free energy approximation.
The independent unknowns are pressure, temperature, and mass fraction.
All species are defined to have standard states that depend upon both the temperature and the pressure. The Redlich-Kister approximation assumes symmetric standard states, where all of the standard state assume that the species are in pure component states at the temperature and pressure of the solution. I don't think it prevents, however, some species from being dilute in the solution.
The molar excess Gibbs free energy is given by the following formula which is a sum over interactions i. Each of the interactions are binary interactions involving two of the species in the phase, denoted, Ai and Bi. This is the generalization of the Redlich-Kister formulation for a phase that has more than 2 species.
\[ G^E = \sum_{i} G^E_{i} \]
where
\[ G^E_{i} = n X_{Ai} X_{Bi} \sum_m \left( A^{i}_m {\left( X_{Ai} - X_{Bi} \right)}^m \right) \]
and where we can break down the Gibbs free energy contributions into enthalpy and entropy contributions
\[ H^E_i = n X_{Ai} X_{Bi} \sum_m \left( H^{i}_m {\left( X_{Ai} - X_{Bi} \right)}^m \right) \]
\[ S^E_i = n X_{Ai} X_{Bi} \sum_m \left( S^{i}_m {\left( X_{Ai} - X_{Bi} \right)}^m \right) \]
where n is the total moles in the solution. The activity of a species defined in the phase is given by an excess Gibbs free energy formulation.
\[ a_k = \gamma_k X_k \]
where
\[ R T \ln( \gamma_k )= \frac{d(n G^E)}{d(n_k)}\Bigg|_{n_i} \]
Taking the derivatives results in the following expression
\[ R T \ln( \gamma_k )= \sum_i \delta_{Ai,k} (1 - X_{Ai}) X_{Bi} \sum_m \left( A^{i}_m {\left( X_{Ai} - X_{Bi} \right)}^m \right) + \sum_i \delta_{Ai,k} X_{Ai} X_{Bi} \sum_m \left( A^{i}_0 + A^{i}_m {\left( X_{Ai} - X_{Bi} \right)}^{m-1} (1 - X_{Ai} + X_{Bi}) \right) \]
This object inherits from the class VPStandardStateTP. Therefore, the specification and calculation of all standard state and reference state values are handled at that level. Various functional forms for the standard state are permissible. The chemical potential for species k is equal to
\[ \mu_k(T,P) = \mu^o_k(T, P) + R T \ln(\gamma_k X_k) \]
The partial molar entropy for species k is given by the following relation,
\[ \tilde{s}_k(T,P) = s^o_k(T,P) - R \ln( \gamma_k X_k ) - R T \frac{d \ln(\gamma_k) }{dT} \]
The partial molar enthalpy for species k is given by
\[ \tilde{h}_k(T,P) = h^o_k(T,P) - R T^2 \frac{d \ln(\gamma_k)}{dT} \]
The partial molar volume for species k is
\[ \tilde V_k(T,P) = V^o_k(T,P) + R T \frac{d \ln(\gamma_k) }{dP} \]
The partial molar Heat Capacity for species k is
\[ \tilde{C}_{p,k}(T,P) = C^o_{p,k}(T,P) - 2 R T \frac{d \ln( \gamma_k )}{dT} - R T^2 \frac{d^2 \ln(\gamma_k) }{{dT}^2} \]
\( C^a_k\) are defined such that \( a_k = C^a_k / C^s_k, \) where \( C^s_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. The activity concentration, \( C^a_k \),is given by the following expression.
\[ C^a_k = C^s_k X_k = \frac{P}{R T} X_k \]
The standard concentration for species k is independent of k and equal to
\[ C^s_k = C^s = \frac{P}{R T} \]
For example, a bulk-phase binary gas reaction between species j and k, producing a new gas species l would have the following equation for its rate of progress variable, \( R^1 \), which has units of kmol m-3 s-1.
\[ R^1 = k^1 C_j^a C_k^a = k^1 (C^s a_j) (C^s a_k) \]
where
\[ C_j^a = C^s a_j \mbox{\quad and \quad} C_k^a = C^s a_k \]
\( C_j^a \) is the activity concentration of species j, and \( C_k^a \) is the activity concentration of species k. \( C^s \) is the standard concentration. \( a_j \) is the activity of species j which is equal to the mole fraction of j.
The reverse rate constant can then be obtained from the law of microscopic reversibility and the equilibrium expression for the system.
\[ \frac{a_j a_k}{ a_l} = K_a^{o,1} = \exp(\frac{\mu^o_l - \mu^o_j - \mu^o_k}{R T} ) \]
\( K_a^{o,1} \) is the dimensionless form of the equilibrium constant, associated with the pressure dependent standard states \( \mu^o_l(T,P) \) and their associated activities, \( a_l \), repeated here:
\[ \mu_l(T,P) = \mu^o_l(T, P) + R T \log(a_l) \]
We can switch over to expressing the equilibrium constant in terms of the reference state chemical potentials
\[ K_a^{o,1} = \exp(\frac{\mu^{ref}_l - \mu^{ref}_j - \mu^{ref}_k}{R T} ) * \frac{P_{ref}}{P} \]
The concentration equilibrium constant, \( K_c \), may be obtained by changing over to activity concentrations. When this is done:
\[ \frac{C^a_j C^a_k}{ C^a_l} = C^o K_a^{o,1} = K_c^1 = \exp(\frac{\mu^{ref}_l - \mu^{ref}_j - \mu^{ref}_k}{R T} ) * \frac{P_{ref}}{RT} \]
Kinetics managers will calculate the concentration equilibrium constant, \( K_c \), using the second and third part of the above expression as a definition for the concentration equilibrium constant.
For completeness, the pressure equilibrium constant may be obtained as well
\[ \frac{P_j P_k}{ P_l P_{ref}} = K_p^1 = \exp(\frac{\mu^{ref}_l - \mu^{ref}_j - \mu^{ref}_k}{R T} ) \]
\( K_p \) is the simplest form of the equilibrium constant for ideal gases. However, it isn't necessarily the simplest form of the equilibrium constant for other types of phases; \( K_c \) is used instead because it is completely general.
The reverse rate of progress may be written down as
\[ R^{-1} = k^{-1} C_l^a = k^{-1} (C^o a_l) \]
where we can use the concept of microscopic reversibility to write the reverse rate constant in terms of the forward rate constant and the concentration equilibrium constant, \( K_c \).
\[ k^{-1} = k^1 K^1_c \]
\(k^{-1} \) has units of s-1.
Definition at line 216 of file RedlichKisterVPSSTP.h.
Constructor.
This doesn't do much more than initialize constants with default values.
Definition at line 21 of file RedlichKisterVPSSTP.cpp.
RedlichKisterVPSSTP | ( | const std::string & | inputFile, |
const std::string & | id = "" |
||
) |
Construct a RedlichKisterVPSSTP object from an input file.
inputFile | Name of the input file containing the phase definition |
id | name (ID) of the phase in the input file. If empty, the first phase definition in the input file will be used. |
Definition at line 28 of file RedlichKisterVPSSTP.cpp.
References ThermoPhase::initThermoFile().
RedlichKisterVPSSTP | ( | XML_Node & | phaseRef, |
const std::string & | id = "" |
||
) |
Construct and initialize a RedlichKisterVPSSTP ThermoPhase object directly from an XML database.
phaseRef | XML phase node containing the description of the phase |
id | id attribute containing the name of the phase. (default is the empty string) |
Definition at line 37 of file RedlichKisterVPSSTP.cpp.
References Cantera::importPhase().
|
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.
Definition at line 245 of file RedlichKisterVPSSTP.h.
|
virtual |
Molar enthalpy. Units: J/kmol.
Reimplemented from ThermoPhase.
Definition at line 74 of file RedlichKisterVPSSTP.cpp.
References RedlichKisterVPSSTP::getPartialMolarEnthalpies(), Phase::m_kk, and GibbsExcessVPSSTP::moleFractions_.
|
virtual |
Molar entropy. Units: J/kmol/K.
Reimplemented from ThermoPhase.
Definition at line 85 of file RedlichKisterVPSSTP.cpp.
References RedlichKisterVPSSTP::getPartialMolarEntropies(), Phase::m_kk, and GibbsExcessVPSSTP::moleFractions_.
|
virtual |
Molar heat capacity at constant pressure. Units: J/kmol/K.
Reimplemented from ThermoPhase.
Definition at line 96 of file RedlichKisterVPSSTP.cpp.
References RedlichKisterVPSSTP::getPartialMolarCp(), Phase::m_kk, and GibbsExcessVPSSTP::moleFractions_.
Referenced by RedlichKisterVPSSTP::cv_mole().
|
virtual |
Molar heat capacity at constant volume. Units: J/kmol/K.
Reimplemented from ThermoPhase.
Definition at line 107 of file RedlichKisterVPSSTP.cpp.
References RedlichKisterVPSSTP::cp_mole(), and Cantera::GasConstant.
|
virtual |
Get the array of non-dimensional molar-based ln activity coefficients at the current solution temperature, pressure, and solution concentration.
lnac | Output vector of ln activity coefficients. Length: m_kk. |
Reimplemented from ThermoPhase.
Definition at line 48 of file RedlichKisterVPSSTP.cpp.
References GibbsExcessVPSSTP::lnActCoeff_Scaled_, Phase::m_kk, and RedlichKisterVPSSTP::s_update_lnActCoeff().
|
virtual |
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.
mu | Output vector of species chemical potentials. Length: m_kk. Units: J/kmol |
Reimplemented from ThermoPhase.
Definition at line 60 of file RedlichKisterVPSSTP.cpp.
References VPStandardStateTP::getStandardChemPotentials(), GibbsExcessVPSSTP::lnActCoeff_Scaled_, Phase::m_kk, GibbsExcessVPSSTP::moleFractions_, ThermoPhase::RT(), RedlichKisterVPSSTP::s_update_lnActCoeff(), and Cantera::SmallNumber.
|
virtual |
Returns an array of partial molar enthalpies for the species in the mixture.
Units (J/kmol)
For this phase, the partial molar enthalpies are equal to the standard state enthalpies modified by the derivative of the molality-based activity coefficient wrt temperature
\[ \bar h_k(T,P) = h^o_k(T,P) - R T^2 \frac{d \ln(\gamma_k)}{dT} \]
hbar | Vector of returned partial molar enthalpies (length m_kk, units = J/kmol) |
Reimplemented from ThermoPhase.
Definition at line 112 of file RedlichKisterVPSSTP.cpp.
References VPStandardStateTP::getEnthalpy_RT().
Referenced by RedlichKisterVPSSTP::enthalpy_mole().
|
virtual |
Returns an array of partial molar entropies for the species in the mixture.
Units (J/kmol)
For this phase, the partial molar enthalpies are equal to the standard state enthalpies modified by the derivative of the activity coefficient wrt temperature
\[ \bar s_k(T,P) = s^o_k(T,P) - R T^2 \frac{d \ln(\gamma_k)}{dT} - R \ln( \gamma_k X_k) - R T \frac{d \ln(\gamma_k) }{dT} \]
sbar | Vector of returned partial molar entropies (length m_kk, units = J/kmol/K) |
Reimplemented from ThermoPhase.
Definition at line 150 of file RedlichKisterVPSSTP.cpp.
References VPStandardStateTP::getEntropy_R().
Referenced by RedlichKisterVPSSTP::entropy_mole().
|
virtual |
Returns an array of partial molar entropies for the species in the mixture.
Units (J/kmol)
For this phase, the partial molar enthalpies are equal to the standard state enthalpies modified by the derivative of the activity coefficient wrt temperature
\[ ??????????????? \bar s_k(T,P) = s^o_k(T,P) - R T^2 \frac{d \ln(\gamma_k)}{dT} - R \ln( \gamma_k X_k) - R T \frac{d \ln(\gamma_k) }{dT} ??????????????? \]
cpbar | Vector of returned partial molar heat capacities (length m_kk, units = J/kmol/K) |
Reimplemented from ThermoPhase.
Definition at line 131 of file RedlichKisterVPSSTP.cpp.
References VPStandardStateTP::getCp_R().
Referenced by RedlichKisterVPSSTP::cp_mole().
|
virtual |
Return an array of partial molar volumes for the species in the mixture.
Units: m^3/kmol.
Frequently, for this class of thermodynamics representations, the excess Volume due to mixing is zero. Here, we set it as a default. It may be overridden in derived classes.
vbar | Output vector of species partial molar volumes. Length = m_kk. units are m^3/kmol. |
Reimplemented from GibbsExcessVPSSTP.
Definition at line 171 of file RedlichKisterVPSSTP.cpp.
References VPStandardStateTP::getStandardVolumes(), and Phase::m_kk.
|
virtual |
Get the array of temperature second derivatives of the log activity coefficients.
units = 1/Kelvin
d2lnActCoeffdT2 | Output vector of temperature 2nd derivatives of the log Activity Coefficients. length = m_kk |
Definition at line 342 of file RedlichKisterVPSSTP.cpp.
References GibbsExcessVPSSTP::d2lnActCoeffdT2_Scaled_, Phase::m_kk, and RedlichKisterVPSSTP::s_update_dlnActCoeff_dT().
|
virtual |
Get the array of temperature derivatives of the log activity coefficients.
This function is virtual, and first appears in GibbsExcessVPSSTP.
units = 1/Kelvin
dlnActCoeffdT | Output vector of temperature derivatives of the log Activity Coefficients. length = m_kk |
Reimplemented from GibbsExcessVPSSTP.
Definition at line 334 of file RedlichKisterVPSSTP.cpp.
References GibbsExcessVPSSTP::dlnActCoeffdT_Scaled_, Phase::m_kk, and RedlichKisterVPSSTP::s_update_dlnActCoeff_dT().
|
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 VPStandardStateTP.
Definition at line 180 of file RedlichKisterVPSSTP.cpp.
References RedlichKisterVPSSTP::addBinaryInteraction(), AnyMap::hasKey(), RedlichKisterVPSSTP::initLengths(), VPStandardStateTP::initThermo(), ThermoPhase::m_input, and Phase::species().
|
virtual |
Import and initialize a ThermoPhase object using an XML tree.
Here we read extra information about the XML description of a phase. Regular information about elements and species and their reference state thermodynamic information have already been read at this point. For example, we do not need to call this function for ideal gas equations of state. This function is called from importPhase() after the elements and the species are initialized with default ideal solution level data.
The default implementation in ThermoPhase calls the virtual function initThermo() and then sets the "state" of the phase by looking for an XML element named "state", and then interpreting its contents by calling the virtual function setStateFromXML().
phaseNode | This object must be the phase node of a complete XML tree description of the phase, including all of the species data. In other words while "phase" must point to an XML phase object, it must have sibling nodes "speciesData" that describe the species in the phase. |
id | ID of the phase. If nonnull, a check is done to see if phaseNode is pointing to the phase with the correct id. |
Reimplemented from ThermoPhase.
Definition at line 201 of file RedlichKisterVPSSTP.cpp.
References Cantera::caseInsensitiveEquals(), XML_Node::child(), XML_Node::hasChild(), XML_Node::id(), ThermoPhase::initThermoXML(), XML_Node::name(), XML_Node::nChildren(), and RedlichKisterVPSSTP::readXMLBinarySpecies().
void addBinaryInteraction | ( | const std::string & | speciesA, |
const std::string & | speciesB, | ||
const double * | excess_enthalpy, | ||
size_t | n_enthalpy, | ||
const double * | excess_entropy, | ||
size_t | n_entropy | ||
) |
Add a binary species interaction with the specified parameters.
speciesA | name of the first species |
speciesB | name of the second species |
excess_enthalpy | coefficients of the excess enthalpy polynomial |
n_enthalpy | number of excess enthalpy polynomial coefficients |
excess_entropy | coefficients of the excess entropy polynomial |
n_entropy | number of excess entropy polynomial coefficients |
Definition at line 551 of file RedlichKisterVPSSTP.cpp.
References Phase::charge(), RedlichKisterVPSSTP::dlnActCoeff_dX_, RedlichKisterVPSSTP::m_HE_m_ij, RedlichKisterVPSSTP::m_N_ij, RedlichKisterVPSSTP::m_pSpecies_A_ij, RedlichKisterVPSSTP::m_pSpecies_B_ij, RedlichKisterVPSSTP::m_SE_m_ij, Cantera::npos, RedlichKisterVPSSTP::numBinaryInteractions_, Array2D::resize(), and Phase::speciesIndex().
Referenced by RedlichKisterVPSSTP::initThermo(), and RedlichKisterVPSSTP::readXMLBinarySpecies().
|
virtual |
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.
dTds | Input of temperature change along the path |
dXds | Input vector of changes in mole fraction along the path. length = m_kk Along the path length it must be the case that the mole fractions sum to one. |
dlnActCoeffds | Output vector of the directional derivatives of the log Activity Coefficients along the path. length = m_kk units are 1/units(s). if s is a physical coordinate then the units are 1/m. |
Reimplemented from ThermoPhase.
Definition at line 463 of file RedlichKisterVPSSTP.cpp.
References RedlichKisterVPSSTP::dlnActCoeff_dX_, GibbsExcessVPSSTP::dlnActCoeffdT_Scaled_, Phase::m_kk, RedlichKisterVPSSTP::s_update_dlnActCoeff_dT(), and RedlichKisterVPSSTP::s_update_dlnActCoeff_dX_().
|
virtual |
Get the array of ln mole fraction derivatives of the log activity coefficients - diagonal component only.
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 mole fraction variable that represents the standard state. This quantity is to be used in conjunction with derivatives of that mole fraction variable when the derivative of the chemical potential is taken.
units = dimensionless
dlnActCoeffdlnX_diag | Output vector of derivatives of the log Activity Coefficients wrt the mole fractions. length = m_kk |
Reimplemented from ThermoPhase.
Definition at line 487 of file RedlichKisterVPSSTP.cpp.
References GibbsExcessVPSSTP::dlnActCoeffdlnX_diag_, Phase::m_kk, and RedlichKisterVPSSTP::s_update_dlnActCoeff_dlnX_diag().
|
virtual |
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 (i.e. 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
dlnActCoeffdlnN_diag | Output vector of derivatives of the log Activity Coefficients. length = m_kk |
Reimplemented from VPStandardStateTP.
Definition at line 476 of file RedlichKisterVPSSTP.cpp.
References RedlichKisterVPSSTP::dlnActCoeff_dX_, Phase::m_kk, GibbsExcessVPSSTP::moleFractions_, and RedlichKisterVPSSTP::s_update_dlnActCoeff_dX_().
|
virtual |
Get the array of derivatives of the log activity coefficients with respect to the log of the species mole numbers.
Implementations should take the derivative of the logarithm of the activity coefficient with respect to a species log mole number (with all other species mole numbers held constant). The default treatment in the ThermoPhase object is to set this vector to zero.
units = 1 / kmol
dlnActCoeffdlnN[ ld * k + m] will contain the derivative of log act_coeff for the m-th species with respect to the number of moles of the k-th species.
\[ \frac{d \ln(\gamma_m) }{d \ln( n_k ) }\Bigg|_{n_i} \]
ld | Number of rows in the matrix |
dlnActCoeffdlnN | Output vector of derivatives of the log Activity Coefficients. length = m_kk * m_kk |
Reimplemented from GibbsExcessVPSSTP.
Definition at line 495 of file RedlichKisterVPSSTP.cpp.
References GibbsExcessVPSSTP::dlnActCoeffdlnN_, Phase::m_kk, and RedlichKisterVPSSTP::s_update_dlnActCoeff_dX_().
|
private |
Process an XML node called "binaryNeutralSpeciesParameters".
This node contains all of the parameters necessary to describe the Redlich-Kister model for a particular binary interaction. This function reads the XML file and writes the coefficients it finds to an internal data structures.
xmlBinarySpecies | Reference to the XML_Node named "binaryNeutralSpeciesParameters" containing the binary interaction |
Definition at line 506 of file RedlichKisterVPSSTP.cpp.
References RedlichKisterVPSSTP::addBinaryInteraction(), XML_Node::attrib(), XML_Node::child(), Cantera::getFloatArray(), XML_Node::name(), XML_Node::nChildren(), Cantera::npos, Phase::speciesIndex(), and Cantera::toLowerCopy().
Referenced by RedlichKisterVPSSTP::initThermoXML().
|
private |
Initialize lengths of local variables after all species have been identified.
Definition at line 196 of file RedlichKisterVPSSTP.cpp.
References GibbsExcessVPSSTP::dlnActCoeffdlnN_, Phase::m_kk, and Array2D::resize().
Referenced by RedlichKisterVPSSTP::initThermo().
|
private |
Update the activity coefficients.
This function will be called to update the internally stored natural logarithm of the activity coefficients
Definition at line 244 of file RedlichKisterVPSSTP.cpp.
Referenced by RedlichKisterVPSSTP::getChemPotentials(), and RedlichKisterVPSSTP::getLnActivityCoefficients().
|
private |
Update the derivative of the log of the activity coefficients wrt T.
This function will be called to update the internally stored derivative of the natural logarithm of the activity coefficients wrt temperature.
Definition at line 292 of file RedlichKisterVPSSTP.cpp.
References GibbsExcessVPSSTP::d2lnActCoeffdT2_Scaled_, GibbsExcessVPSSTP::dlnActCoeffdT_Scaled_, Phase::m_kk, RedlichKisterVPSSTP::m_N_ij, RedlichKisterVPSSTP::m_pSpecies_A_ij, RedlichKisterVPSSTP::m_pSpecies_B_ij, RedlichKisterVPSSTP::m_SE_m_ij, GibbsExcessVPSSTP::moleFractions_, and RedlichKisterVPSSTP::numBinaryInteractions_.
Referenced by RedlichKisterVPSSTP::getd2lnActCoeffdT2(), RedlichKisterVPSSTP::getdlnActCoeffds(), and RedlichKisterVPSSTP::getdlnActCoeffdT().
|
private |
Internal routine that calculates the derivative of the activity coefficients wrt the mole fractions.
This routine calculates the the derivative of the activity coefficients wrt to mole fraction with all other mole fractions held constant. This is strictly not permitted. However, if the resulting matrix is multiplied by a permissible deltaX vector then everything is ok.
This is the natural way to handle concentration derivatives in this routine.
Definition at line 400 of file RedlichKisterVPSSTP.cpp.
Referenced by RedlichKisterVPSSTP::getdlnActCoeffdlnN(), RedlichKisterVPSSTP::getdlnActCoeffdlnN_diag(), and RedlichKisterVPSSTP::getdlnActCoeffds().
|
private |
Internal routine that calculates the total derivative of the activity coefficients with respect to the log of the mole fractions.
This function will be called to update the internally stored vector of the total derivatives (i.e. not assuming other mole fractions are constant) of the natural logarithm of the activity coefficients with respect to the log of the mole fraction.
Definition at line 350 of file RedlichKisterVPSSTP.cpp.
Referenced by RedlichKisterVPSSTP::getdlnActCoeffdlnX_diag().
|
protected |
number of binary interaction expressions
Definition at line 441 of file RedlichKisterVPSSTP.h.
Referenced by RedlichKisterVPSSTP::addBinaryInteraction(), and RedlichKisterVPSSTP::s_update_dlnActCoeff_dT().
|
protected |
vector of species indices representing species A in the interaction
Each Redlich-Kister excess Gibbs free energy term involves two species, A and B. This vector identifies species A.
Definition at line 448 of file RedlichKisterVPSSTP.h.
Referenced by RedlichKisterVPSSTP::addBinaryInteraction(), and RedlichKisterVPSSTP::s_update_dlnActCoeff_dT().
|
protected |
vector of species indices representing species B in the interaction
Each Redlich-Kister excess Gibbs free energy term involves two species, A and B. This vector identifies species B.
Definition at line 455 of file RedlichKisterVPSSTP.h.
Referenced by RedlichKisterVPSSTP::addBinaryInteraction(), and RedlichKisterVPSSTP::s_update_dlnActCoeff_dT().
|
protected |
Vector of the length of the polynomial for the interaction.
Definition at line 458 of file RedlichKisterVPSSTP.h.
Referenced by RedlichKisterVPSSTP::addBinaryInteraction(), and RedlichKisterVPSSTP::s_update_dlnActCoeff_dT().
|
mutableprotected |
Enthalpy term for the binary mole fraction interaction of the excess Gibbs free energy expression.
Definition at line 462 of file RedlichKisterVPSSTP.h.
Referenced by RedlichKisterVPSSTP::addBinaryInteraction().
|
mutableprotected |
Entropy term for the binary mole fraction interaction of the excess Gibbs free energy expression.
Definition at line 466 of file RedlichKisterVPSSTP.h.
Referenced by RedlichKisterVPSSTP::addBinaryInteraction(), and RedlichKisterVPSSTP::s_update_dlnActCoeff_dT().
|
protected |
form of the RedlichKister interaction expression.
Currently there is only one form.
Definition at line 470 of file RedlichKisterVPSSTP.h.
|
protected |
form of the temperature dependence of the Redlich-Kister interaction expression.
Currently there is only one form -> constant wrt temperature.
Definition at line 475 of file RedlichKisterVPSSTP.h.
|
mutableprotected |
Two dimensional array of derivatives of activity coefficients wrt mole fractions.
Definition at line 479 of file RedlichKisterVPSSTP.h.
Referenced by RedlichKisterVPSSTP::addBinaryInteraction(), RedlichKisterVPSSTP::getdlnActCoeffdlnN_diag(), and RedlichKisterVPSSTP::getdlnActCoeffds().