Cantera  2.1.2
Public Member Functions | Protected Attributes | Private Member Functions | List of all members
PhaseCombo_Interaction Class Reference

PhaseCombo_Interaction is a derived class of GibbsExcessVPSSTP that employs the Margules approximation for the excess gibbs free energy while eliminating the entropy of mixing term. More...

#include <PhaseCombo_Interaction.h>

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

Public Member Functions

 PhaseCombo_Interaction ()
 Constructor. More...
 
 PhaseCombo_Interaction (const std::string &inputFile, const std::string &id="")
 Construct and initialize a PhaseCombo_Interaction ThermoPhase object directly from an xml input file. More...
 
 PhaseCombo_Interaction (XML_Node &phaseRef, const std::string &id="")
 Construct and initialize a PhaseCombo_Interaction ThermoPhase object directly from an XML database. More...
 
 PhaseCombo_Interaction (int testProb)
 Special constructor for a hard-coded problem. More...
 
 PhaseCombo_Interaction (const PhaseCombo_Interaction &b)
 Copy constructor. More...
 
PhaseCombo_Interactionoperator= (const PhaseCombo_Interaction &b)
 Assignment operator. More...
 
virtual ThermoPhaseduplMyselfAsThermoPhase () const
 Duplication routine for objects which inherit from ThermoPhase. More...
 
Utilities
virtual int eosType () const
 Equation of state type flag. 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. 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 getActivityCoefficients (doublereal *ac) const
 Get the array of non-dimensional molar-based 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...
 
void getElectrochemPotentials (doublereal *mu) const
 Get the species electrochemical potentials. 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 files importCTML.cpp and ThermoFactory.cpp.

virtual void initThermo ()
 
void initThermoXML (XML_Node &phaseNode, const std::string &id)
 Import and initialize a ThermoPhase object. More...
 
Derivatives of Thermodynamic Variables needed for Applications
virtual void getdlnActCoeffds (const doublereal dTds, const doublereal *const dXds, doublereal *dlnActCoeffds) const
 Get the change in activity coefficients w.r.t. More...
 
virtual void getdlnActCoeffdlnX_diag (doublereal *dlnActCoeffdlnX_diag) const
 Get the array of log concentration-like derivatives of the log activity coefficients - diagonal component. More...
 
virtual void getdlnActCoeffdlnN_diag (doublereal *dlnActCoeffdlnN_diag) const
 Get the array of derivatives of the log activity coefficients wrt mole numbers - diagonal 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 ln species mole numbers. More...
 
- Public Member Functions inherited from GibbsExcessVPSSTP
 GibbsExcessVPSSTP ()
 
 GibbsExcessVPSSTP (const GibbsExcessVPSSTP &b)
 Copy constructor. More...
 
GibbsExcessVPSSTPoperator= (const GibbsExcessVPSSTP &b)
 Assignment operator. 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
 Returns the natural logarithm of the standard concentration of the kth species. More...
 
virtual void getUnitsStandardConc (double *uA, int k=0, int sizeUA=6) const
 Returns the units of the standard and generalized concentrations Note they have the same units, as their ratio is defined to be equal to the activity of the kth species in the solution, which is unitless. 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 getdlnActCoeffdlnX (doublereal *dlnActCoeffdlnX) const
 Get the array of log concentration-like derivatives of the log activity coefficients. More...
 
void getElectrochemPotentials (doublereal *mu) const
 Get the species electrochemical potentials. More...
 
virtual const vector_fpgetPartialMolarVolumes () const
 
virtual void setState_TP (doublereal t, doublereal p)
 Set the temperature (K) and pressure (Pa) More...
 
virtual void setMassFractions (const doublereal *const y)
 Set the mass fractions to the specified values, and then normalize them so that they sum to 1.0. More...
 
virtual void setMassFractions_NoNorm (const doublereal *const y)
 Set the mass fractions to the specified values without normalizing. More...
 
virtual void setMoleFractions (const doublereal *const x)
 Set the mole fractions to the specified values, and then normalize them so that they sum to 1.0. More...
 
virtual void setMoleFractions_NoNorm (const doublereal *const x)
 Set the mole fractions to the specified values without normalizing. More...
 
virtual void setConcentrations (const doublereal *const c)
 Set the concentrations to the specified values within the phase. More...
 
virtual void setPressure (doublereal p)
 Set the internally stored pressure (Pa) at constant temperature and composition. More...
 
- Public Member Functions inherited from VPStandardStateTP
 VPStandardStateTP ()
 Constructor. More...
 
 VPStandardStateTP (const VPStandardStateTP &b)
 Copy Constructor. More...
 
VPStandardStateTPoperator= (const VPStandardStateTP &b)
 Assignment operator. More...
 
virtual ~VPStandardStateTP ()
 Destructor. 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...
 
void getChemPotentials_RT (doublereal *mu) const
 Get the array of non-dimensional species chemical potentials. More...
 
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...
 
void getPureGibbs (doublereal *gpure) const
 Get the standard state Gibbs functions for each species at the current T and P. 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...
 
virtual const vector_fpgetStandardVolumes () const
 
virtual void setTemperature (const doublereal temp)
 Set the temperature of the phase. More...
 
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 void getEnthalpy_RT_ref (doublereal *hrt) const
 Returns the vector of nondimensional enthalpies of the reference state at the current temperature of the solution and the reference pressure for the species. More...
 
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
 
virtual void getEntropy_R_ref (doublereal *er) const
 
virtual void getCp_R_ref (doublereal *cprt) const
 
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...
 
virtual void setParametersFromXML (const XML_Node &eosdata)
 Set equation of state parameter values from XML entries. More...
 
void setVPSSMgr (VPSSMgr *vp_ptr)
 set the VPSS Mgr More...
 
VPSSMgrprovideVPSSMgr ()
 Return a pointer to the VPSSMgr for this phase. More...
 
void createInstallPDSS (size_t k, const XML_Node &s, const XML_Node *phaseNode_ptr)
 
PDSSprovidePDSS (size_t k)
 
const PDSSprovidePDSS (size_t k) const
 
- Public Member Functions inherited from ThermoPhase
 ThermoPhase ()
 Constructor. More...
 
virtual ~ThermoPhase ()
 Destructor. Deletes the species thermo manager. More...
 
 ThermoPhase (const ThermoPhase &right)
 Copy Constructor for the ThermoPhase object. More...
 
ThermoPhaseoperator= (const ThermoPhase &right)
 Assignment operator. More...
 
doublereal _RT () const
 Return the Gas Constant multiplied by the current temperature. 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 int k) const
 Report the 298 K Heat of Formation of the standard state of one species (J kmol-1) More...
 
virtual void modifyOneHf298SS (const int 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 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 cv_vib (int, double) const
 
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 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...
 
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 getdPartialMolarVolumes_dT (doublereal *d_vbar_dT) const
 Return an array of derivatives of partial molar volumes wrt temperature for the species in the mixture. More...
 
virtual void getdPartialMolarVolumes_dP (doublereal *d_vbar_dP) const
 Return an array of derivatives of partial molar volumes wrt pressure for the species in the mixture. More...
 
virtual void getdStandardVolumes_dT (doublereal *d_vol_dT) const
 Get the derivative of the molar volumes of the species standard states wrt temperature at the current T and P of the solution. More...
 
virtual void getdStandardVolumes_dP (doublereal *d_vol_dP) const
 Get the derivative molar volumes of the species standard states wrt pressure at the current T and P of the solution. 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...
 
virtual void setReferenceComposition (const doublereal *const x)
 Sets the reference composition. More...
 
virtual void getReferenceComposition (doublereal *const x) const
 Gets the reference composition. More...
 
doublereal enthalpy_mass () const
 Specific enthalpy. More...
 
doublereal intEnergy_mass () const
 Specific internal energy. More...
 
doublereal entropy_mass () const
 Specific entropy. More...
 
doublereal gibbs_mass () const
 Specific Gibbs function. More...
 
doublereal cp_mass () const
 Specific heat at constant pressure. More...
 
doublereal cv_mass () const
 Specific heat at constant volume. More...
 
virtual void setToEquilState (const doublereal *lambda_RT)
 This method is used by the ChemEquil equilibrium solver. More...
 
void setElementPotentials (const vector_fp &lambda)
 Stores the element potentials in the ThermoPhase object. More...
 
bool getElementPotentials (doublereal *lambda) const
 Returns the element potentials stored in the ThermoPhase object. More...
 
virtual doublereal critTemperature () const
 Critical temperature (K). More...
 
virtual doublereal critPressure () const
 Critical pressure (Pa). 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 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...
 
void setSpeciesThermo (SpeciesThermo *spthermo)
 Install a species thermodynamic property manager. More...
 
virtual SpeciesThermospeciesThermo (int k=-1)
 Return a changeable reference to the calculation manager for species reference-state thermodynamic properties. More...
 
virtual void initThermoFile (const std::string &inputFile, const std::string &id)
 
virtual void installSlavePhases (Cantera::XML_Node *phaseNode)
 Add in species from Slave phases. 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 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) 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...
 
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, 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, 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 (doublereal h, doublereal p, doublereal tol=1.e-4)
 Set the internally stored specific enthalpy (J/kg) and pressure (Pa) of the phase. More...
 
virtual void setState_UV (doublereal u, doublereal v, doublereal tol=1.e-4)
 Set the specific internal energy (J/kg) and specific volume (m^3/kg). More...
 
virtual void setState_SP (doublereal s, doublereal p, doublereal tol=1.e-4)
 Set the specific entropy (J/kg/K) and pressure (Pa). More...
 
virtual void setState_SV (doublereal s, doublereal v, doublereal tol=1.e-4)
 Set the specific entropy (J/kg/K) and specific volume (m^3/kg). More...
 
- Public Member Functions inherited from Phase
 Phase ()
 Default constructor. More...
 
virtual ~Phase ()
 Destructor. More...
 
 Phase (const Phase &right)
 Copy Constructor. More...
 
Phaseoperator= (const Phase &right)
 Assignment operator. More...
 
XML_Nodexml ()
 Returns a reference to the XML_Node stored for the phase. More...
 
void saveState (vector_fp &state) const
 Save the current internal state of the phase Write to vector 'state' the current internal state. More...
 
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...
 
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...
 
doublereal size (size_t k) const
 This routine returns the size of species k. 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 void freezeSpecies ()
 Call when finished adding species. More...
 
bool speciesFrozen ()
 True if freezeSpecies has been called. More...
 
virtual bool ready () const
 
int stateMFNumber () const
 Return the State Mole Fraction Number. 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_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 Throws an exception if m is greater than nElements()-1. More...
 
void checkElementArraySize (size_t mm) const
 Check that an array size is at least nElements() Throws an exception if mm is less than 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 Throws an exception if k is greater than nSpecies()-1. More...
 
void checkSpeciesArraySize (size_t kk) const
 Check that an array size is at least nSpecies() Throws an exception if kk is less than nSpecies(). More...
 
void setMoleFractionsByName (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 (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, 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, 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...
 
void getMoleFractionsByName (compositionMap &x) const
 Get the mole fractions by name. More...
 
doublereal moleFraction (size_t k) const
 Return the mole fraction of a single species. More...
 
doublereal moleFraction (const std::string &name) const
 Return the mole fraction of a single species. More...
 
doublereal massFraction (size_t k) const
 Return the mass fraction of a single species. More...
 
doublereal massFraction (const std::string &name) const
 Return the mass fraction of a single species. More...
 
void getMoleFractions (doublereal *const x) const
 Get the species mole fraction vector. More...
 
void getMassFractions (doublereal *const y) const
 Get the species mass fractions. More...
 
const doublereal * massFractions () const
 Return a const pointer to the mass fraction array. More...
 
void getConcentrations (doublereal *const c) const
 Get the species concentrations (kmol/m^3). More...
 
doublereal concentration (const size_t k) const
 Concentration of species k. More...
 
const doublereal * moleFractdivMMW () const
 Returns a const pointer to the start of the moleFraction/MW array. More...
 
doublereal temperature () const
 Temperature (K). More...
 
virtual doublereal density () const
 Density (kg/m^3). More...
 
doublereal molarDensity () const
 Molar density (kmol/m^3). More...
 
doublereal molarVolume () const
 Molar volume (m^3/kmol). More...
 
virtual void setDensity (const doublereal density_)
 Set the internally stored density (kg/m^3) of the phase Note the density of a phase is an independent variable. More...
 
virtual void setMolarDensity (const doublereal 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_Y (const doublereal *const Q) const
 Evaluate the mass-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...
 
doublereal sum_xlogQ (doublereal *const Q) const
 Evaluate \( \sum_k X_k \log Q_k \). More...
 
void addElement (const std::string &symbol, doublereal weight=-12345.0)
 Add an element. More...
 
void addElement (const XML_Node &e)
 Add an element from an XML specification. More...
 
void addUniqueElement (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, checking for uniqueness The uniqueness is checked by comparing the string symbol. More...
 
void addUniqueElement (const XML_Node &e)
 Add an element, checking for uniqueness The uniqueness is checked by comparing the string symbol. More...
 
void addElementsFromXML (const XML_Node &phase)
 Add all elements referenced in an XML_Node tree. More...
 
void freezeElements ()
 Prohibit addition of more elements, and prepare to add species. More...
 
bool elementsFrozen ()
 True if freezeElements has been called. More...
 
size_t addUniqueElementAfterFreeze (const std::string &symbol, doublereal weight, int atomicNumber, doublereal entropy298=ENTROPY298_UNKNOWN, int elem_type=CT_ELEM_TYPE_ABSPOS)
 Add an element after elements have been frozen, checking for uniqueness The uniqueness is checked by comparing the string symbol. More...
 
void addSpecies (const std::string &name, const doublereal *comp, doublereal charge=0.0, doublereal size=1.0)
 
void addUniqueSpecies (const std::string &name, const doublereal *comp, doublereal charge=0.0, doublereal size=1.0)
 Add a species to the phase, checking for uniqueness of the name This routine checks for uniqueness of the string name. More...
 

Protected Attributes

size_t numBinaryInteractions_
 number of binary interaction expressions More...
 
vector_fp m_HE_b_ij
 Enthalpy term for the binary mole fraction interaction of the excess gibbs free energy expression. More...
 
vector_fp m_HE_c_ij
 Enthalpy term for the ternary mole fraction interaction of the excess gibbs free energy expression. More...
 
vector_fp m_HE_d_ij
 Enthalpy term for the quaternary mole fraction interaction of the excess gibbs free energy expression. More...
 
vector_fp m_SE_b_ij
 Entropy term for the binary mole fraction interaction of the excess gibbs free energy expression. More...
 
vector_fp m_SE_c_ij
 Entropy term for the ternary mole fraction interaction of the excess gibbs free energy expression. More...
 
vector_fp m_SE_d_ij
 Entropy term for the quaternary mole fraction interaction of the excess gibbs free energy expression. More...
 
vector_fp m_VHE_b_ij
 Enthalpy term for the binary mole fraction interaction of the excess gibbs free energy expression. More...
 
vector_fp m_VHE_c_ij
 Enthalpy term for the ternary mole fraction interaction of the excess gibbs free energy expression. More...
 
vector_fp m_VHE_d_ij
 Enthalpy term for the quaternary mole fraction interaction of the excess gibbs free energy expression. More...
 
vector_fp m_VSE_b_ij
 Entropy term for the binary mole fraction interaction of the excess gibbs free energy expression. More...
 
vector_fp m_VSE_c_ij
 Entropy term for the ternary mole fraction interaction of the excess gibbs free energy expression. More...
 
vector_fp m_VSE_d_ij
 Entropy term for the quaternary mole fraction interaction of the excess gibbs free energy expression. 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...
 
int formMargules_
 form of the Margules interaction expression More...
 
int formTempModel_
 form of the temperature dependence of the Margules interaction expression More...
 
- Protected Attributes inherited from GibbsExcessVPSSTP
std::vector< doublereal > moleFractions_
 Storage for the current values of the mole fractions of the species. More...
 
std::vector< doublereal > lnActCoeff_Scaled_
 Storage for the current values of the activity coefficients of the species. More...
 
std::vector< doublereal > 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...
 
std::vector< doublereal > 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...
 
std::vector< doublereal > 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...
 
std::vector< doublereal > 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...
 
std::vector< doublereal > m_pp
 Temporary storage space that is fair game. More...
 
- Protected Attributes inherited from VPStandardStateTP
doublereal m_Pcurrent
 Current value of the pressure - state variable. More...
 
doublereal m_Tlast_ss
 The last temperature at which the standard statethermodynamic properties were calculated at. More...
 
doublereal m_Plast_ss
 The last pressure at which the Standard State thermodynamic properties were calculated at. More...
 
doublereal m_P0
 
VPSSMgrm_VPSS_ptr
 Pointer to the VPSS manager that calculates all of the standard state info efficiently. More...
 
std::vector< PDSS * > m_PDSS_storage
 Storage for the PDSS objects for the species. More...
 
- Protected Attributes inherited from ThermoPhase
SpeciesThermom_spthermo
 Pointer to the calculation manager for species reference-state thermodynamic properties. 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. More...
 
vector_fp m_lambdaRRT
 Vector of element potentials. More...
 
bool m_hasElementPotentials
 Boolean indicating whether there is a valid set of saved element potentials for this phase. More...
 
bool m_chargeNeutralityNecessary
 Boolean indicating whether a charge neutrality condition is a necessity. More...
 
int m_ssConvention
 Contains the standard state convention. More...
 
std::vector< doublereal > xMol_Ref
 Reference Mole Fraction Composition. More...
 
- Protected Attributes inherited from Phase
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_speciesSize
 Vector of species sizes. More...
 
vector_fp m_speciesCharge
 Vector of species charges. length m_kk. More...
 

Private Member Functions

void readXMLBinarySpecies (XML_Node &xmlBinarySpecies)
 Process an XML node called "binaryNeutralSpeciesParameters". More...
 
void resizeNumInteractions (const size_t num)
 Resize internal arrays within the object that depend upon the number of binary Margules interaction terms. 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_dlnX_diag () const
 Update the derivative of the log of the activity coefficients wrt log(mole fraction) More...
 
void s_update_dlnActCoeff_dlnN_diag () const
 Update the derivative of the log of the activity coefficients wrt log(moles) - diagonal only. More...
 
void s_update_dlnActCoeff_dlnN () const
 Update the derivative of the log of the activity coefficients wrt log(moles_m) More...
 
doublereal err (const std::string &msg) const
 Error function. More...
 

Additional Inherited Members

- Protected Member Functions inherited from GibbsExcessVPSSTP
double checkMFSum (const doublereal *const x) const
 utility routine to check mole fraction sum More...
 
void calcDensity ()
 Calculate the density of the mixture using the partial molar volumes and mole fractions as input. More...
 
- Protected Member Functions inherited from VPStandardStateTP
virtual void _updateStandardStateThermo () const
 Updates the standard state thermodynamic functions at the current T and P of the solution. More...
 
const vector_fpGibbs_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 init (const vector_fp &mw)
 
void setMolecularWeight (const int k, const double mw)
 Set the molecular weight of a single species to a given value. More...
 

Detailed Description

PhaseCombo_Interaction is a derived class of GibbsExcessVPSSTP that employs the Margules approximation for the excess gibbs free energy while eliminating the entropy of mixing term.

PhaseCombo_Interaction derives from class GibbsExcessVPSSTP which is derived from VPStandardStateTP, and overloads the virtual methods defined there with ones that use expressions appropriate for the Margules Excess gibbs free energy approximation. The reader should refer to the MargulesVPSSTP class for information on that class. This class in addition adds a term to the activity coefficient that eliminates the ideal solution mixing term within the chemical potential. This is a very radical thing to do, but it is supported by experimental evidence under some conditions.

The independent unknowns are pressure, temperature, and mass fraction.

Several concepts are introduced. The first concept is that 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. Currently, these variables and the calculation method are handled by the VPSSMgr class, for which VPStandardStateTP owns a pointer to.

To support the above functionality, pressure and temperature variables, m_plast_ss and m_tlast_ss, are kept which store the last pressure and temperature used in the evaluation of standard state properties.

This class is introduced to represent specific conditions observed in thermal batteries. HOwever, it may be physically motivated to represent conditions where there may be a mixture of compounds that are not "mixed" at the molecular level. Therefore, there is no mixing term.

The lack of a mixing term has profound effects. First, the mole fraction of a species can now be identically zero due to thermodynamic considerations. The phase behaves more like a series of phases. That's why we named it PhaseCombo.


Specification of Species Standard State Properties


All species are defined to have standard states that depend upon both the temperature and the pressure. The Margules 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.


Specification of Solution Thermodynamic Properties


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 Margules formulation for a phase that has more than 2 species. The second term in the excess gibbs free energy is a negation of the ideal solution's mixing term.

\[ G^E = \sum_i \left( H_{Ei} - T S_{Ei} \right) - \sum_i \left( n_i R T \ln{X_i} \right) \]

\[ H^E_i = n X_{Ai} X_{Bi} \left( h_{o,i} + h_{1,i} X_{Bi} \right) \]

\[ S^E_i = n X_{Ai} X_{Bi} \left( s_{o,i} + s_{1,i} X_{Bi} \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 \left( \left( \delta_{Ai,k} X_{Bi} + \delta_{Bi,k} X_{Ai} - X_{Ai} X_{Bi} \right) \left( g^E_{o,i} + g^E_{1,i} X_{Bi} \right) + \left( \delta_{Bi,k} - X_{Bi} \right) X_{Ai} X_{Bi} g^E_{1,i} \right) - RT \ln{X_k} \]

where \( g^E_{o,i} = h_{o,i} - T s_{o,i} \) and \( g^E_{1,i} = h_{1,i} - T s_{1,i} \) and where \( X_k \) is the mole fraction of species k.

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} \]


Application within Kinetics Managers


\( 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 reate constant and the concentration equilibrium constant, \( K_c \).

\[ k^{-1} = k^1 K^1_c \]

\(k^{-1} \) has units of s-1.


Instantiation of the Class


The constructor for this phase is located in the default ThermoFactory for Cantera. A new PhaseCombo_Interaction object may be created by the following code snippet:

XML_Node *xc = get_XML_File("LiFeS_X_combo.xml");
XML_Node * const xs = xc->findNameID("phase", "LiFeS_X");
ThermoPhase *l_tp = newPhase(*xs);
PhaseCombo_Interaction *LiFeS_X_solid = dynamic_cast <PhaseCombo_Interaction *>(l_tp);

or by the following code

std::string id = "LiFeS_X";
Cantera::ThermoPhase *LiFeS_X_Phase = Cantera::newPhase("LiFeS_X_combo.xml", id);
PhaseCombo_Interaction *LiFeS_X_solid = dynamic_cast <PhaseCombo_Interaction *>(l_tp);

or by the following constructor:

XML_Node *xc = get_XML_File("LiFeS_X_combo.xml");
XML_Node * const xs = xc->findNameID("phase", "LiFeS_X");

XML Example


An example of an XML Element named phase setting up a PhaseCombo_Interaction object named LiFeS_X is given below.

<phase dim="3" id="LiFeS_X">
<elementArray datasrc="elements.xml">
Li Fe S
</elementArray>
<speciesArray datasrc="#species_LiFeS">
LiTFe1S2(S) Li2Fe1S2(S)
</speciesArray>
<thermo model="PhaseCombo_Interaction">
<activityCoefficients model="Margules" TempModel="constant">
<binaryNeutralSpeciesParameters speciesA="LiTFe1S2(S)" speciesB="Li2Fe1S2(S)">
<excessEnthalpy model="poly_Xb" terms="2" units="kJ/mol">
84.67069219, -269.1959421
</excessEnthalpy>
<excessEntropy model="poly_Xb" terms="2" units="J/mol/K">
100.7511565, -361.4222659
</excessEntropy>
<excessVolume_Enthalpy model="poly_Xb" terms="2" units="ml/mol">
0, 0
</excessVolume_Enthalpy>
<excessVolume_Entropy model="poly_Xb" terms="2" units="ml/mol/K">
0, 0
</excessVolume_Entropy>
</binaryNeutralSpeciesParameters>
</activityCoefficients>
</thermo>
<transport model="none"/>
<kinetics model="none"/>
</phase>

The model attribute "PhaseCombo_Interaction" of the thermo XML element identifies the phase as being of the type handled by the PhaseCombo_Interaction object.

Definition at line 333 of file PhaseCombo_Interaction.h.

Constructor & Destructor Documentation

Constructor.

This doesn't do much more than initialize constants with default values for water at 25C. Water molecular weight comes from the default elements.xml file. It actually differs slightly from the IAPWS95 value of 18.015268. However, density conservation and therefore element conservation is the more important principle to follow.

Definition at line 21 of file PhaseCombo_Interaction.cpp.

Referenced by PhaseCombo_Interaction::duplMyselfAsThermoPhase().

PhaseCombo_Interaction ( const std::string &  inputFile,
const std::string &  id = "" 
)

Construct and initialize a PhaseCombo_Interaction ThermoPhase object directly from an xml input file.

Parameters
inputFileName of the input file containing the phase XML data to set up the object
idID of the phase in the input file. Defaults to the empty string.

Definition at line 29 of file PhaseCombo_Interaction.cpp.

References ThermoPhase::initThermoFile().

PhaseCombo_Interaction ( XML_Node phaseRef,
const std::string &  id = "" 
)

Construct and initialize a PhaseCombo_Interaction ThermoPhase object directly from an XML database.

Parameters
phaseRefXML phase node containing the description of the phase
idid attribute containing the name of the phase. (default is the empty string)

Definition at line 39 of file PhaseCombo_Interaction.cpp.

References Cantera::findXMLPhase(), and Cantera::importPhase().

PhaseCombo_Interaction ( int  testProb)

Copy constructor.

Parameters
bclass to be copied

Definition at line 49 of file PhaseCombo_Interaction.cpp.

References PhaseCombo_Interaction::operator=().

Member Function Documentation

PhaseCombo_Interaction & operator= ( const PhaseCombo_Interaction b)
ThermoPhase * duplMyselfAsThermoPhase ( ) const
virtual

Duplication routine for objects which inherit from ThermoPhase.

This virtual routine can be used to duplicate thermophase objects inherited from ThermoPhase even if the application only has a pointer to ThermoPhase to work with.

Reimplemented from GibbsExcessVPSSTP.

Definition at line 86 of file PhaseCombo_Interaction.cpp.

References PhaseCombo_Interaction::PhaseCombo_Interaction().

int eosType ( ) const
virtual

Equation of state type flag.

The ThermoPhase base class returns zero. Subclasses should define this to return a unique non-zero value. Known constants defined for this purpose are listed in mix_defs.h.

Reimplemented from GibbsExcessVPSSTP.

Definition at line 155 of file PhaseCombo_Interaction.cpp.

Referenced by PhaseCombo_Interaction::err().

doublereal enthalpy_mole ( ) const
virtual

Molar enthalpy. Units: J/kmol.

Reimplemented from ThermoPhase.

Definition at line 216 of file PhaseCombo_Interaction.cpp.

References PhaseCombo_Interaction::getPartialMolarEnthalpies(), GibbsExcessVPSSTP::moleFractions_, and Phase::nSpecies().

doublereal entropy_mole ( ) const
virtual

Molar entropy. Units: J/kmol.

Reimplemented from ThermoPhase.

Definition at line 228 of file PhaseCombo_Interaction.cpp.

References PhaseCombo_Interaction::getPartialMolarEntropies(), GibbsExcessVPSSTP::moleFractions_, and Phase::nSpecies().

doublereal cp_mole ( ) const
virtual

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

Reimplemented from ThermoPhase.

Definition at line 240 of file PhaseCombo_Interaction.cpp.

References PhaseCombo_Interaction::getPartialMolarCp(), GibbsExcessVPSSTP::moleFractions_, and Phase::nSpecies().

Referenced by PhaseCombo_Interaction::cv_mole().

doublereal cv_mole ( ) const
virtual

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

Reimplemented from ThermoPhase.

Definition at line 252 of file PhaseCombo_Interaction.cpp.

References PhaseCombo_Interaction::cp_mole(), and Cantera::GasConstant.

void getActivityCoefficients ( doublereal *  ac) const
virtual

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

Parameters
acOutput vector of activity coefficients. Length: m_kk.

Reimplemented from GibbsExcessVPSSTP.

Definition at line 164 of file PhaseCombo_Interaction.cpp.

References GibbsExcessVPSSTP::lnActCoeff_Scaled_, Phase::m_kk, and PhaseCombo_Interaction::s_update_lnActCoeff().

void getChemPotentials ( doublereal *  mu) const
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.

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

Reimplemented from ThermoPhase.

Definition at line 192 of file PhaseCombo_Interaction.cpp.

References Cantera::GasConstant, VPStandardStateTP::getStandardChemPotentials(), GibbsExcessVPSSTP::lnActCoeff_Scaled_, Phase::m_kk, GibbsExcessVPSSTP::moleFractions_, PhaseCombo_Interaction::s_update_lnActCoeff(), Cantera::SmallNumber, and Phase::temperature().

Referenced by PhaseCombo_Interaction::getElectrochemPotentials().

void getPartialMolarEnthalpies ( doublereal *  hbar) const
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} \]

Parameters
hbarVector of returned partial molar enthalpies (length m_kk, units = J/kmol)

Reimplemented from ThermoPhase.

Definition at line 257 of file PhaseCombo_Interaction.cpp.

References GibbsExcessVPSSTP::dlnActCoeffdT_Scaled_, Cantera::GasConstant, VPStandardStateTP::getEnthalpy_RT(), Phase::m_kk, PhaseCombo_Interaction::s_update_dlnActCoeff_dT(), PhaseCombo_Interaction::s_update_lnActCoeff(), and Phase::temperature().

Referenced by PhaseCombo_Interaction::enthalpy_mole().

void getPartialMolarEntropies ( doublereal *  sbar) const
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} \]

Parameters
sbarVector of returned partial molar entropies (length m_kk, units = J/kmol/K)

Reimplemented from ThermoPhase.

Definition at line 308 of file PhaseCombo_Interaction.cpp.

References GibbsExcessVPSSTP::dlnActCoeffdT_Scaled_, Cantera::GasConstant, VPStandardStateTP::getEntropy_R(), GibbsExcessVPSSTP::lnActCoeff_Scaled_, Phase::m_kk, GibbsExcessVPSSTP::moleFractions_, PhaseCombo_Interaction::s_update_dlnActCoeff_dT(), PhaseCombo_Interaction::s_update_lnActCoeff(), Cantera::SmallNumber, and Phase::temperature().

Referenced by PhaseCombo_Interaction::entropy_mole().

void getPartialMolarCp ( doublereal *  cpbar) const
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} ??????????????? \]

Parameters
cpbarVector of returned partial molar heat capacities (length m_kk, units = J/kmol/K)

Reimplemented from ThermoPhase.

Definition at line 283 of file PhaseCombo_Interaction.cpp.

References GibbsExcessVPSSTP::d2lnActCoeffdT2_Scaled_, GibbsExcessVPSSTP::dlnActCoeffdT_Scaled_, Cantera::GasConstant, VPStandardStateTP::getCp_R(), Phase::m_kk, PhaseCombo_Interaction::s_update_dlnActCoeff_dT(), PhaseCombo_Interaction::s_update_lnActCoeff(), and Phase::temperature().

Referenced by PhaseCombo_Interaction::cp_mole().

void getPartialMolarVolumes ( doublereal *  vbar) const
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.

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

Reimplemented from GibbsExcessVPSSTP.

Definition at line 335 of file PhaseCombo_Interaction.cpp.

References VPStandardStateTP::getStandardVolumes(), Phase::m_kk, PhaseCombo_Interaction::m_pSpecies_A_ij, PhaseCombo_Interaction::m_pSpecies_B_ij, PhaseCombo_Interaction::m_VHE_b_ij, PhaseCombo_Interaction::m_VHE_c_ij, PhaseCombo_Interaction::m_VSE_b_ij, PhaseCombo_Interaction::m_VSE_c_ij, GibbsExcessVPSSTP::moleFractions_, PhaseCombo_Interaction::numBinaryInteractions_, and Phase::temperature().

void getElectrochemPotentials ( doublereal *  mu) const

Get the species electrochemical potentials.

These are partial molar quantities. This method adds a term \( Fz_k \phi_k \) to the to each chemical potential.

Units: J/kmol

Parameters
muoutput vector containing the species electrochemical potentials. Length: m_kk., units = J/kmol

Definition at line 183 of file PhaseCombo_Interaction.cpp.

References Phase::charge(), ThermoPhase::electricPotential(), PhaseCombo_Interaction::getChemPotentials(), and Phase::m_kk.

void getd2lnActCoeffdT2 ( doublereal *  d2lnActCoeffdT2) const
virtual

Get the array of temperature second derivatives of the log activity coefficients.

This function is a virtual class, but it first appears in GibbsExcessVPSSTP class and derived classes from GibbsExcessVPSSTP.

units = 1/Kelvin

Parameters
d2lnActCoeffdT2Output vector of temperature 2nd derivatives of the log Activity Coefficients. length = m_kk

Definition at line 536 of file PhaseCombo_Interaction.cpp.

References GibbsExcessVPSSTP::d2lnActCoeffdT2_Scaled_, Phase::m_kk, and PhaseCombo_Interaction::s_update_dlnActCoeff_dT().

void getdlnActCoeffdT ( doublereal *  dlnActCoeffdT) const
virtual

Get the array of temperature derivatives of the log activity coefficients.

This function is a virtual class, but it first appears in GibbsExcessVPSSTP class and derived classes from GibbsExcessVPSSTP.

units = 1/Kelvin

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

Reimplemented from GibbsExcessVPSSTP.

Definition at line 528 of file PhaseCombo_Interaction.cpp.

References GibbsExcessVPSSTP::dlnActCoeffdT_Scaled_, Phase::m_kk, and PhaseCombo_Interaction::s_update_dlnActCoeff_dT().

void initThermo ( )
virtual

Initialize. This method is provided to allow subclasses to perform any initialization required after all species have been added. For example, it might be used to resize internal work arrays that must have an entry for each species. The base class implementation does nothing, and subclasses that do not require initialization do not need to overload this method. When importing a CTML phase description, this method is called just prior to returning from function importPhase.

See Also
importCTML.cpp

Reimplemented from GibbsExcessVPSSTP.

Definition at line 378 of file PhaseCombo_Interaction.cpp.

References PhaseCombo_Interaction::initLengths(), and GibbsExcessVPSSTP::initThermo().

void initThermoXML ( XML_Node phaseNode,
const std::string &  id 
)
virtual

Import and initialize a ThermoPhase object.

Parameters
phaseNodeThis 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.
idID of the phase. If nonnull, a check is done to see if phaseNode is pointing to the phase with the correct id.

Reimplemented from VPStandardStateTP.

Definition at line 390 of file PhaseCombo_Interaction.cpp.

References XML_Node::attrib(), XML_Node::child(), XML_Node::hasChild(), XML_Node::id(), VPStandardStateTP::initThermoXML(), Cantera::lowercase(), XML_Node::name(), XML_Node::nChildren(), PhaseCombo_Interaction::readXMLBinarySpecies(), and Phase::size().

void getdlnActCoeffds ( const doublereal  dTds,
const doublereal *const  dXds,
doublereal *  dlnActCoeffds 
) const
virtual

Get the change in activity coefficients w.r.t.

change in state (temp, mole fraction, etc.) along a line in parameter space or along a line in physical space

Parameters
dTdsInput of temperature change along the path
dXdsInput 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.
dlnActCoeffdsOutput 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 544 of file PhaseCombo_Interaction.cpp.

References GibbsExcessVPSSTP::dlnActCoeffdT_Scaled_, Cantera::GasConstant, PhaseCombo_Interaction::m_HE_b_ij, PhaseCombo_Interaction::m_HE_c_ij, Phase::m_kk, PhaseCombo_Interaction::m_pSpecies_A_ij, PhaseCombo_Interaction::m_pSpecies_B_ij, PhaseCombo_Interaction::m_SE_b_ij, PhaseCombo_Interaction::m_SE_c_ij, GibbsExcessVPSSTP::moleFractions_, PhaseCombo_Interaction::numBinaryInteractions_, PhaseCombo_Interaction::s_update_dlnActCoeff_dT(), Cantera::SmallNumber, and Phase::temperature().

void getdlnActCoeffdlnX_diag ( doublereal *  dlnActCoeffdlnX_diag) const
virtual

Get the array of log concentration-like derivatives of the log activity coefficients - diagonal component.

This function is a virtual method. 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.

units = dimensionless

Parameters
dlnActCoeffdlnX_diagOutput vector of the diagonal component of the log(mole fraction) derivatives of the log Activity Coefficients. length = m_kk

Reimplemented from ThermoPhase.

Definition at line 746 of file PhaseCombo_Interaction.cpp.

References GibbsExcessVPSSTP::dlnActCoeffdlnX_diag_, Phase::m_kk, and PhaseCombo_Interaction::s_update_dlnActCoeff_dlnX_diag().

void getdlnActCoeffdlnN_diag ( doublereal *  dlnActCoeffdlnN_diag) const
virtual

Get the array of derivatives of the log activity coefficients wrt mole numbers - diagonal only.

This function is a virtual method. 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. mole fraction, molality, etc.) that represents the standard state.

units = dimensionless

Parameters
dlnActCoeffdlnN_diagOutput vector of the diagonal entries for the log(mole fraction) derivatives of the log Activity Coefficients. length = m_kk

Reimplemented from VPStandardStateTP.

Definition at line 738 of file PhaseCombo_Interaction.cpp.

References GibbsExcessVPSSTP::dlnActCoeffdlnN_diag_, Phase::m_kk, and PhaseCombo_Interaction::s_update_dlnActCoeff_dlnN_diag().

void getdlnActCoeffdlnN ( const size_t  ld,
doublereal *const  dlnActCoeffdlnN 
)
virtual

Get the array of derivatives of the log activity coefficients with respect to the ln species mole numbers.

Implementations should take the derivative of the logarithm of the activity coefficient with respect to a log of a species mole number (with all other species mole numbers held constant)

units = 1 / kmol

dlnActCoeffdlnN[ ld * k + m] will contain the derivative of log act_coeff for the mth species with respect to the number of moles of the kth species.

\[ \frac{d \ln(\gamma_m) }{d \ln( n_k ) }\Bigg|_{n_i} \]

Parameters
ldNumber of rows in the matrix
dlnActCoeffdlnNOutput vector of derivatives of the log Activity Coefficients. length = m_kk * m_kk

Reimplemented from GibbsExcessVPSSTP.

Definition at line 754 of file PhaseCombo_Interaction.cpp.

References GibbsExcessVPSSTP::dlnActCoeffdlnN_, Phase::m_kk, and PhaseCombo_Interaction::s_update_dlnActCoeff_dlnN().

void readXMLBinarySpecies ( XML_Node xmlBinarySpecies)
private

Process an XML node called "binaryNeutralSpeciesParameters".

This node contains all of the parameters necessary to describe the Margules model for a particular binary interaction. This function reads the XML file and writes the coefficients it finds to an internal data structures.

Parameters
xmlBinarySpeciesReference to the XML_Node named "binaryNeutralSpeciesParameters" containing the binary interaction

Definition at line 785 of file PhaseCombo_Interaction.cpp.

References XML_Node::attrib(), Phase::charge(), XML_Node::child(), ctml::getFloatArray(), Cantera::lowercase(), PhaseCombo_Interaction::m_HE_b_ij, PhaseCombo_Interaction::m_HE_c_ij, PhaseCombo_Interaction::m_pSpecies_A_ij, PhaseCombo_Interaction::m_pSpecies_B_ij, PhaseCombo_Interaction::m_SE_b_ij, PhaseCombo_Interaction::m_SE_c_ij, PhaseCombo_Interaction::m_VHE_b_ij, PhaseCombo_Interaction::m_VHE_c_ij, PhaseCombo_Interaction::m_VSE_b_ij, PhaseCombo_Interaction::m_VSE_c_ij, XML_Node::name(), XML_Node::nChildren(), Cantera::npos, PhaseCombo_Interaction::numBinaryInteractions_, PhaseCombo_Interaction::resizeNumInteractions(), Phase::speciesIndex(), and Phase::speciesName().

Referenced by PhaseCombo_Interaction::initThermoXML().

void resizeNumInteractions ( const size_t  num)
private
void initLengths ( )
private

Initialize lengths of local variables after all species have been identified.

Definition at line 384 of file PhaseCombo_Interaction.cpp.

References GibbsExcessVPSSTP::dlnActCoeffdlnN_, Phase::m_kk, Phase::nSpecies(), and Array2D::resize().

Referenced by PhaseCombo_Interaction::initThermo().

void s_update_lnActCoeff ( ) const
private
void s_update_dlnActCoeff_dT ( ) const
private
void s_update_dlnActCoeff_dlnX_diag ( ) const
private

Update the derivative of the log of the activity coefficients wrt log(mole fraction)

This function will be called to update the internally stored derivative of the natural logarithm of the activity coefficients wrt logarithm of the mole fractions.

Definition at line 714 of file PhaseCombo_Interaction.cpp.

References GibbsExcessVPSSTP::dlnActCoeffdlnX_diag_, Cantera::GasConstant, PhaseCombo_Interaction::m_HE_b_ij, PhaseCombo_Interaction::m_HE_c_ij, Phase::m_kk, PhaseCombo_Interaction::m_pSpecies_A_ij, PhaseCombo_Interaction::m_pSpecies_B_ij, PhaseCombo_Interaction::m_SE_b_ij, PhaseCombo_Interaction::m_SE_c_ij, GibbsExcessVPSSTP::moleFractions_, PhaseCombo_Interaction::numBinaryInteractions_, and Phase::temperature().

Referenced by PhaseCombo_Interaction::getdlnActCoeffdlnX_diag().

void s_update_dlnActCoeff_dlnN_diag ( ) const
private

Update the derivative of the log of the activity coefficients wrt log(moles) - diagonal only.

This function will be called to update the internally stored diagonal entries for the derivative of the natural logarithm of the activity coefficients wrt logarithm of the moles.

Definition at line 596 of file PhaseCombo_Interaction.cpp.

References GibbsExcessVPSSTP::dlnActCoeffdlnN_diag_, Cantera::GasConstant, PhaseCombo_Interaction::m_HE_b_ij, PhaseCombo_Interaction::m_HE_c_ij, Phase::m_kk, PhaseCombo_Interaction::m_pSpecies_A_ij, PhaseCombo_Interaction::m_pSpecies_B_ij, PhaseCombo_Interaction::m_SE_b_ij, PhaseCombo_Interaction::m_SE_c_ij, GibbsExcessVPSSTP::moleFractions_, PhaseCombo_Interaction::numBinaryInteractions_, Cantera::SmallNumber, and Phase::temperature().

Referenced by PhaseCombo_Interaction::getdlnActCoeffdlnN_diag().

void s_update_dlnActCoeff_dlnN ( ) const
private

Update the derivative of the log of the activity coefficients wrt log(moles_m)

This function will be called to update the internally stored derivative of the natural logarithm of the activity coefficients wrt logarithm of the mole number of species

Definition at line 647 of file PhaseCombo_Interaction.cpp.

References GibbsExcessVPSSTP::dlnActCoeffdlnN_, Cantera::GasConstant, PhaseCombo_Interaction::m_HE_b_ij, PhaseCombo_Interaction::m_HE_c_ij, Phase::m_kk, PhaseCombo_Interaction::m_pSpecies_A_ij, PhaseCombo_Interaction::m_pSpecies_B_ij, PhaseCombo_Interaction::m_SE_b_ij, PhaseCombo_Interaction::m_SE_c_ij, GibbsExcessVPSSTP::moleFractions_, PhaseCombo_Interaction::numBinaryInteractions_, Cantera::SmallNumber, Phase::temperature(), and Array2D::zero().

Referenced by PhaseCombo_Interaction::getdlnActCoeffdlnN().

doublereal err ( const std::string &  msg) const
private

Error function.

Print an error string and exit

Parameters
msgMessage to be printed

Definition at line 371 of file PhaseCombo_Interaction.cpp.

References PhaseCombo_Interaction::eosType(), and Cantera::int2str().

Member Data Documentation

size_t numBinaryInteractions_
protected
vector_fp m_HE_b_ij
mutableprotected
vector_fp m_HE_c_ij
mutableprotected
vector_fp m_HE_d_ij
mutableprotected

Enthalpy term for the quaternary mole fraction interaction of the excess gibbs free energy expression.

Definition at line 768 of file PhaseCombo_Interaction.h.

Referenced by PhaseCombo_Interaction::operator=(), PhaseCombo_Interaction::PhaseCombo_Interaction(), and PhaseCombo_Interaction::resizeNumInteractions().

vector_fp m_SE_b_ij
mutableprotected
vector_fp m_SE_c_ij
mutableprotected
vector_fp m_SE_d_ij
mutableprotected

Entropy term for the quaternary mole fraction interaction of the excess gibbs free energy expression.

Definition at line 780 of file PhaseCombo_Interaction.h.

Referenced by PhaseCombo_Interaction::operator=(), PhaseCombo_Interaction::PhaseCombo_Interaction(), and PhaseCombo_Interaction::resizeNumInteractions().

vector_fp m_VHE_b_ij
mutableprotected
vector_fp m_VHE_c_ij
mutableprotected
vector_fp m_VHE_d_ij
mutableprotected

Enthalpy term for the quaternary mole fraction interaction of the excess gibbs free energy expression.

Definition at line 792 of file PhaseCombo_Interaction.h.

Referenced by PhaseCombo_Interaction::operator=(), PhaseCombo_Interaction::PhaseCombo_Interaction(), and PhaseCombo_Interaction::resizeNumInteractions().

vector_fp m_VSE_b_ij
mutableprotected
vector_fp m_VSE_c_ij
mutableprotected
vector_fp m_VSE_d_ij
mutableprotected

Entropy term for the quaternary mole fraction interaction of the excess gibbs free energy expression.

Definition at line 804 of file PhaseCombo_Interaction.h.

Referenced by PhaseCombo_Interaction::operator=(), PhaseCombo_Interaction::PhaseCombo_Interaction(), and PhaseCombo_Interaction::resizeNumInteractions().

std::vector<size_t> m_pSpecies_A_ij
protected
std::vector<size_t> m_pSpecies_B_ij
protected
int formMargules_
protected

form of the Margules interaction expression

Currently there is only one form.

Definition at line 824 of file PhaseCombo_Interaction.h.

Referenced by PhaseCombo_Interaction::operator=().

int formTempModel_
protected

form of the temperature dependence of the Margules interaction expression

Currently there is only one form -> constant wrt temperature.

Definition at line 830 of file PhaseCombo_Interaction.h.

Referenced by PhaseCombo_Interaction::operator=().


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