Cantera  2.3.0
Public Member Functions | Public Attributes | Private Member Functions | Private Attributes | List of all members
vcs_VolPhase Class Reference

Phase information and Phase calculations for vcs. More...

#include <vcs_VolPhase.h>

Collaboration diagram for vcs_VolPhase:
[legend]

Public Member Functions

 vcs_VolPhase (VCS_SOLVE *owningSolverObject=0)
 
 vcs_VolPhase (const vcs_VolPhase &b)
 
vcs_VolPhaseoperator= (const vcs_VolPhase &b)
 
void resize (const size_t phaseNum, const size_t numSpecies, const size_t numElem, const char *const phaseName, const double molesInert=0.0)
 The resize() function fills in all of the initial information if it is not given in the constructor. More...
 
void elemResize (const size_t numElemConstraints)
 
double AC_calc_one (size_t kspec) const
 Evaluate activity coefficients and return the kspec coefficient. More...
 
void setMoleFractionsState (const double molNum, const double *const moleFracVec, const int vcsStateStatus)
 Set the moles and/or mole fractions within the phase. More...
 
void setMolesFromVCS (const int stateCalc, const double *molesSpeciesVCS=0)
 Set the moles within the phase. More...
 
void setMolesFromVCSCheck (const int vcsStateStatus, const double *molesSpeciesVCS, const double *const TPhMoles)
 Set the moles within the phase. More...
 
void updateFromVCS_MoleNumbers (const int stateCalc)
 Update the moles within the phase, if necessary. More...
 
void sendToVCS_ActCoeff (const int stateCalc, double *const AC)
 Fill in an activity coefficients vector within a VCS_SOLVE object. More...
 
void setElectricPotential (const double phi)
 set the electric potential of the phase More...
 
double electricPotential () const
 Returns the electric field of the phase. More...
 
double GStar_calc_one (size_t kspec) const
 Gibbs free energy calculation for standard state of one species. More...
 
double G0_calc_one (size_t kspec) const
 Gibbs free energy calculation at a temperature for the reference state of a species, return a value for one species. More...
 
double VolStar_calc_one (size_t kspec) const
 Molar volume calculation for standard state of one species. More...
 
double sendToVCS_VolPM (double *const VolPM) const
 Fill in the partial molar volume vector for VCS. More...
 
void sendToVCS_GStar (double *const gstar) const
 Fill in the partial molar volume vector for VCS. More...
 
void setState_TP (const double temperature_Kelvin, const double pressure_PA)
 Sets the temperature and pressure in this object and underlying ThermoPhase objects. More...
 
void setState_T (const double temperature_Kelvin)
 Sets the temperature in this object and underlying ThermoPhase objects. More...
 
void sendToVCS_LnActCoeffJac (Array2D &LnACJac_VCS)
 
void setPtrThermoPhase (ThermoPhase *tp_ptr)
 Set the pointer for Cantera's ThermoPhase parameter. More...
 
const ThermoPhaseptrThermoPhase () const
 Return a const ThermoPhase pointer corresponding to this phase. More...
 
double totalMoles () const
 Return the total moles in the phase. More...
 
double molefraction (size_t kspec) const
 Returns the mole fraction of the kspec species. More...
 
void setTotalMoles (const double totalMols)
 Sets the total moles in the phase. More...
 
void setMolesOutOfDate (int stateCalc=-1)
 Sets the mole flag within the object to out of date. More...
 
void setMolesCurrent (int vcsStateStatus)
 Sets the mole flag within the object to be current. More...
 
const vector_fpmoleFractions () const
 Return a const reference to the mole fractions stored in the object. More...
 
double moleFraction (size_t klocal) const
 
void setCreationMoleNumbers (const double *const n_k, const std::vector< size_t > &creationGlobalRxnNumbers)
 Sets the creationMoleNum's within the phase object. More...
 
const vector_fpcreationMoleNumbers (std::vector< size_t > &creationGlobalRxnNumbers) const
 Return a const reference to the creationMoleNumbers stored in the object. More...
 
bool isIdealSoln () const
 Returns whether the phase is an ideal solution phase. More...
 
size_t phiVarIndex () const
 Return the index of the species that represents the the voltage of the phase. More...
 
void setPhiVarIndex (size_t phiVarIndex)
 
vcs_SpeciesPropertiesspeciesProperty (const size_t kindex)
 Retrieve the kth Species structure for the species belonging to this phase. More...
 
int exists () const
 int indicating whether the phase exists or not More...
 
void setExistence (const int existence)
 Set the existence flag in the object. More...
 
size_t spGlobalIndexVCS (const size_t spIndex) const
 Return the Global VCS index of the kth species in the phase. More...
 
void setSpGlobalIndexVCS (const size_t spIndex, const size_t spGlobalIndex)
 set the Global VCS index of the kth species in the phase More...
 
void setTotalMolesInert (const double tMolesInert)
 Sets the total moles of inert in the phase. More...
 
double totalMolesInert () const
 Returns the value of the total kmol of inert in the phase. More...
 
size_t elemGlobalIndex (const size_t e) const
 Returns the global index of the local element index for the phase. More...
 
void setElemGlobalIndex (const size_t eLocal, const size_t eGlobal)
 sets a local phase element to a global index value More...
 
size_t nElemConstraints () const
 Returns the number of element constraints. More...
 
std::string elementName (const size_t e) const
 Name of the element constraint with index e. More...
 
int elementType (const size_t e) const
 Type of the element constraint with index e. More...
 
void setElementType (const size_t e, const int eType)
 Set the element Type of the element constraint with index e. More...
 
size_t transferElementsFM (const ThermoPhase *const tPhase)
 Transfer all of the element information from the ThermoPhase object to the vcs_VolPhase object. More...
 
const Array2DgetFormulaMatrix () const
 Get a constant form of the Species Formula Matrix. More...
 
int speciesUnknownType (const size_t k) const
 Returns the type of the species unknown. More...
 
int elementActive (const size_t e) const
 
size_t nSpecies () const
 Return the number of species in the phase. More...
 
std::string eos_name () const
 Return the name corresponding to the equation of state. More...
 

Public Attributes

size_t VP_ID_
 Original ID of the phase in the problem. More...
 
bool m_singleSpecies
 If true, this phase consists of a single species. More...
 
bool m_gasPhase
 If true, this phase is a gas-phase like phase. More...
 
int m_eqnState
 Type of the equation of state. More...
 
size_t ChargeNeutralityElement
 This is the element number for the charge neutrality condition of the phase. More...
 
int p_activityConvention
 Convention for the activity formulation. More...
 
std::string PhaseName
 String name for the phase. More...
 

Private Member Functions

void setMoleFractions (const double *const xmol)
 Set the mole fractions from a conventional mole fraction vector. More...
 
void _updateActCoeff () const
 Evaluate the activity coefficients at the current conditions. More...
 
void _updateGStar () const
 Gibbs free energy calculation for standard states. More...
 
void _updateG0 () const
 Gibbs free energy calculation at a temperature for the reference state of each species. More...
 
void _updateVolStar () const
 Molar volume calculation for standard states. More...
 
double _updateVolPM () const
 Calculate the partial molar volumes of all species and return the total volume. More...
 
void _updateLnActCoeffJac ()
 Evaluation of Activity Coefficient Jacobians. More...
 
void _updateMoleFractionDependencies ()
 Updates the mole fraction dependencies. More...
 

Private Attributes

VCS_SOLVEm_owningSolverObject
 Backtrack value of VCS_SOLVE *. More...
 
size_t m_numElemConstraints
 Number of element constraints within the problem. More...
 
std::vector< std::string > m_elementNames
 vector of strings containing the element constraint names More...
 
vector_int m_elementActive
 boolean indicating whether an element constraint is active for the current problem More...
 
vector_int m_elementType
 Type of the element constraint. More...
 
Array2D m_formulaMatrix
 Formula Matrix for the phase. More...
 
vector_int m_speciesUnknownType
 Type of the species unknown. More...
 
std::vector< size_t > m_elemGlobalIndex
 Index of the element number in the global list of elements stored in VCS_PROB or VCS_SOLVE. More...
 
size_t m_numSpecies
 Number of species in the phase. More...
 
double m_totalMolesInert
 Total moles of inert in the phase. More...
 
bool m_isIdealSoln
 Boolean indicating whether the phase is an ideal solution and therefore its molar-based activity coefficients are uniformly equal to one. More...
 
int m_existence
 Current state of existence: More...
 
int m_MFStartIndex
 
std::vector< size_t > IndSpecies
 Index into the species vectors. More...
 
std::vector< vcs_SpeciesProperties * > ListSpeciesPtr
 Vector of Species structures for the species belonging to this phase. More...
 
ThermoPhaseTP_ptr
 If we are using Cantera, this is the pointer to the ThermoPhase object. More...
 
double v_totalMoles
 Total mols in the phase. units are kmol. More...
 
vector_fp Xmol_
 Vector of the current mole fractions for species in the phase. More...
 
vector_fp creationMoleNumbers_
 Vector of current creationMoleNumbers_. More...
 
std::vector< size_t > creationGlobalRxnNumbers_
 Vector of creation global reaction numbers for the phase stability problem. More...
 
size_t m_phiVarIndex
 If the potential is a solution variable in VCS, it acts as a species. More...
 
double m_totalVol
 Total Volume of the phase. Units are m**3. More...
 
vector_fp SS0ChemicalPotential
 Vector of calculated SS0 chemical potentials for the current Temperature. More...
 
vector_fp StarChemicalPotential
 Vector of calculated Star chemical potentials for the current Temperature and pressure. More...
 
vector_fp StarMolarVol
 Vector of the Star molar Volumes of the species. units m3 / kmol. More...
 
vector_fp PartialMolarVol
 Vector of the Partial molar Volumes of the species. units m3 / kmol. More...
 
vector_fp ActCoeff
 Vector of calculated activity coefficients for the current state. More...
 
Array2D np_dLnActCoeffdMolNumber
 Vector of the derivatives of the ln activity coefficient wrt to the current mole number multiplied by the current phase moles. More...
 
int m_vcsStateStatus
 Status. More...
 
double m_phi
 Value of the potential for the phase (Volts) More...
 
bool m_UpToDate
 Boolean indicating whether the object has an up-to-date mole number vector and potential with respect to the current vcs state calc status. More...
 
bool m_UpToDate_AC
 Boolean indicating whether activity coefficients are up to date. More...
 
bool m_UpToDate_VolStar
 Boolean indicating whether Star volumes are up to date. More...
 
bool m_UpToDate_VolPM
 Boolean indicating whether partial molar volumes are up to date. More...
 
bool m_UpToDate_GStar
 Boolean indicating whether GStar is up to date. More...
 
bool m_UpToDate_G0
 Boolean indicating whether G0 is up to date. More...
 
double Temp_
 Current value of the temperature for this object, and underlying objects. More...
 
double Pres_
 Current value of the pressure for this object, and underlying objects. More...
 

Detailed Description

Phase information and Phase calculations for vcs.

Each phase in a vcs calculation has a vcs_VolPhase object associated with it. This object helps to coordinate property evaluations for species within the phase. Usually these evaluations must be carried out on a per phase basis. However, vcs frequently needs per species quantities. Therefore, we need an interface layer between vcs and Cantera's ThermoPhase.

The species stay in the same ordering within this structure. The vcs algorithm will change the ordering of species in the global species list. However, the indexing of species in this list stays the same. This structure contains structures that point to the species belonging to this phase in the global vcs species list.

This object is considered not to own the underlying Cantera ThermoPhase object for the phase.

This object contains an idea of the temperature and pressure. It checks to see if if the temperature and pressure has changed before calling underlying property evaluation routines.

The object contains values for the electric potential of a phase. It coordinates the evaluation of properties wrt when the electric potential of a phase has changed.

The object knows about the mole fractions of the phase. It controls the values of mole fractions, and coordinates the property evaluation wrt to changes in the mole fractions. It also will keep track of the likely values of mole fractions in multicomponent phases even when the phase doesn't actually exist within the thermo program.

The object knows about the total moles of a phase. It checks to see if the phase currently exists or not, and modifies its behavior accordingly.

Activity coefficients and volume calculations are lagged. They are only called when they are needed (and when the state has changed so that they need to be recalculated).

Definition at line 81 of file vcs_VolPhase.h.

Member Function Documentation

◆ resize()

void resize ( const size_t  phaseNum,
const size_t  numSpecies,
const size_t  numElem,
const char *const  phaseName,
const double  molesInert = 0.0 
)

The resize() function fills in all of the initial information if it is not given in the constructor.

Parameters
phaseNumindex of the phase in the vcs problem
numSpeciesNumber of species in the phase
phaseNameString name for the phase
molesInertkmoles of inert in the phase (defaults to zero)

Definition at line 182 of file vcs_VolPhase.cpp.

References AssertThrowMsg, Cantera::npos, VCS_SPECIES_TYPE_MOLNUM, and VCS_STATECALC_OLD.

Referenced by Cantera::vcs_Cantera_to_vprob().

◆ AC_calc_one()

double AC_calc_one ( size_t  kspec) const

Evaluate activity coefficients and return the kspec coefficient.

We carry out a calculation whenever m_UpToDate_AC is false. Specifically whenever a phase goes zero, we do not carry out calculations on it.

Parameters
kspecspecies number
Deprecated:
Unused. To be removed after Cantera 2.3.

Definition at line 288 of file vcs_VolPhase.cpp.

References Cantera::warn_deprecated().

Referenced by VCS_SPECIES_THERMO::eval_ac().

◆ setMoleFractionsState()

void setMoleFractionsState ( const double  molNum,
const double *const  moleFracVec,
const int  vcsStateStatus 
)

Set the moles and/or mole fractions within the phase.

Parameters
molNumtotal moles in the phase
moleFracVecVector of input mole fractions
vcsStateStatusStatus flag for this update

Definition at line 364 of file vcs_VolPhase.cpp.

References VCS_PHASE_EXIST_NO, VCS_PHASE_EXIST_YES, VCS_PHASE_EXIST_ZEROEDPHASE, and VCS_STATECALC_TMP.

Referenced by VCS_SOLVE::vcs_prob_update().

◆ setMolesFromVCS()

void setMolesFromVCS ( const int  stateCalc,
const double *  molesSpeciesVCS = 0 
)

Set the moles within the phase.

This function takes as input the mole numbers in vcs format, and then updates this object with their values. This is essentially a gather routine.

Parameters
molesSpeciesVCSArray of mole numbers. Note, the indices for species in this array may not be contiguous. IndSpecies[] is needed to gather the species into the local contiguous vector format.

Definition at line 414 of file vcs_VolPhase.cpp.

References AssertThrowMsg, Cantera::npos, VCS_PHASE_EXIST_ALWAYS, VCS_PHASE_EXIST_NO, VCS_PHASE_EXIST_YES, VCS_SPECIES_TYPE_INTERFACIALVOLTAGE, VCS_STATECALC_NEW, and VCS_STATECALC_OLD.

Referenced by VCS_SOLVE::vcs_CalcLnActCoeffJac(), Cantera::vcs_Cantera_update_vprob(), VCS_SOLVE::vcs_chemPotPhase(), and VCS_SOLVE::vcs_VolTotal().

◆ setMolesFromVCSCheck()

void setMolesFromVCSCheck ( const int  vcsStateStatus,
const double *  molesSpeciesVCS,
const double *const  TPhMoles 
)

Set the moles within the phase.

This function takes as input the mole numbers in vcs format, and then updates this object with their values. This is essentially a gather routine.

Additionally it checks to see that the total moles value in TPhMoles[iplace] is equal to the internally computed value. If this isn't the case, an error exit is carried out.

Parameters
vcsStateStatusState calc value either VCS_STATECALC_OLD or VCS_STATECALC_NEW. With any other value nothing is done.
molesSpeciesVCSarray of mole numbers. Note, the indices for species in this array may not be contiguous. IndSpecies[] is needed to gather the species into the local contiguous vector format.
TPhMolesVCS's array containing the number of moles in each phase.

Definition at line 495 of file vcs_VolPhase.cpp.

References Cantera::vcs_doubleEqual().

Referenced by VCS_SOLVE::vcs_delete_species(), VCS_SOLVE::vcs_reinsert_deleted(), and VCS_SOLVE::vcs_updateVP().

◆ updateFromVCS_MoleNumbers()

void updateFromVCS_MoleNumbers ( const int  stateCalc)

Update the moles within the phase, if necessary.

This function takes as input the stateCalc value, which determines where within VCS_SOLVE to fetch the mole numbers. It then updates this object with their values. This is essentially a gather routine.

Parameters
stateCalcState calc value either VCS_STATECALC_OLD or VCS_STATECALC_NEW. With any other value nothing is done.

Definition at line 513 of file vcs_VolPhase.cpp.

References VCS_STATECALC_NEW, and VCS_STATECALC_OLD.

Referenced by VCS_SOLVE::vcs_dfe().

◆ sendToVCS_ActCoeff()

void sendToVCS_ActCoeff ( const int  stateCalc,
double *const  AC 
)

Fill in an activity coefficients vector within a VCS_SOLVE object.

This routine will calculate the activity coefficients for the current phase, and fill in the corresponding entries in the VCS activity coefficients vector.

Parameters
ACvector of activity coefficients for all of the species in all of the phases in a VCS problem. Only the entries for the current phase are filled in.

Definition at line 521 of file vcs_VolPhase.cpp.

Referenced by VCS_SOLVE::vcs_add_all_deleted(), VCS_SOLVE::vcs_chemPotPhase(), VCS_SOLVE::vcs_dfe(), and VCS_SOLVE::vcs_phaseStabilityTest().

◆ setElectricPotential()

void setElectricPotential ( const double  phi)

set the electric potential of the phase

Parameters
phielectric potential (volts)

Definition at line 557 of file vcs_VolPhase.cpp.

Referenced by Cantera::vcs_Cantera_to_vprob(), Cantera::vcs_Cantera_update_vprob(), and VCS_SOLVE::vcs_prob_update().

◆ electricPotential()

double electricPotential ( ) const

Returns the electric field of the phase.

Units are potential

Definition at line 568 of file vcs_VolPhase.cpp.

Referenced by VCS_SOLVE::vcs_chemPotPhase(), VCS_SOLVE::vcs_dfe(), and VCS_SOLVE::vcs_prob_update().

◆ GStar_calc_one()

double GStar_calc_one ( size_t  kspec) const

Gibbs free energy calculation for standard state of one species.

Calculate the Gibbs free energies for the standard state of the kth species. The results are held internally within the object.

Parameters
kspecSpecies number (within the phase)
Returns
the Gibbs free energy for the standard state of the kth species.

Definition at line 318 of file vcs_VolPhase.cpp.

Referenced by VCS_SPECIES_THERMO::GStar_R_calc(), and VCS_PROB::prob_report().

◆ G0_calc_one()

double G0_calc_one ( size_t  kspec) const

Gibbs free energy calculation at a temperature for the reference state of a species, return a value for one species.

Parameters
kspecspecies index
Returns
value of the Gibbs free energy

Definition at line 304 of file vcs_VolPhase.cpp.

Referenced by VCS_PROB::prob_report().

◆ VolStar_calc_one()

double VolStar_calc_one ( size_t  kspec) const

Molar volume calculation for standard state of one species.

Calculate the molar volume for the standard states. The results are held internally within the object.

Parameters
kspecSpecies number (within the phase)
Returns
molar volume of the kspec species's standard state (m**3/kmol)

Definition at line 600 of file vcs_VolPhase.cpp.

Referenced by VCS_SPECIES_THERMO::VolStar_calc().

◆ sendToVCS_VolPM()

double sendToVCS_VolPM ( double *const  VolPM) const

Fill in the partial molar volume vector for VCS.

This routine will calculate the partial molar volumes for the current phase (if needed), and fill in the corresponding entries in the VCS partial molar volumes vector.

Parameters
VolPMvector of partial molar volumes for all of the species in all of the phases in a VCS problem. Only the entries for the current phase are filled in.

Definition at line 534 of file vcs_VolPhase.cpp.

Referenced by VCS_SOLVE::vcs_VolTotal().

◆ sendToVCS_GStar()

void sendToVCS_GStar ( double *const  gstar) const

Fill in the partial molar volume vector for VCS.

This routine will calculate the partial molar volumes for the current phase (if needed), and fill in the corresponding entries in the VCS partial molar volumes vector.

Parameters
VolPMvector of partial molar volumes for all of the species in all of the phases in a VCS problem. Only the entries for the current phase are filled in.
Todo:
This function's documentation is incorrect.

Definition at line 546 of file vcs_VolPhase.cpp.

Referenced by VCS_SOLVE::vcs_evalSS_TP().

◆ setState_TP()

void setState_TP ( const double  temperature_Kelvin,
const double  pressure_PA 
)

Sets the temperature and pressure in this object and underlying ThermoPhase objects.

Parameters
temperature_Kelvin(Kelvin)
pressure_PAPressure (MKS units - Pascal)

Definition at line 573 of file vcs_VolPhase.cpp.

Referenced by VCS_SPECIES_THERMO::GStar_R_calc(), VCS_PROB::prob_report(), Cantera::vcs_Cantera_update_vprob(), VCS_SOLVE::vcs_evalSS_TP(), VCS_SOLVE::vcs_VolTotal(), and VCS_SPECIES_THERMO::VolStar_calc().

◆ setState_T()

void setState_T ( const double  temperature_Kelvin)

Sets the temperature in this object and underlying ThermoPhase objects.

Parameters
temperature_Kelvin(Kelvin)

Definition at line 589 of file vcs_VolPhase.cpp.

◆ setPtrThermoPhase()

void setPtrThermoPhase ( ThermoPhase tp_ptr)

Set the pointer for Cantera's ThermoPhase parameter.

When we first initialize the ThermoPhase object, we read the state of the ThermoPhase into vcs_VolPhase object.

Parameters
tp_ptrPointer to the ThermoPhase object corresponding to this phase.

Definition at line 710 of file vcs_VolPhase.cpp.

References plogf, and Phase::temperature().

Referenced by Cantera::vcs_Cantera_to_vprob().

◆ ptrThermoPhase()

const ThermoPhase * ptrThermoPhase ( ) const

Return a const ThermoPhase pointer corresponding to this phase.

Returns
pointer to the ThermoPhase.

Definition at line 745 of file vcs_VolPhase.cpp.

◆ totalMoles()

double totalMoles ( ) const

◆ molefraction()

double molefraction ( size_t  kspec) const

Returns the mole fraction of the kspec species.

Parameters
kspecIndex of the species in the phase
Returns
Value of the mole fraction

Definition at line 755 of file vcs_VolPhase.cpp.

Referenced by VCS_SOLVE::vcs_prob_update().

◆ setTotalMoles()

void setTotalMoles ( const double  totalMols)

Sets the total moles in the phase.

We don't have to flag the internal state as changing here because we have just changed the total moles.

Parameters
totalMolsTotal moles in the phase (kmol)

Definition at line 775 of file vcs_VolPhase.cpp.

References AssertThrowMsg, VCS_PHASE_EXIST_ALWAYS, VCS_PHASE_EXIST_NO, and VCS_PHASE_EXIST_YES.

Referenced by Cantera::vcs_Cantera_to_vprob(), VCS_SOLVE::vcs_delete_multiphase(), VCS_SOLVE::vcs_prob_update(), and VCS_SOLVE::vcs_tmoles().

◆ setMolesOutOfDate()

void setMolesOutOfDate ( int  stateCalc = -1)

Sets the mole flag within the object to out of date.

This will trigger the object to go get the current mole numbers when it needs it.

Definition at line 797 of file vcs_VolPhase.cpp.

◆ setMolesCurrent()

void setMolesCurrent ( int  vcsStateStatus)

Sets the mole flag within the object to be current.

Definition at line 805 of file vcs_VolPhase.cpp.

◆ setMoleFractions()

void setMoleFractions ( const double *const  xmol)
private

Set the mole fractions from a conventional mole fraction vector.

Parameters
xmolValue of the mole fractions for the species in the phase. These are contiguous.

Definition at line 326 of file vcs_VolPhase.cpp.

◆ moleFractions()

const vector_fp & moleFractions ( ) const

Return a const reference to the mole fractions stored in the object.

Definition at line 354 of file vcs_VolPhase.cpp.

Referenced by VCS_SOLVE::vcs_prob_update().

◆ setCreationMoleNumbers()

void setCreationMoleNumbers ( const double *const  n_k,
const std::vector< size_t > &  creationGlobalRxnNumbers 
)

Sets the creationMoleNum's within the phase object.

Parameters
F_kPointer to a vector of n_k's

Definition at line 760 of file vcs_VolPhase.cpp.

◆ creationMoleNumbers()

const vector_fp & creationMoleNumbers ( std::vector< size_t > &  creationGlobalRxnNumbers) const

Return a const reference to the creationMoleNumbers stored in the object.

Returns
a const reference to the vector of creationMoleNumbers

Definition at line 769 of file vcs_VolPhase.cpp.

Referenced by VCS_SOLVE::vcs_phaseStabilityTest(), and VCS_SOLVE::vcs_popPhaseRxnStepSizes().

◆ isIdealSoln()

bool isIdealSoln ( ) const

Returns whether the phase is an ideal solution phase.

Definition at line 847 of file vcs_VolPhase.cpp.

Referenced by VCS_SOLVE::vcs_CalcLnActCoeffJac().

◆ phiVarIndex()

size_t phiVarIndex ( ) const

Return the index of the species that represents the the voltage of the phase.

Definition at line 852 of file vcs_VolPhase.cpp.

Referenced by Cantera::vcs_Cantera_update_vprob(), and VCS_SOLVE::vcs_prob_update().

◆ speciesProperty()

vcs_SpeciesProperties * speciesProperty ( const size_t  kindex)

Retrieve the kth Species structure for the species belonging to this phase.

The index into this vector is the species index within the phase.

Parameters
kindexkth species index.

Definition at line 866 of file vcs_VolPhase.cpp.

Referenced by VCS_SOLVE::vcs_prep_oneTime().

◆ exists()

int exists ( ) const

int indicating whether the phase exists or not

returns the m_existence int for the phase

  • VCS_PHASE_EXIST_ZEROEDPHASE = -6: Set to not exist by fiat from a higher level. This is used in phase stability boundary calculations
  • VCS_PHASE_EXIST_NO = 0: Doesn't exist currently
  • VCS_PHASE_EXIST_MINORCONC = 1: Exists, but the concentration is so low that an alternate method is used to calculate the total phase concentrations.
  • VCS_PHASE_EXIST_YES = 2 : Does exist currently
  • VCS_PHASE_EXIST_ALWAYS = 3: Always exists because it contains inerts which can't exist in any other phase. Or, the phase exists always because it consists of a single species, which is identified with the voltage, i.e., it's an electron metal phase.

Definition at line 871 of file vcs_VolPhase.cpp.

Referenced by VCS_SOLVE::recheck_deleted_phase(), VCS_SOLVE::vcs_delete_species(), VCS_SOLVE::vcs_phasePopDeterminePossibleList(), VCS_SOLVE::vcs_popPhaseID(), VCS_SOLVE::vcs_popPhasePossible(), VCS_SOLVE::vcs_popPhaseRxnStepSizes(), VCS_SOLVE::vcs_reinsert_deleted(), and VCS_SOLVE::vcs_RxnStepSizes().

◆ setExistence()

void setExistence ( const int  existence)

Set the existence flag in the object.

Note the total moles of the phase must have been set appropriately before calling this routine.

Parameters
existencePhase existence flag
Note
try to eliminate this routine

Definition at line 876 of file vcs_VolPhase.cpp.

References VCS_PHASE_EXIST_NO, and VCS_PHASE_EXIST_ZEROEDPHASE.

Referenced by Cantera::vcs_Cantera_update_vprob(), VCS_SOLVE::vcs_reinsert_deleted(), and VCS_SOLVE::vcs_SSPhase().

◆ spGlobalIndexVCS()

size_t spGlobalIndexVCS ( const size_t  spIndex) const

Return the Global VCS index of the kth species in the phase.

Parameters
spIndexlocal species index (0 to the number of species in the phase)
Returns
the VCS_SOLVE species index of the species. This changes as rearrangements are carried out.

Definition at line 896 of file vcs_VolPhase.cpp.

Referenced by VCS_PROB::prob_report(), VCS_SOLVE::recheck_deleted_phase(), Cantera::vcs_Cantera_update_vprob(), VCS_SOLVE::vcs_chemPotPhase(), VCS_SOLVE::vcs_deltag(), VCS_SOLVE::vcs_deltag_Phase(), VCS_SOLVE::vcs_phasePopDeterminePossibleList(), VCS_SOLVE::vcs_phaseStabilityTest(), VCS_SOLVE::vcs_popPhaseID(), VCS_SOLVE::vcs_popPhasePossible(), VCS_SOLVE::vcs_popPhaseRxnStepSizes(), VCS_SOLVE::vcs_prob_update(), and VCS_SOLVE::vcs_switch_pos().

◆ setSpGlobalIndexVCS()

void setSpGlobalIndexVCS ( const size_t  spIndex,
const size_t  spGlobalIndex 
)

set the Global VCS index of the kth species in the phase

Parameters
spIndexlocal species index (0 to the number of species in the phase)
Returns
the VCS_SOLVE species index of the that species This changes as rearrangements are carried out.

Definition at line 901 of file vcs_VolPhase.cpp.

Referenced by VCS_PROB::addOnePhaseSpecies(), and VCS_SOLVE::vcs_switch_pos().

◆ setTotalMolesInert()

void setTotalMolesInert ( const double  tMolesInert)

Sets the total moles of inert in the phase.

Parameters
tMolesInertValue of the total kmols of inert species in the phase.

Definition at line 910 of file vcs_VolPhase.cpp.

References VCS_PHASE_EXIST_ALWAYS, VCS_PHASE_EXIST_NO, and VCS_PHASE_EXIST_YES.

Referenced by VCS_SOLVE::vcs_prob_update(), and VCS_SOLVE::vcs_redim_TP().

◆ totalMolesInert()

double totalMolesInert ( ) const

Returns the value of the total kmol of inert in the phase.

Definition at line 935 of file vcs_VolPhase.cpp.

Referenced by VCS_PROB::prob_report(), Cantera::vcs_Cantera_update_vprob(), and VCS_SOLVE::vcs_prob_update().

◆ elemGlobalIndex()

size_t elemGlobalIndex ( const size_t  e) const

Returns the global index of the local element index for the phase.

Definition at line 940 of file vcs_VolPhase.cpp.

References AssertThrow.

Referenced by VCS_PROB::addOnePhaseSpecies(), and VCS_SOLVE::vcs_switch_elem_pos().

◆ setElemGlobalIndex()

void setElemGlobalIndex ( const size_t  eLocal,
const size_t  eGlobal 
)

sets a local phase element to a global index value

Parameters
eLocalLocal phase element index
eGlobalGlobal phase element index

Definition at line 946 of file vcs_VolPhase.cpp.

References AssertThrow.

Referenced by VCS_PROB::addPhaseElements(), and VCS_SOLVE::vcs_switch_elem_pos().

◆ nElemConstraints()

size_t nElemConstraints ( ) const

Returns the number of element constraints.

Definition at line 953 of file vcs_VolPhase.cpp.

Referenced by VCS_PROB::addOnePhaseSpecies(), VCS_PROB::addPhaseElements(), and VCS_SOLVE::vcs_switch_elem_pos().

◆ elementName()

std::string elementName ( const size_t  e) const

Name of the element constraint with index e.

Parameters
eElement index.

Definition at line 958 of file vcs_VolPhase.cpp.

Referenced by VCS_PROB::addPhaseElements().

◆ elementType()

int elementType ( const size_t  e) const

Type of the element constraint with index e.

Parameters
eElement index.

Definition at line 1091 of file vcs_VolPhase.cpp.

Referenced by VCS_PROB::addPhaseElements().

◆ setElementType()

void setElementType ( const size_t  e,
const int  eType 
)

Set the element Type of the element constraint with index e.

Parameters
eElement index
eTypetype of the element.
Deprecated:
Unused. To be removed after Cantera 2.3.

Definition at line 1096 of file vcs_VolPhase.cpp.

References Cantera::warn_deprecated().

◆ transferElementsFM()

size_t transferElementsFM ( const ThermoPhase *const  tPhase)

Transfer all of the element information from the ThermoPhase object to the vcs_VolPhase object.

Also decide whether we need a new charge neutrality element in the phase to enforce a charge neutrality constraint.

Parameters
tPhasePointer to the ThermoPhase object

Definition at line 987 of file vcs_VolPhase.cpp.

References Phase::charge(), Cantera::chargeNeutralityElement(), Phase::elementName(), Phase::elementType(), Cantera::hasChargedSpecies(), Phase::id(), Phase::nAtoms(), Phase::nElements(), Cantera::npos, Phase::nSpecies(), VCS_ELEM_TYPE_CHARGENEUTRALITY, VCS_ELEM_TYPE_ELECTRONCHARGE, VCS_SPECIES_TYPE_INTERFACIALVOLTAGE, and VCS_SPECIES_TYPE_MOLNUM.

◆ getFormulaMatrix()

const Array2D & getFormulaMatrix ( ) const

Get a constant form of the Species Formula Matrix.

Returns a double** pointer such that fm[e][f] is the formula matrix entry for element e for species k

Definition at line 1103 of file vcs_VolPhase.cpp.

Referenced by VCS_PROB::addOnePhaseSpecies().

◆ speciesUnknownType()

int speciesUnknownType ( const size_t  k) const

Returns the type of the species unknown.

Parameters
kspecies index
Returns
the SpeciesUnknownType[k] = type of species
  • Normal -> VCS_SPECIES_TYPE_MOLUNK (unknown is the mole number in the phase)
  • metal electron -> VCS_SPECIES_INTERFACIALVOLTAGE (unknown is the interfacial voltage (volts))

Definition at line 1108 of file vcs_VolPhase.cpp.

Referenced by VCS_SOLVE::vcs_prob_update().

◆ nSpecies()

size_t nSpecies ( ) const

◆ eos_name()

std::string eos_name ( ) const

Return the name corresponding to the equation of state.

Definition at line 1123 of file vcs_VolPhase.cpp.

Referenced by VCS_PROB::prob_report(), and Cantera::vcs_Cantera_update_vprob().

◆ _updateActCoeff()

void _updateActCoeff ( ) const
private

Evaluate the activity coefficients at the current conditions.

We carry out a calculation whenever m_UpToDate_AC is false. Specifically whenever a phase goes zero, we do not carry out calculations on it.

Definition at line 278 of file vcs_VolPhase.cpp.

◆ _updateGStar()

void _updateGStar ( ) const
private

Gibbs free energy calculation for standard states.

Calculate the Gibbs free energies for the standard states The results are held internally within the object.

Definition at line 312 of file vcs_VolPhase.cpp.

◆ _updateG0()

void _updateG0 ( ) const
private

Gibbs free energy calculation at a temperature for the reference state of each species.

Definition at line 298 of file vcs_VolPhase.cpp.

◆ _updateVolStar()

void _updateVolStar ( ) const
private

Molar volume calculation for standard states.

Calculate the molar volume for the standard states. The results are held internally within the object. Units are in m**3/kmol.

Definition at line 594 of file vcs_VolPhase.cpp.

◆ _updateVolPM()

double _updateVolPM ( ) const
private

Calculate the partial molar volumes of all species and return the total volume.

Calculates these quantities internally and then stores them

Returns
total volume [m^3]

Definition at line 608 of file vcs_VolPhase.cpp.

◆ _updateLnActCoeffJac()

void _updateLnActCoeffJac ( )
private

Evaluation of Activity Coefficient Jacobians.

This is the derivative of the ln of the activity coefficient with respect to mole number of jth species. (temp, pressure, and other mole numbers held constant)

We employ a finite difference derivative approach here. Because we have to change the mole numbers, this is not a const function, even though the paradigm would say that it should be.

Definition at line 629 of file vcs_VolPhase.cpp.

◆ _updateMoleFractionDependencies()

void _updateMoleFractionDependencies ( )
private

Updates the mole fraction dependencies.

Whenever the mole fractions change, this routine should be called.

Definition at line 343 of file vcs_VolPhase.cpp.

Member Data Documentation

◆ m_owningSolverObject

VCS_SOLVE* m_owningSolverObject
private

Backtrack value of VCS_SOLVE *.

Note the default for this is 0. That's a valid value too, since VCS_PROB also uses vcs_VolPhase objects.

Definition at line 555 of file vcs_VolPhase.h.

◆ VP_ID_

size_t VP_ID_

Original ID of the phase in the problem.

If a non-ideal phase splits into two due to a miscibility gap, these numbers will stay the same after the split.

Definition at line 563 of file vcs_VolPhase.h.

Referenced by VCS_PROB::prob_report(), Cantera::vcs_Cantera_update_vprob(), and VCS_SOLVE::vcs_prob_specify().

◆ m_singleSpecies

bool m_singleSpecies

◆ m_gasPhase

bool m_gasPhase

If true, this phase is a gas-phase like phase.

A RTlog(p/1atm) term is added onto the chemical potential for inert species if this is true.

Definition at line 573 of file vcs_VolPhase.h.

Referenced by VCS_PROB::prob_report(), Cantera::vcs_Cantera_to_vprob(), Cantera::vcs_Cantera_update_vprob(), VCS_SOLVE::vcs_GibbsPhase(), and VCS_SOLVE::vcs_Total_Gibbs().

◆ m_eqnState

int m_eqnState

Type of the equation of state.

The known types are listed at the top of this file.

Definition at line 579 of file vcs_VolPhase.h.

Referenced by Cantera::vcs_Cantera_to_vprob().

◆ ChargeNeutralityElement

size_t ChargeNeutralityElement

This is the element number for the charge neutrality condition of the phase.

If it has one. If it does not have a charge neutrality constraint, then this value is equal to -1

Definition at line 587 of file vcs_VolPhase.h.

◆ p_activityConvention

int p_activityConvention

Convention for the activity formulation.

  • 0 = molar based activities (default)
  • 1 = Molality based activities, mu = mu_0 + ln a_molality. Standard state is based on unity molality

Definition at line 595 of file vcs_VolPhase.h.

Referenced by Cantera::vcs_Cantera_to_vprob().

◆ m_numElemConstraints

size_t m_numElemConstraints
private

Number of element constraints within the problem.

This is usually equal to the number of elements.

Definition at line 602 of file vcs_VolPhase.h.

◆ m_elementNames

std::vector<std::string> m_elementNames
private

vector of strings containing the element constraint names

Length = nElemConstraints

Definition at line 608 of file vcs_VolPhase.h.

◆ m_elementActive

vector_int m_elementActive
private

boolean indicating whether an element constraint is active for the current problem

Definition at line 612 of file vcs_VolPhase.h.

◆ m_elementType

vector_int m_elementType
private

Type of the element constraint.

m_elType[j] = type of the element:

  • 0 VCS_ELEM_TYPE_ABSPOS Normal element that is positive or zero in all species.
  • 1 VCS_ELEM_TYPE_ELECTRONCHARGE element dof that corresponds to the charge DOF.
  • 2 VCS_ELEM_TYPE_OTHERCONSTRAINT Other constraint which may mean that a species has neg 0 or pos value of that constraint (other than charge)

Definition at line 625 of file vcs_VolPhase.h.

◆ m_formulaMatrix

Array2D m_formulaMatrix
private

Formula Matrix for the phase.

FormulaMatrix(kspec,j) = Formula Matrix for the species Number of elements, j, in the kspec species

Definition at line 632 of file vcs_VolPhase.h.

◆ m_speciesUnknownType

vector_int m_speciesUnknownType
private

Type of the species unknown.

SpeciesUnknownType[k] = type of species

  • Normal -> VCS_SPECIES_TYPE_MOLUNK. (unknown is the mole number in the phase)
  • metal electron -> VCS_SPECIES_INTERFACIALVOLTAGE. (unknown is the interfacial voltage (volts))

Definition at line 642 of file vcs_VolPhase.h.

◆ m_elemGlobalIndex

std::vector<size_t> m_elemGlobalIndex
private

Index of the element number in the global list of elements stored in VCS_PROB or VCS_SOLVE.

Definition at line 646 of file vcs_VolPhase.h.

◆ m_numSpecies

size_t m_numSpecies
private

Number of species in the phase.

Definition at line 649 of file vcs_VolPhase.h.

◆ PhaseName

std::string PhaseName

◆ m_totalMolesInert

double m_totalMolesInert
private

Total moles of inert in the phase.

Definition at line 657 of file vcs_VolPhase.h.

◆ m_isIdealSoln

bool m_isIdealSoln
private

Boolean indicating whether the phase is an ideal solution and therefore its molar-based activity coefficients are uniformly equal to one.

Definition at line 662 of file vcs_VolPhase.h.

◆ m_existence

int m_existence
private

Current state of existence:

  • VCS_PHASE_EXIST_ZEROEDPHASE = -6: Set to not exist by fiat from a higher level. This is used in phase stability boundary calculations
  • VCS_PHASE_EXIST_NO = 0: Doesn't exist currently
  • VCS_PHASE_EXIST_MINORCONC = 1: Exists, but the concentration is so low that an alternate method is used to calculate the total phase concentrations.
  • VCS_PHASE_EXIST_YES = 2 : Does exist currently
  • VCS_PHASE_EXIST_ALWAYS = 3: Always exists because it contains inerts which can't exist in any other phase. Or, the phase exists always because it consists of a single species, which is identified with the voltage, i.e., its an electron metal phase.

Definition at line 678 of file vcs_VolPhase.h.

◆ m_MFStartIndex

int m_MFStartIndex
private

This is always equal to zero. Am anticipating the case where the phase potential is species # 0, for multiphase phases. Right now we have the phase potential equal to 0 for single species phases, where we set by hand the mole fraction of species 0 to one.

Definition at line 688 of file vcs_VolPhase.h.

◆ IndSpecies

std::vector<size_t> IndSpecies
private

Index into the species vectors.

Maps the phase species number into the global species number. Note, as part of the vcs algorithm, the order of the species vector is changed during the algorithm

Definition at line 696 of file vcs_VolPhase.h.

◆ ListSpeciesPtr

std::vector<vcs_SpeciesProperties*> ListSpeciesPtr
private

Vector of Species structures for the species belonging to this phase.

The index into this vector is the species index within the phase.

Definition at line 702 of file vcs_VolPhase.h.

◆ TP_ptr

ThermoPhase* TP_ptr
private

If we are using Cantera, this is the pointer to the ThermoPhase object.

If not, this is null.

Definition at line 708 of file vcs_VolPhase.h.

◆ v_totalMoles

double v_totalMoles
private

Total mols in the phase. units are kmol.

Definition at line 711 of file vcs_VolPhase.h.

◆ Xmol_

vector_fp Xmol_
private

Vector of the current mole fractions for species in the phase.

Definition at line 714 of file vcs_VolPhase.h.

◆ creationMoleNumbers_

vector_fp creationMoleNumbers_
private

Vector of current creationMoleNumbers_.

These are the actual unknowns in the phase stability problem

Definition at line 720 of file vcs_VolPhase.h.

◆ creationGlobalRxnNumbers_

std::vector<size_t> creationGlobalRxnNumbers_
private

Vector of creation global reaction numbers for the phase stability problem.

The phase stability problem requires a global reaction number for each species in the phase. Usually this is the krxn = kglob - M for species in the phase that are not components. For component species, the choice of the reaction is one which maximizes the chance that the phase pops into (or remains in) existence.

The index here is the local phase species index. the value of the variable is the global vcs reaction number. Note, that the global reaction number will go out of order when the species positions are swapped. So, this number has to be recalculated.

Length = number of species in phase

Definition at line 737 of file vcs_VolPhase.h.

◆ m_phiVarIndex

size_t m_phiVarIndex
private

If the potential is a solution variable in VCS, it acts as a species.

This is the species index in the phase for the potential

Definition at line 741 of file vcs_VolPhase.h.

◆ m_totalVol

double m_totalVol
mutableprivate

Total Volume of the phase. Units are m**3.

Definition at line 744 of file vcs_VolPhase.h.

◆ SS0ChemicalPotential

vector_fp SS0ChemicalPotential
mutableprivate

Vector of calculated SS0 chemical potentials for the current Temperature.

Note, This is the chemical potential derived strictly from the polynomial in temperature. Pressure effects have to be added in to get to the standard state. Units are J/kmol.

Definition at line 753 of file vcs_VolPhase.h.

◆ StarChemicalPotential

vector_fp StarChemicalPotential
mutableprivate

Vector of calculated Star chemical potentials for the current Temperature and pressure.

Note, This is the chemical potential at unit activity. Thus, we can call it the standard state chemical potential as well. Units are J/kmol.

Definition at line 761 of file vcs_VolPhase.h.

◆ StarMolarVol

vector_fp StarMolarVol
mutableprivate

Vector of the Star molar Volumes of the species. units m3 / kmol.

Definition at line 764 of file vcs_VolPhase.h.

◆ PartialMolarVol

vector_fp PartialMolarVol
mutableprivate

Vector of the Partial molar Volumes of the species. units m3 / kmol.

Definition at line 767 of file vcs_VolPhase.h.

◆ ActCoeff

vector_fp ActCoeff
mutableprivate

Vector of calculated activity coefficients for the current state.

Whether or not this vector is current is determined by the bool m_UpToDate_AC.

Definition at line 774 of file vcs_VolPhase.h.

◆ np_dLnActCoeffdMolNumber

Array2D np_dLnActCoeffdMolNumber
mutableprivate

Vector of the derivatives of the ln activity coefficient wrt to the current mole number multiplied by the current phase moles.

np_dLnActCoeffdMolNumber(k,j);

  • j = id of the species mole number
  • k = id of the species activity coefficient

Definition at line 783 of file vcs_VolPhase.h.

◆ m_vcsStateStatus

int m_vcsStateStatus
private

Status.

valid values are

  • VCS_STATECALC_OLD
  • VCS_STATECALC_NEW
  • VCS_STATECALC_TMP

Definition at line 792 of file vcs_VolPhase.h.

◆ m_phi

double m_phi
private

Value of the potential for the phase (Volts)

Definition at line 795 of file vcs_VolPhase.h.

◆ m_UpToDate

bool m_UpToDate
private

Boolean indicating whether the object has an up-to-date mole number vector and potential with respect to the current vcs state calc status.

Definition at line 799 of file vcs_VolPhase.h.

◆ m_UpToDate_AC

bool m_UpToDate_AC
mutableprivate

Boolean indicating whether activity coefficients are up to date.

Activity coefficients and volume calculations are lagged. They are only called when they are needed (and when the state has changed so that they need to be recalculated).

Definition at line 807 of file vcs_VolPhase.h.

◆ m_UpToDate_VolStar

bool m_UpToDate_VolStar
mutableprivate

Boolean indicating whether Star volumes are up to date.

Activity coefficients and volume calculations are lagged. They are only called when they are needed (and when the state has changed so that they need to be recalculated). Star volumes are sensitive to temperature and pressure

Definition at line 816 of file vcs_VolPhase.h.

◆ m_UpToDate_VolPM

bool m_UpToDate_VolPM
mutableprivate

Boolean indicating whether partial molar volumes are up to date.

Activity coefficients and volume calculations are lagged. They are only called when they are needed (and when the state has changed so that they need to be recalculated). partial molar volumes are sensitive to everything

Definition at line 825 of file vcs_VolPhase.h.

◆ m_UpToDate_GStar

bool m_UpToDate_GStar
mutableprivate

Boolean indicating whether GStar is up to date.

GStar is sensitive to the temperature and the pressure, only

Definition at line 831 of file vcs_VolPhase.h.

◆ m_UpToDate_G0

bool m_UpToDate_G0
mutableprivate

Boolean indicating whether G0 is up to date.

G0 is sensitive to the temperature and the pressure, only

Definition at line 837 of file vcs_VolPhase.h.

◆ Temp_

double Temp_
private

Current value of the temperature for this object, and underlying objects.

Definition at line 840 of file vcs_VolPhase.h.

◆ Pres_

double Pres_
private

Current value of the pressure for this object, and underlying objects.

Definition at line 843 of file vcs_VolPhase.h.


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