Cantera  2.5.1
Public Member Functions | Private Member Functions | Private Attributes | List of all members
MultiPhase Class Reference

A class for multiphase mixtures. More...

#include <MultiPhase.h>

Collaboration diagram for MultiPhase:
[legend]

Public Member Functions

 MultiPhase ()
 Constructor. More...
 
virtual ~MultiPhase ()
 Destructor. More...
 
void addPhases (std::vector< ThermoPhase * > &phases, const vector_fp &phaseMoles)
 Add a vector of phases to the mixture. More...
 
void addPhases (MultiPhase &mix)
 Add all phases present in 'mix' to this mixture. More...
 
void addPhase (ThermoPhase *p, doublereal moles)
 Add a phase to the mixture. More...
 
size_t nElements () const
 Number of elements. More...
 
void checkElementIndex (size_t m) const
 Check that the specified element index is in range. More...
 
void checkElementArraySize (size_t mm) const
 Check that an array size is at least nElements(). More...
 
std::string elementName (size_t m) const
 Returns the name of the global element m. More...
 
size_t elementIndex (const std::string &name) const
 Returns the index of the element with name name. More...
 
size_t nSpecies () const
 Number of species, summed over all phases. More...
 
void checkSpeciesIndex (size_t k) const
 Check that the specified species index is in range. More...
 
void checkSpeciesArraySize (size_t kk) const
 Check that an array size is at least nSpecies(). More...
 
std::string speciesName (const size_t kGlob) const
 Name of species with global index kGlob. More...
 
doublereal nAtoms (const size_t kGlob, const size_t mGlob) const
 Returns the Number of atoms of global element mGlob in global species kGlob. More...
 
void getMoleFractions (doublereal *const x) const
 Returns the global Species mole fractions. More...
 
void init ()
 Process phases and build atomic composition array. More...
 
std::string phaseName (const size_t iph) const
 Returns the name of the n'th phase. More...
 
int phaseIndex (const std::string &pName) const
 Returns the index, given the phase name. More...
 
doublereal phaseMoles (const size_t n) const
 Return the number of moles in phase n. More...
 
void setPhaseMoles (const size_t n, const doublereal moles)
 Set the number of moles of phase with index n. More...
 
thermo_tphase (size_t n)
 Return a reference to phase n. More...
 
void checkPhaseIndex (size_t m) const
 Check that the specified phase index is in range Throws an exception if m is greater than nPhases() More...
 
void checkPhaseArraySize (size_t mm) const
 Check that an array size is at least nPhases() Throws an exception if mm is less than nPhases(). More...
 
doublereal speciesMoles (size_t kGlob) const
 Returns the moles of global species k. units = kmol. More...
 
size_t speciesIndex (size_t k, size_t p) const
 Return the global index of the species belonging to phase number p with local index k within the phase. More...
 
size_t speciesIndex (const std::string &speciesName, const std::string &phaseName)
 Return the global index of the species belonging to phase name phaseName with species name speciesName. More...
 
doublereal minTemp () const
 Minimum temperature for which all solution phases have valid thermo data. More...
 
doublereal maxTemp () const
 Maximum temperature for which all solution phases have valid thermo data. More...
 
doublereal charge () const
 Total charge summed over all phases (Coulombs). More...
 
doublereal phaseCharge (size_t p) const
 Charge (Coulombs) of phase with index p. More...
 
doublereal elementMoles (size_t m) const
 Total moles of global element m, summed over all phases. More...
 
void getChemPotentials (doublereal *mu) const
 Returns a vector of Chemical potentials. More...
 
void getValidChemPotentials (doublereal not_mu, doublereal *mu, bool standard=false) const
 Returns a vector of Valid chemical potentials. More...
 
doublereal temperature () const
 Temperature [K]. More...
 
void equilibrate (const std::string &XY, const std::string &solver="auto", double rtol=1e-9, int max_steps=50000, int max_iter=100, int estimate_equil=0, int log_level=0)
 Equilibrate a MultiPhase object. More...
 
void setTemperature (const doublereal T)
 Set the temperature [K]. More...
 
void setState_TP (const doublereal T, const doublereal Pres)
 Set the state of the underlying ThermoPhase objects in one call. More...
 
void setState_TPMoles (const doublereal T, const doublereal Pres, const doublereal *Moles)
 Set the state of the underlying ThermoPhase objects in one call. More...
 
doublereal pressure () const
 Pressure [Pa]. More...
 
doublereal volume () const
 The total mixture volume [m^3]. More...
 
void setPressure (doublereal P)
 Set the pressure [Pa]. More...
 
doublereal enthalpy () const
 The enthalpy of the mixture [J]. More...
 
doublereal IntEnergy () const
 The internal energy of the mixture [J]. More...
 
doublereal entropy () const
 The entropy of the mixture [J/K]. More...
 
doublereal gibbs () const
 The Gibbs function of the mixture [J]. More...
 
doublereal cp () const
 Heat capacity at constant pressure [J/K]. More...
 
size_t nPhases () const
 Number of phases. More...
 
bool solutionSpecies (size_t kGlob) const
 Return true is species kGlob is a species in a multicomponent solution phase. More...
 
size_t speciesPhaseIndex (const size_t kGlob) const
 Returns the phase index of the Kth "global" species. More...
 
doublereal moleFraction (const size_t kGlob) const
 Returns the mole fraction of global species k. More...
 
void setPhaseMoleFractions (const size_t n, const doublereal *const x)
 Set the Mole fractions of the nth phase. More...
 
void setMolesByName (const compositionMap &xMap)
 Set the number of moles of species in the mixture. More...
 
void setMolesByName (const std::string &x)
 Set the moles via a string containing their names. More...
 
void getMoles (doublereal *molNum) const
 Get the mole numbers of all species in the multiphase object. More...
 
void setMoles (const doublereal *n)
 Sets all of the global species mole numbers. More...
 
void addSpeciesMoles (const int indexS, const doublereal addedMoles)
 Adds moles of a certain species to the mixture. More...
 
void getElemAbundances (doublereal *elemAbundances) const
 Retrieves a vector of element abundances. More...
 
bool tempOK (size_t p) const
 Return true if the phase p has valid thermo data for the current temperature. More...
 
void uploadMoleFractionsFromPhases ()
 Update the locally-stored composition within this object to match the current compositions of the phase objects. More...
 
void updatePhases () const
 Set the states of the phase objects to the locally-stored state within this MultiPhase object. More...
 

Private Member Functions

void calcElemAbundances () const
 Calculate the element abundance vector. More...
 
double equilibrate_MultiPhaseEquil (int XY, doublereal err, int maxsteps, int maxiter, int loglevel)
 Set the mixture to a state of chemical equilibrium using the MultiPhaseEquil solver. More...
 

Private Attributes

vector_fp m_moles
 Vector of the number of moles in each phase. More...
 
std::vector< ThermoPhase * > m_phase
 Vector of the ThermoPhase pointers. More...
 
DenseMatrix m_atoms
 Global Stoichiometric Coefficient array. More...
 
vector_fp m_moleFractions
 Locally stored vector of mole fractions of all species comprising the MultiPhase object. More...
 
std::vector< size_t > m_spphase
 Mapping between the global species number and the phase ID. More...
 
std::vector< size_t > m_spstart
 Vector of ints containing of first species index in the global list of species for each phase. More...
 
std::vector< std::string > m_enames
 String names of the global elements. More...
 
vector_int m_atomicNumber
 Atomic number of each global element. More...
 
std::vector< std::string > m_snames
 Vector of species names in the problem. More...
 
std::map< std::string, size_t > m_enamemap
 Returns the global element index, given the element string name. More...
 
doublereal m_temp
 Current value of the temperature (kelvin) More...
 
doublereal m_press
 Current value of the pressure (Pa) More...
 
size_t m_nel
 Number of distinct elements in all of the phases. More...
 
size_t m_nsp
 Number of distinct species in all of the phases. More...
 
bool m_init
 True if the init() routine has been called, and the MultiPhase frozen. More...
 
size_t m_eloc
 Global ID of the element corresponding to the electronic charge. More...
 
std::vector< bool > m_temp_OK
 Vector of bools indicating whether temperatures are ok for phases. More...
 
doublereal m_Tmin
 Minimum temperature for which thermo parameterizations are valid. More...
 
doublereal m_Tmax
 Minimum temperature for which thermo parameterizations are valid. More...
 
vector_fp m_elemAbundances
 Vector of element abundances. More...
 

Detailed Description

A class for multiphase mixtures.

The mixture can contain any number of phases of any type.

This object is the basic tool used by Cantera for use in Multiphase equilibrium calculations.

It is a container for a set of phases. Each phase has a given number of kmoles. Therefore, MultiPhase may be considered an "extrinsic" thermodynamic object, in contrast to the ThermoPhase object, which is an "intrinsic" thermodynamic object.

MultiPhase may be considered to be "upstream" of the ThermoPhase objects in the sense that setting a property within MultiPhase, such as temperature, pressure, or species mole number, affects the underlying ThermoPhase object, but not the other way around.

All phases have the same temperature and pressure, and a specified number of moles for each phase. The phases do not need to have the same elements. For example, a mixture might consist of a gaseous phase with elements (H, C, O, N), a solid carbon phase containing only element C, etc. A master element set will be constructed for the mixture that is the intersection of the elements of each phase.

Below, reference is made to global species and global elements. These refer to the collective species and elements encompassing all of the phases tracked by the object.

The global element list kept by this object is an intersection of the element lists of all the phases that comprise the MultiPhase.

The global species list kept by this object is a concatenated list of all of the species in all the phases that comprise the MultiPhase. The ordering of species is contiguous with respect to the phase id.

Definition at line 57 of file MultiPhase.h.

Constructor & Destructor Documentation

◆ MultiPhase()

Constructor.

The constructor takes no arguments, since phases are added using method addPhase().

Definition at line 21 of file MultiPhase.cpp.

◆ ~MultiPhase()

virtual ~MultiPhase ( )
inlinevirtual

Destructor.

Does nothing. Class MultiPhase does not take "ownership" (i.e. responsibility for destroying) the phase objects.

Definition at line 69 of file MultiPhase.h.

Member Function Documentation

◆ addPhases() [1/2]

void addPhases ( std::vector< ThermoPhase * > &  phases,
const vector_fp phaseMoles 
)

Add a vector of phases to the mixture.

See the single addPhases command. This just does a bunch of phases at one time

Parameters
phasesVector of pointers to phases
phaseMolesVector of mole numbers in each phase (kmol)

Definition at line 40 of file MultiPhase.cpp.

References MultiPhase::addPhase(), MultiPhase::init(), and MultiPhase::phaseMoles().

◆ addPhases() [2/2]

void addPhases ( MultiPhase mix)

Add all phases present in 'mix' to this mixture.

Parameters
mixAdd all of the phases in another MultiPhase object to the current object.

Definition at line 33 of file MultiPhase.cpp.

References MultiPhase::addPhase(), MultiPhase::m_moles, MultiPhase::m_phase, and MultiPhase::nPhases().

◆ addPhase()

void addPhase ( ThermoPhase p,
doublereal  moles 
)

◆ nElements()

size_t nElements ( ) const
inline

Number of elements.

Definition at line 98 of file MultiPhase.h.

References MultiPhase::m_nel.

Referenced by Cantera::BasisOptimize(), Cantera::ElemRearrange(), and MultiPhaseEquil::MultiPhaseEquil().

◆ checkElementIndex()

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

Definition at line 716 of file MultiPhase.cpp.

References MultiPhase::m_nel.

◆ checkElementArraySize()

void checkElementArraySize ( size_t  mm) const

Check that an array size is at least nElements().

Throws an exception if mm is less than nElements(). Used before calls which take an array pointer.

Definition at line 723 of file MultiPhase.cpp.

References MultiPhase::m_nel.

◆ elementName()

std::string elementName ( size_t  m) const

Returns the name of the global element m.

Parameters
mindex of the global element

Definition at line 730 of file MultiPhase.cpp.

References MultiPhase::m_enames.

Referenced by MultiPhaseEquil::MultiPhaseEquil().

◆ elementIndex()

size_t elementIndex ( const std::string &  name) const

Returns the index of the element with name name.

Parameters
nameString name of the global element

Definition at line 735 of file MultiPhase.cpp.

References MultiPhase::m_enames, MultiPhase::m_nel, and Cantera::npos.

◆ nSpecies()

size_t nSpecies ( ) const
inline

◆ checkSpeciesIndex()

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

Definition at line 745 of file MultiPhase.cpp.

References MultiPhase::m_nsp.

◆ checkSpeciesArraySize()

void checkSpeciesArraySize ( size_t  kk) const

Check that an array size is at least nSpecies().

Throws an exception if kk is less than nSpecies(). Used before calls which take an array pointer.

Definition at line 752 of file MultiPhase.cpp.

References MultiPhase::m_nsp.

◆ speciesName()

std::string speciesName ( const size_t  kGlob) const

Name of species with global index kGlob.

Parameters
kGlobglobal species index

Definition at line 759 of file MultiPhase.cpp.

References MultiPhase::m_snames.

Referenced by vcs_MultiPhaseEquil::equilibrate_TP(), MultiPhaseEquil::MultiPhaseEquil(), MultiPhase::setMolesByName(), and MultiPhase::speciesIndex().

◆ nAtoms()

doublereal nAtoms ( const size_t  kGlob,
const size_t  mGlob 
) const

Returns the Number of atoms of global element mGlob in global species kGlob.

Parameters
kGlobglobal species index
mGlobglobal element index
Returns
the number of atoms.

Definition at line 764 of file MultiPhase.cpp.

References MultiPhase::m_atoms.

Referenced by MultiPhaseEquil::getComponents(), and MultiPhaseEquil::MultiPhaseEquil().

◆ getMoleFractions()

void getMoleFractions ( doublereal *const  x) const

Returns the global Species mole fractions.

Write the array of species mole fractions into array x. The mole fractions are normalized to sum to one in each phase.

Parameters
xvector of mole fractions. Length = number of global species.

Definition at line 769 of file MultiPhase.cpp.

References MultiPhase::m_moleFractions.

◆ init()

void init ( )

◆ phaseName()

std::string phaseName ( const size_t  iph) const

Returns the name of the n'th phase.

Parameters
iphphase Index

Definition at line 774 of file MultiPhase.cpp.

References Phase::id(), and MultiPhase::m_phase.

Referenced by MultiPhase::speciesIndex().

◆ phaseIndex()

int phaseIndex ( const std::string &  pName) const

Returns the index, given the phase name.

Parameters
pNameName of the phase
Returns
the index. A value of -1 means the phase isn't in the object.

Definition at line 780 of file MultiPhase.cpp.

References MultiPhase::m_phase, and MultiPhase::nPhases().

Referenced by MultiPhase::speciesIndex().

◆ phaseMoles()

doublereal phaseMoles ( const size_t  n) const

Return the number of moles in phase n.

Parameters
nIndex of the phase.

Definition at line 790 of file MultiPhase.cpp.

References MultiPhase::m_moles.

Referenced by MultiPhase::addPhases(), MultiPhaseEquil::computeReactionSteps(), vcs_MultiPhaseEquil::equilibrate_HP(), vcs_MultiPhaseEquil::equilibrate_SP(), and Cantera::operator<<().

◆ setPhaseMoles()

void setPhaseMoles ( const size_t  n,
const doublereal  moles 
)

Set the number of moles of phase with index n.

Parameters
nIndex of the phase
molesNumber of moles in the phase (kmol)

Definition at line 795 of file MultiPhase.cpp.

References MultiPhase::m_moles.

Referenced by VCS_SOLVE::vcs_prob_update().

◆ phase()

ThermoPhase & phase ( size_t  n)

Return a reference to phase n.

The state of phase n is also updated to match the state stored locally in the mixture object.

Parameters
nPhase Index
Returns
Reference to the ThermoPhase object for the phase

Definition at line 159 of file MultiPhase.cpp.

References MultiPhase::init(), MultiPhase::m_init, MultiPhase::m_moleFractions, MultiPhase::m_phase, MultiPhase::m_press, MultiPhase::m_spstart, and MultiPhase::m_temp.

Referenced by MultiPhaseEquil::computeReactionSteps(), Cantera::operator<<(), vcs_MultiPhaseEquil::reportCSV(), VCS_SOLVE::vcs_prob_specifyFully(), and VCS_SOLVE::VCS_SOLVE().

◆ checkPhaseIndex()

void checkPhaseIndex ( size_t  m) const

Check that the specified phase index is in range Throws an exception if m is greater than nPhases()

Definition at line 170 of file MultiPhase.cpp.

References MultiPhase::nPhases().

◆ checkPhaseArraySize()

void checkPhaseArraySize ( size_t  mm) const

Check that an array size is at least nPhases() Throws an exception if mm is less than nPhases().

Used before calls which take an array pointer.

Definition at line 177 of file MultiPhase.cpp.

References MultiPhase::nPhases().

◆ speciesMoles()

doublereal speciesMoles ( size_t  kGlob) const

Returns the moles of global species k. units = kmol.

Parameters
kGlobGlobal species index k

Definition at line 184 of file MultiPhase.cpp.

References MultiPhase::m_moleFractions, MultiPhase::m_moles, and MultiPhase::m_spphase.

Referenced by MultiPhaseEquil::computeReactionSteps(), vcs_MultiPhaseEquil::equilibrate_TP(), and MultiPhaseEquil::MultiPhaseEquil().

◆ speciesIndex() [1/2]

size_t speciesIndex ( size_t  k,
size_t  p 
) const
inline

Return the global index of the species belonging to phase number p with local index k within the phase.

Parameters
klocal index of the species within the phase
pindex of the phase

Definition at line 227 of file MultiPhase.h.

References MultiPhase::m_spstart.

Referenced by MultiPhase::elementMoles(), and MultiPhase::phaseCharge().

◆ speciesIndex() [2/2]

size_t speciesIndex ( const std::string &  speciesName,
const std::string &  phaseName 
)

Return the global index of the species belonging to phase name phaseName with species name speciesName.

Parameters
speciesNameSpecies Name
phaseNamePhase Name
Returns
the global index

If the species or phase name is not recognized, this routine throws a CanteraError.

Definition at line 214 of file MultiPhase.cpp.

References MultiPhase::init(), MultiPhase::m_init, MultiPhase::m_phase, MultiPhase::m_spstart, Cantera::npos, MultiPhase::phaseIndex(), MultiPhase::phaseName(), and MultiPhase::speciesName().

◆ minTemp()

doublereal minTemp ( ) const
inline

Minimum temperature for which all solution phases have valid thermo data.

Stoichiometric phases are not considered, since they may have thermo data only valid for conditions for which they are stable.

Definition at line 247 of file MultiPhase.h.

References MultiPhase::m_Tmin.

Referenced by vcs_MultiPhaseEquil::equilibrate(), vcs_MultiPhaseEquil::equilibrate_HP(), vcs_MultiPhaseEquil::equilibrate_SP(), vcs_MultiPhaseEquil::equilibrate_TV(), and MultiPhase::updatePhases().

◆ maxTemp()

doublereal maxTemp ( ) const
inline

Maximum temperature for which all solution phases have valid thermo data.

Stoichiometric phases are not considered, since they may have thermo data only valid for conditions for which they are stable.

Definition at line 254 of file MultiPhase.h.

References MultiPhase::m_Tmax.

Referenced by vcs_MultiPhaseEquil::equilibrate(), vcs_MultiPhaseEquil::equilibrate_HP(), vcs_MultiPhaseEquil::equilibrate_SP(), and vcs_MultiPhaseEquil::equilibrate_TV().

◆ charge()

doublereal charge ( ) const

Total charge summed over all phases (Coulombs).

Definition at line 205 of file MultiPhase.cpp.

References MultiPhase::nPhases(), and MultiPhase::phaseCharge().

◆ phaseCharge()

doublereal phaseCharge ( size_t  p) const

Charge (Coulombs) of phase with index p.

The net charge is computed as

\[ Q_p = N_p \sum_k F z_k X_k \]

where the sum runs only over species in phase p.

Parameters
pindex of the phase for which the charge is desired.

Definition at line 230 of file MultiPhase.cpp.

References Cantera::Faraday, MultiPhase::m_moleFractions, MultiPhase::m_moles, MultiPhase::m_phase, and MultiPhase::speciesIndex().

Referenced by MultiPhase::charge().

◆ elementMoles()

doublereal elementMoles ( size_t  m) const

Total moles of global element m, summed over all phases.

Parameters
mIndex of the global element

Definition at line 190 of file MultiPhase.cpp.

References MultiPhase::m_atoms, MultiPhase::m_moleFractions, MultiPhase::m_moles, MultiPhase::m_phase, MultiPhase::nPhases(), and MultiPhase::speciesIndex().

Referenced by MultiPhaseEquil::MultiPhaseEquil().

◆ getChemPotentials()

void getChemPotentials ( doublereal *  mu) const

Returns a vector of Chemical potentials.

Write into array mu the chemical potentials of all species [J/kmol]. The chemical potentials are related to the activities by

\( \mu_k = \mu_k^0(T, P) + RT \ln a_k. \).

Parameters
muChemical potential vector. Length = num global species. Units = J/kmol.

Definition at line 241 of file MultiPhase.cpp.

References MultiPhase::m_phase, MultiPhase::nPhases(), and MultiPhase::updatePhases().

Referenced by vcs_MultiPhaseEquil::equilibrate_TP().

◆ getValidChemPotentials()

void getValidChemPotentials ( doublereal  not_mu,
doublereal *  mu,
bool  standard = false 
) const

Returns a vector of Valid chemical potentials.

Write into array mu the chemical potentials of all species with thermo data valid for the current temperature [J/kmol]. For other species, set the chemical potential to the value not_mu. If standard is set to true, then the values returned are standard chemical potentials.

This method is designed for use in computing chemical equilibrium by Gibbs minimization. For solution phases (more than one species), this does the same thing as getChemPotentials. But for stoichiometric phases, this writes into array mu the user-specified value not_mu instead of the chemical potential if the temperature is outside the range for which the thermo data for the one species in the phase are valid. The need for this arises since many condensed phases have thermo data fit only for the temperature range for which they are stable. For example, in the NASA database, the fits for H2O(s) are only done up to 0 C, the fits for H2O(L) are only done from 0 C to 100 C, etc. Using the polynomial fits outside the range for which the fits were done can result in spurious chemical potentials, and can lead to condensed phases appearing when in fact they should be absent.

By setting not_mu to a large positive value, it is possible to force routines which seek to minimize the Gibbs free energy of the mixture to zero out any phases outside the temperature range for which their thermo data are valid.

Parameters
not_muValue of the chemical potential to set species in phases, for which the thermo data is not valid
muVector of chemical potentials. length = Global species, units = J kmol-1
standardIf this method is called with standard set to true, then the composition-independent standard chemical potentials are returned instead of the composition- dependent chemical potentials.

Definition at line 251 of file MultiPhase.cpp.

References MultiPhase::m_phase, MultiPhase::nPhases(), MultiPhase::nSpecies(), MultiPhase::tempOK(), and MultiPhase::updatePhases().

Referenced by MultiPhaseEquil::computeReactionSteps(), MultiPhaseEquil::setInitialMoles(), and MultiPhaseEquil::stepComposition().

◆ temperature()

doublereal temperature ( ) const
inline

◆ setTemperature()

void setTemperature ( const doublereal  T)

Set the temperature [K].

Parameters
Tvalue of the temperature (Kelvin)

Definition at line 707 of file MultiPhase.cpp.

Referenced by vcs_MultiPhaseEquil::equilibrate_HP(), MultiPhase::equilibrate_MultiPhaseEquil(), vcs_MultiPhaseEquil::equilibrate_SP(), and vcs_MultiPhaseEquil::equilibrate_TV().

◆ setState_TP()

void setState_TP ( const doublereal  T,
const doublereal  Pres 
)

Set the state of the underlying ThermoPhase objects in one call.

Parameters
TTemperature of the system (kelvin)
Prespressure of the system (pascal)

Definition at line 424 of file MultiPhase.cpp.

◆ setState_TPMoles()

void setState_TPMoles ( const doublereal  T,
const doublereal  Pres,
const doublereal *  Moles 
)

Set the state of the underlying ThermoPhase objects in one call.

Parameters
TTemperature of the system (kelvin)
Prespressure of the system (pascal)
MolesVector of mole numbers of all the species in all the phases (kmol)

Definition at line 434 of file MultiPhase.cpp.

◆ pressure()

doublereal pressure ( ) const
inline

◆ volume()

doublereal volume ( ) const

The total mixture volume [m^3].

Returns the cumulative sum of the volumes of all the phases in the mixture.

Definition at line 472 of file MultiPhase.cpp.

References MultiPhase::m_moles, MultiPhase::m_phase, and MultiPhase::nPhases().

Referenced by vcs_MultiPhaseEquil::equilibrate_TV(), and VCS_SOLVE::vcs_prob_specifyFully().

◆ setPressure()

void setPressure ( doublereal  P)
inline

Set the pressure [Pa].

Parameters
PSet the pressure in the MultiPhase object (Pa)

Definition at line 404 of file MultiPhase.h.

Referenced by vcs_MultiPhaseEquil::equilibrate_TV().

◆ enthalpy()

doublereal enthalpy ( ) const

◆ IntEnergy()

doublereal IntEnergy ( ) const

The internal energy of the mixture [J].

Definition at line 304 of file MultiPhase.cpp.

References MultiPhase::m_moles, MultiPhase::m_phase, MultiPhase::nPhases(), and MultiPhase::updatePhases().

Referenced by vcs_MultiPhaseEquil::equilibrate(), and vcs_MultiPhaseEquil::equilibrate_HP().

◆ entropy()

doublereal entropy ( ) const

◆ gibbs()

doublereal gibbs ( ) const

The Gibbs function of the mixture [J].

Definition at line 280 of file MultiPhase.cpp.

References MultiPhase::m_moles, MultiPhase::m_phase, MultiPhase::nPhases(), and MultiPhase::updatePhases().

◆ cp()

doublereal cp ( ) const

Heat capacity at constant pressure [J/K].

Note that this does not account for changes in composition of the mixture with temperature.

Definition at line 328 of file MultiPhase.cpp.

References MultiPhase::m_moles, MultiPhase::m_phase, MultiPhase::nPhases(), and MultiPhase::updatePhases().

◆ nPhases()

size_t nPhases ( ) const
inline

◆ solutionSpecies()

bool solutionSpecies ( size_t  kGlob) const

Return true is species kGlob is a species in a multicomponent solution phase.

Parameters
kGlobindex of the global species

Definition at line 271 of file MultiPhase.cpp.

References MultiPhase::m_phase, MultiPhase::m_spphase, and MultiPhase::nSpecies().

Referenced by MultiPhaseEquil::getComponents(), and MultiPhaseEquil::MultiPhaseEquil().

◆ speciesPhaseIndex()

size_t speciesPhaseIndex ( const size_t  kGlob) const

Returns the phase index of the Kth "global" species.

Parameters
kGlobGlobal species index.
Returns
the index of the owning phase.

Definition at line 800 of file MultiPhase.cpp.

References MultiPhase::m_spphase.

Referenced by MultiPhaseEquil::computeReactionSteps(), and MultiPhaseEquil::MultiPhaseEquil().

◆ moleFraction()

doublereal moleFraction ( const size_t  kGlob) const

Returns the mole fraction of global species k.

Parameters
kGlobIndex of the global species.

Definition at line 805 of file MultiPhase.cpp.

References MultiPhase::m_moleFractions.

Referenced by vcs_MultiPhaseEquil::equilibrate_TP().

◆ setPhaseMoleFractions()

void setPhaseMoleFractions ( const size_t  n,
const doublereal *const  x 
)

Set the Mole fractions of the nth phase.

This function sets the mole fractions of the nth phase. Note, the mole number of the phase stays constant

Parameters
nindex of the phase
xVector of input mole fractions.

Definition at line 340 of file MultiPhase.cpp.

References MultiPhase::init(), MultiPhase::m_init, MultiPhase::m_moleFractions, MultiPhase::m_phase, MultiPhase::m_press, MultiPhase::m_spstart, MultiPhase::m_temp, Phase::nSpecies(), and ThermoPhase::setState_TPX().

◆ setMolesByName() [1/2]

void setMolesByName ( const compositionMap xMap)

Set the number of moles of species in the mixture.

Parameters
xMapCompositionMap of the species with nonzero mole numbers. Mole numbers that are less than or equal to zero will be set to zero. units = kmol.

Definition at line 353 of file MultiPhase.cpp.

References Cantera::getValue(), MultiPhase::nSpecies(), MultiPhase::setMoles(), and MultiPhase::speciesName().

Referenced by MultiPhase::setMolesByName().

◆ setMolesByName() [2/2]

void setMolesByName ( const std::string &  x)

Set the moles via a string containing their names.

The string x is in the form of a composition map. Species which are not listed are set to zero.

Parameters
xstring x in the form of a composition map where values are the moles of the species.

Definition at line 363 of file MultiPhase.cpp.

References MultiPhase::m_snames, Cantera::parseCompString(), and MultiPhase::setMolesByName().

◆ getMoles()

void getMoles ( doublereal *  molNum) const

Get the mole numbers of all species in the multiphase object.

Parameters
[out]molNumVector of doubles of length nSpecies containing the global mole numbers (kmol).

Definition at line 370 of file MultiPhase.cpp.

References MultiPhase::m_moleFractions, MultiPhase::m_moles, MultiPhase::m_phase, MultiPhase::nPhases(), and Phase::nSpecies().

Referenced by MultiPhase::addSpeciesMoles().

◆ setMoles()

void setMoles ( const doublereal *  n)

Sets all of the global species mole numbers.

The state of each phase object is also updated to have the specified composition and the mixture temperature and pressure.

Parameters
nVector of doubles of length nSpecies containing the global mole numbers (kmol).

Definition at line 385 of file MultiPhase.cpp.

References Phase::getMoleFractions(), MultiPhase::init(), MultiPhase::m_init, MultiPhase::m_moleFractions, MultiPhase::m_moles, MultiPhase::m_phase, MultiPhase::m_press, MultiPhase::m_temp, MultiPhase::nPhases(), Phase::nSpecies(), and ThermoPhase::setState_TPX().

Referenced by MultiPhase::addSpeciesMoles(), MultiPhaseEquil::finish(), and MultiPhase::setMolesByName().

◆ addSpeciesMoles()

void addSpeciesMoles ( const int  indexS,
const doublereal  addedMoles 
)

Adds moles of a certain species to the mixture.

Parameters
indexSIndex of the species in the MultiPhase object
addedMolesValue of the moles that are added to the species.

Definition at line 415 of file MultiPhase.cpp.

References MultiPhase::getMoles(), MultiPhase::m_nsp, and MultiPhase::setMoles().

◆ getElemAbundances()

void getElemAbundances ( doublereal *  elemAbundances) const

Retrieves a vector of element abundances.

Parameters
elemAbundancesVector of element abundances Length = number of elements in the MultiPhase object. Index is the global element index. Units is in kmol.

Definition at line 442 of file MultiPhase.cpp.

References MultiPhase::calcElemAbundances(), MultiPhase::m_elemAbundances, and MultiPhase::m_nel.

◆ tempOK()

bool tempOK ( size_t  p) const

Return true if the phase p has valid thermo data for the current temperature.

Parameters
pIndex of the phase.

Definition at line 810 of file MultiPhase.cpp.

References MultiPhase::m_temp_OK.

Referenced by MultiPhase::getValidChemPotentials(), and MultiPhaseEquil::MultiPhaseEquil().

◆ uploadMoleFractionsFromPhases()

void uploadMoleFractionsFromPhases ( )

Update the locally-stored composition within this object to match the current compositions of the phase objects.

Query the underlying ThermoPhase objects for their mole fractions and fill in the mole fraction vector of this current object. Adjust element compositions within this object to match.

This is an upload operation in the sense that we are taking downstream information (ThermoPhase object info) and applying it to an upstream object (MultiPhase object).

Definition at line 815 of file MultiPhase.cpp.

References MultiPhase::calcElemAbundances(), Phase::getMoleFractions(), MultiPhase::m_moleFractions, MultiPhase::m_phase, MultiPhase::nPhases(), and Phase::nSpecies().

Referenced by MultiPhase::init(), and VCS_SOLVE::vcs_prob_update().

◆ updatePhases()

void updatePhases ( ) const

Set the states of the phase objects to the locally-stored state within this MultiPhase object.

This method sets each phase to the mixture temperature and pressure, and sets the phase mole fractions based on the mixture mole numbers.

This is an download operation in the sense that we are taking upstream object information (MultiPhase object) and applying it to downstream objects (ThermoPhase object information)

Therefore, the term, "update", is appropriate for a downstream operation.

Definition at line 826 of file MultiPhase.cpp.

References MultiPhase::m_moleFractions, MultiPhase::m_phase, MultiPhase::m_press, MultiPhase::m_temp, MultiPhase::m_temp_OK, MultiPhase::minTemp(), and MultiPhase::nPhases().

Referenced by MultiPhase::cp(), MultiPhase::enthalpy(), MultiPhase::entropy(), MultiPhase::equilibrate(), MultiPhase::getChemPotentials(), MultiPhase::getValidChemPotentials(), MultiPhase::gibbs(), MultiPhase::init(), MultiPhase::IntEnergy(), and Cantera::operator<<().

◆ calcElemAbundances()

void calcElemAbundances ( ) const
private

◆ equilibrate_MultiPhaseEquil()

double equilibrate_MultiPhaseEquil ( int  XY,
doublereal  err,
int  maxsteps,
int  maxiter,
int  loglevel 
)
private

Set the mixture to a state of chemical equilibrium using the MultiPhaseEquil solver.

Parameters
XYInteger flag specifying properties to hold fixed.
errError tolerance for \(\Delta \mu/RT \) for all reactions. Also used as the relative error tolerance for the outer loop.
maxstepsMaximum number of steps to take in solving the fixed TP problem.
maxiterMaximum number of "outer" iterations for problems holding fixed something other than (T,P).
loglevelLevel of diagnostic output

Definition at line 482 of file MultiPhase.cpp.

References MultiPhase::enthalpy(), MultiPhase::init(), MultiPhase::m_init, MultiPhase::m_temp, MultiPhase::m_Tmax, MultiPhase::m_Tmin, MultiPhase::setTemperature(), and Cantera::Undef.

Referenced by MultiPhase::equilibrate().

Member Data Documentation

◆ m_moles

vector_fp m_moles
private

◆ m_phase

std::vector<ThermoPhase*> m_phase
private

◆ m_atoms

DenseMatrix m_atoms
private

Global Stoichiometric Coefficient array.

This is a two dimensional array m_atoms(m, k). The first index is the global element index. The second index, k, is the global species index. The value is the number of atoms of type m in species k.

Definition at line 580 of file MultiPhase.h.

Referenced by MultiPhase::calcElemAbundances(), MultiPhase::elementMoles(), MultiPhase::init(), and MultiPhase::nAtoms().

◆ m_moleFractions

vector_fp m_moleFractions
private

◆ m_spphase

std::vector<size_t> m_spphase
private

Mapping between the global species number and the phase ID.

m_spphase[kGlobal] = iPhase Length = number of global species

Definition at line 591 of file MultiPhase.h.

Referenced by MultiPhase::init(), MultiPhase::solutionSpecies(), MultiPhase::speciesMoles(), and MultiPhase::speciesPhaseIndex().

◆ m_spstart

std::vector<size_t> m_spstart
private

Vector of ints containing of first species index in the global list of species for each phase.

kfirst = m_spstart[ip], kfirst is the index of the first species in the ip'th phase.

Definition at line 599 of file MultiPhase.h.

Referenced by MultiPhase::init(), MultiPhase::phase(), MultiPhase::setPhaseMoleFractions(), and MultiPhase::speciesIndex().

◆ m_enames

std::vector<std::string> m_enames
private

String names of the global elements.

This has a length equal to the number of global elements.

Definition at line 603 of file MultiPhase.h.

Referenced by MultiPhase::addPhase(), MultiPhase::elementIndex(), MultiPhase::elementName(), and MultiPhase::init().

◆ m_atomicNumber

vector_int m_atomicNumber
private

Atomic number of each global element.

Definition at line 606 of file MultiPhase.h.

Referenced by MultiPhase::addPhase().

◆ m_snames

std::vector<std::string> m_snames
private

Vector of species names in the problem.

Vector is over all species defined in the object, the global species index.

Definition at line 610 of file MultiPhase.h.

Referenced by MultiPhase::init(), MultiPhase::setMolesByName(), and MultiPhase::speciesName().

◆ m_enamemap

std::map<std::string, size_t> m_enamemap
private

Returns the global element index, given the element string name.

-> used in the construction. However, wonder if it needs to be global.

Definition at line 616 of file MultiPhase.h.

Referenced by MultiPhase::addPhase().

◆ m_temp

doublereal m_temp
private

◆ m_press

doublereal m_press
private

◆ m_nel

size_t m_nel
private

◆ m_nsp

size_t m_nsp
private

◆ m_init

bool m_init
private

◆ m_eloc

size_t m_eloc
private

Global ID of the element corresponding to the electronic charge.

If there is none, then this is equal to -1

Definition at line 635 of file MultiPhase.h.

Referenced by MultiPhase::addPhase().

◆ m_temp_OK

std::vector<bool> m_temp_OK
mutableprivate

Vector of bools indicating whether temperatures are ok for phases.

If the current temperature is outside the range of valid temperatures for the phase thermodynamics, the phase flag is set to false.

Definition at line 642 of file MultiPhase.h.

Referenced by MultiPhase::addPhase(), MultiPhase::tempOK(), and MultiPhase::updatePhases().

◆ m_Tmin

doublereal m_Tmin
private

Minimum temperature for which thermo parameterizations are valid.

Stoichiometric phases are ignored in this determination. units Kelvin

Definition at line 646 of file MultiPhase.h.

Referenced by MultiPhase::addPhase(), MultiPhase::equilibrate_MultiPhaseEquil(), and MultiPhase::minTemp().

◆ m_Tmax

doublereal m_Tmax
private

Minimum temperature for which thermo parameterizations are valid.

Stoichiometric phases are ignored in this determination. units Kelvin

Definition at line 650 of file MultiPhase.h.

Referenced by MultiPhase::addPhase(), MultiPhase::equilibrate_MultiPhaseEquil(), and MultiPhase::maxTemp().

◆ m_elemAbundances

vector_fp m_elemAbundances
mutableprivate

Vector of element abundances.

m_elemAbundances[mGlobal] = kmol of element mGlobal summed over all species in all phases.

Definition at line 657 of file MultiPhase.h.

Referenced by MultiPhase::calcElemAbundances(), MultiPhase::getElemAbundances(), and MultiPhase::init().


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