Cantera
2.0
|
A class for multiphase mixtures. More...
#include <MultiPhase.h>
Public Types | |
typedef size_t | index_t |
Shorthand for an index variable that can't be negative. | |
Public Member Functions | |
MultiPhase () | |
Constructor. | |
MultiPhase (const MultiPhase &right) | |
Copy Constructor. | |
virtual | ~MultiPhase () |
Destructor. | |
MultiPhase & | operator= (const MultiPhase &right) |
Assignment operator. | |
void | addPhases (std::vector< ThermoPhase * > &phases, const vector_fp &phaseMoles) |
Add a vector of phases to the mixture. | |
void | addPhases (MultiPhase &mix) |
Add all phases present in 'mix' to this mixture. | |
void | addPhase (ThermoPhase *p, doublereal moles) |
Add a phase to the mixture. | |
size_t | nElements () const |
Number of elements. | |
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. | |
void | checkElementArraySize (size_t mm) const |
Check that an array size is at least nElements() Throws an exception if mm is less than nElements(). | |
std::string | elementName (size_t m) const |
Returns the string name of the global element m. | |
size_t | elementIndex (std::string name) const |
Returns the index of the element with name name. | |
size_t | nSpecies () const |
Number of species, summed over all phases. | |
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. | |
void | checkSpeciesArraySize (size_t kk) const |
Check that an array size is at least nSpecies() Throws an exception if kk is less than nSpecies(). | |
std::string | speciesName (const size_t kGlob) const |
Name of species with global index kGlob. | |
doublereal | nAtoms (const size_t kGlob, const size_t mGlob) const |
Returns the Number of atoms of global element mGlob in global species kGlob. | |
void | getMoleFractions (doublereal *const x) const |
Returns the global Species mole fractions. | |
void | init () |
Process phases and build atomic composition array. | |
std::string | phaseName (const index_t iph) const |
Returns the name of the n'th phase. | |
int | phaseIndex (const std::string &pName) const |
Returns the index, given the phase name. | |
doublereal | phaseMoles (const index_t n) const |
Return the number of moles in phase n. | |
void | setPhaseMoles (const index_t n, const doublereal moles) |
Set the number of moles of phase with index n. | |
ThermoPhase & | phase (index_t n) |
Return a ThermoPhase reference to phase n. | |
void | checkPhaseIndex (size_t m) const |
Check that the specified phase index is in range Throws an exception if m is greater than nPhases() | |
void | checkPhaseArraySize (size_t mm) const |
Check that an array size is at least nPhases() Throws an exception if mm is less than nPhases(). | |
doublereal | speciesMoles (index_t kGlob) const |
Returns the moles of global species k . | |
size_t | speciesIndex (index_t k, index_t p) const |
Return the global index of the species belonging to phase number p with local index k within the phase. | |
size_t | speciesIndex (std::string speciesName, std::string phaseName) |
Return the global index of the species belonging to phase name phaseName with species name speciesName . | |
doublereal | minTemp () const |
Minimum temperature for which all solution phases have valid thermo data. | |
doublereal | maxTemp () const |
Maximum temperature for which all solution phases have valid thermo data. | |
doublereal | charge () const |
Total charge summed over all phases (Coulombs). | |
doublereal | phaseCharge (index_t p) const |
Charge (Coulombs) of phase with index p. | |
doublereal | elementMoles (index_t m) const |
Total moles of global element m, summed over all phases. | |
void | getChemPotentials (doublereal *mu) const |
Returns a vector of Chemical potentials. | |
void | getValidChemPotentials (doublereal not_mu, doublereal *mu, bool standard=false) const |
Returns a vector of Valid chemical potentials. | |
doublereal | temperature () const |
Temperature [K]. | |
doublereal | equilibrate (int XY, doublereal err=1.0e-9, int maxsteps=1000, int maxiter=200, int loglevel=-99) |
Set the mixture to a state of chemical equilibrium. | |
void | setTemperature (const doublereal T) |
Set the temperature [K]. | |
void | setState_TP (const doublereal T, const doublereal Pres) |
Set the state of the underlying ThermoPhase objects in one call. | |
void | setState_TPMoles (const doublereal T, const doublereal Pres, const doublereal *Moles) |
Set the state of the underlying ThermoPhase objects in one call. | |
doublereal | pressure () const |
Pressure [Pa]. | |
doublereal | volume () const |
Volume [m^3]. | |
void | setPressure (doublereal P) |
Set the pressure [Pa]. | |
doublereal | enthalpy () const |
Enthalpy [J]. | |
doublereal | IntEnergy () const |
Enthalpy [J]. | |
doublereal | entropy () const |
Entropy [J/K]. | |
doublereal | gibbs () const |
Gibbs function [J]. | |
doublereal | cp () const |
Heat capacity at constant pressure [J/K]. | |
index_t | nPhases () const |
Number of phases. | |
bool | solutionSpecies (index_t kGlob) const |
Return true is species kGlob is a species in a multicomponent solution phase. | |
size_t | speciesPhaseIndex (const index_t kGlob) const |
Returns the phase index of the Kth "global" species. | |
doublereal | moleFraction (const index_t kGlob) const |
Returns the mole fraction of global species k. | |
void | setPhaseMoleFractions (const index_t n, const doublereal *const x) |
Set the Mole fractions of the nth phase. | |
void | setMolesByName (compositionMap &xMap) |
Set the number numbers of species in the MultiPhase. | |
void | setMolesByName (const std::string &x) |
Set the Moles via a string containing their names. | |
void | getMoles (doublereal *molNum) const |
Return a vector of global species mole numbers. | |
void | setMoles (const doublereal *n) |
Sets all of the global species mole numbers. | |
void | addSpeciesMoles (const int indexS, const doublereal addedMoles) |
Adds moles of a certain species to the mixture. | |
void | getElemAbundances (doublereal *elemAbundances) const |
Retrieves a vector of element abundances. | |
bool | tempOK (index_t p) const |
Return true if the phase p has valid thermo data for the current temperature. | |
void | updateMoleFractions () |
Update the locally-stored composition within this object to match the current compositions of the phase objects. | |
void | uploadMoleFractionsFromPhases () |
Update the locally-stored composition within this object to match the current compositions of the phase objects. | |
void | updatePhases () const |
Set the states of the phase objects to the locally-stored state within this MultiPhase object. | |
Private Member Functions | |
void | calcElemAbundances () const |
Calculate the element abundance vector. | |
Private Attributes | |
vector_fp | m_moles |
Vector of the number of moles in each phase. | |
std::vector< ThermoPhase * > | m_phase |
Vector of the ThermoPhase Pointers. | |
DenseMatrix | m_atoms |
Global Stoichiometric Coefficient array. | |
vector_fp | m_moleFractions |
Locally stored vector of mole fractions of all species comprising the MultiPhase object. | |
std::vector< size_t > | m_spphase |
Mapping between the global species number and the phase ID. | |
std::vector< size_t > | m_spstart |
Vector of ints containing of first species index in the global list of species for each phase. | |
std::vector< std::string > | m_enames |
String names of the global elements. | |
vector_int | m_atomicNumber |
Atomic number of each element. | |
std::vector< std::string > | m_snames |
Vector of species names in the problem. | |
std::map< std::string, size_t > | m_enamemap |
Returns the global element index, given the element string name. | |
index_t | m_np |
Number of phases in the MultiPhase object. | |
doublereal | m_temp |
Current value of the temperature (kelvin) | |
doublereal | m_press |
Current value of the pressure (Pa) | |
index_t | m_nel |
Number of distinct elements in all of the phases. | |
index_t | m_nsp |
Number of distinct species in all of the phases. | |
bool | m_init |
True if the init() routine has been called, and the MultiPhase frozen. | |
size_t | m_eloc |
Global ID of the element corresponding to the electronic charge. | |
std::vector< bool > | m_temp_OK |
Vector of bools indicating whether temperatures are ok for phases. | |
doublereal | m_Tmin |
Minimum temperature for which thermo parameterizations are valid. | |
doublereal | m_Tmax |
Minimum temperature for which thermo parameterizations are valid. | |
vector_fp | m_elemAbundances |
Vector of element abundances. | |
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 60 of file MultiPhase.h.
typedef size_t index_t |
Shorthand for an index variable that can't be negative.
Definition at line 65 of file MultiPhase.h.
MultiPhase | ( | ) |
Constructor.
The constructor takes no arguments, since phases are added using method addPhase().
Definition at line 21 of file MultiPhase.cpp.
MultiPhase | ( | const MultiPhase & | right | ) |
Copy Constructor.
right | Object to be copied |
Definition at line 38 of file MultiPhase.cpp.
References MultiPhase::operator=().
|
virtual |
Destructor.
Does nothing. Class MultiPhase does not take "ownership" (i.e. responsibility for destroying) the phase objects.
Definition at line 58 of file MultiPhase.cpp.
MultiPhase & operator= | ( | const MultiPhase & | right | ) |
Assignment operator.
right | Object to be copied |
Definition at line 66 of file MultiPhase.cpp.
References MultiPhase::m_atoms, MultiPhase::m_elemAbundances, MultiPhase::m_eloc, MultiPhase::m_enamemap, MultiPhase::m_enames, MultiPhase::m_init, MultiPhase::m_moleFractions, MultiPhase::m_moles, MultiPhase::m_nel, MultiPhase::m_np, MultiPhase::m_nsp, MultiPhase::m_phase, MultiPhase::m_press, MultiPhase::m_spphase, MultiPhase::m_spstart, MultiPhase::m_temp, MultiPhase::m_temp_OK, MultiPhase::m_Tmax, and MultiPhase::m_Tmin.
Referenced by MultiPhase::MultiPhase().
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
phases | Vector of pointers to phases |
phaseMoles | Vector of mole numbers in each phase (kmol) |
Definition at line 103 of file MultiPhase.cpp.
References MultiPhase::addPhase(), and MultiPhase::init().
void addPhases | ( | MultiPhase & | mix | ) |
Add all phases present in 'mix' to this mixture.
mix | Add all of the phases in another MultiPhase object to the current object. |
Definition at line 94 of file MultiPhase.cpp.
References MultiPhase::addPhase(), MultiPhase::m_moles, MultiPhase::m_np, and MultiPhase::m_phase.
void addPhase | ( | ThermoPhase * | p, |
doublereal | moles | ||
) |
Add a phase to the mixture.
This function must be called before the init() function is called, which serves to freeze the MultiPhase.
p | pointer to the phase object |
moles | total number of moles of all species in this phase |
Definition at line 114 of file MultiPhase.cpp.
References Phase::atomicNumber(), Phase::elementName(), MultiPhase::m_atomicNumber, MultiPhase::m_eloc, MultiPhase::m_enamemap, MultiPhase::m_enames, MultiPhase::m_init, MultiPhase::m_moles, MultiPhase::m_nel, MultiPhase::m_np, MultiPhase::m_nsp, MultiPhase::m_phase, MultiPhase::m_press, MultiPhase::m_temp, MultiPhase::m_temp_OK, MultiPhase::m_Tmax, MultiPhase::m_Tmin, ThermoPhase::maxTemp(), ThermoPhase::minTemp(), Phase::nElements(), Phase::nSpecies(), ThermoPhase::pressure(), and Phase::temperature().
Referenced by Crystal::addLattice(), MultiPhase::addPhases(), Cantera::equilibrate(), ChemEquil::estimateElementPotentials(), ChemEquil::setInitialMoles(), and Cantera::vcs_equilibrate().
|
inline |
Number of elements.
Definition at line 121 of file MultiPhase.h.
References MultiPhase::m_nel.
Referenced by vcs_MultiPhaseEquil::determine_PhaseStability(), vcs_MultiPhaseEquil::equilibrate_TP(), MultiPhaseEquil::MultiPhaseEquil(), and vcs_MultiPhaseEquil::vcs_MultiPhaseEquil().
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 1026 of file MultiPhase.cpp.
References MultiPhase::m_nel.
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 1033 of file MultiPhase.cpp.
References MultiPhase::m_nel.
std::string elementName | ( | size_t | m | ) | const |
Returns the string name of the global element m.
m | index of the global element |
Definition at line 1042 of file MultiPhase.cpp.
References MultiPhase::m_enames.
Referenced by MultiPhaseEquil::MultiPhaseEquil().
size_t elementIndex | ( | std::string | name | ) | const |
Returns the index of the element with name name.
name | String name of the global element |
Definition at line 1048 of file MultiPhase.cpp.
References MultiPhase::m_enames, MultiPhase::m_nel, and Cantera::npos.
|
inline |
Number of species, summed over all phases.
Definition at line 147 of file MultiPhase.h.
References MultiPhase::m_nsp.
Referenced by vcs_MultiPhaseEquil::determine_PhaseStability(), vcs_MultiPhaseEquil::equilibrate_TP(), MultiPhase::getValidChemPotentials(), MultiPhaseEquil::MultiPhaseEquil(), MultiPhase::setMolesByName(), MultiPhase::solutionSpecies(), and vcs_MultiPhaseEquil::vcs_MultiPhaseEquil().
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 1058 of file MultiPhase.cpp.
References MultiPhase::m_nsp.
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 1065 of file MultiPhase.cpp.
References MultiPhase::m_nsp.
std::string speciesName | ( | const size_t | kGlob | ) | const |
Name of species with global index kGlob.
kGlob | global species index |
Definition at line 1074 of file MultiPhase.cpp.
References MultiPhase::m_snames.
Referenced by MultiPhaseEquil::MultiPhaseEquil(), MultiPhaseEquil::reactionString(), MultiPhase::setMolesByName(), and MultiPhaseEquil::stepComposition().
doublereal nAtoms | ( | const size_t | kGlob, |
const size_t | mGlob | ||
) | const |
Returns the Number of atoms of global element mGlob in global species kGlob.
kGlob | global species index |
mGlob | global element index |
Definition at line 1079 of file MultiPhase.cpp.
References MultiPhase::m_atoms.
Referenced by MultiPhaseEquil::getComponents(), and MultiPhaseEquil::MultiPhaseEquil().
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.
x | vector of mole fractions. Length = number of global species. |
Definition at line 1084 of file MultiPhase.cpp.
References MultiPhase::m_moleFractions.
void init | ( | ) |
Process phases and build atomic composition array.
This method must be called after all phases are added, before doing anything else with the mixture. After init() has been called, no more phases may be added.
Definition at line 193 of file MultiPhase.cpp.
References Phase::elementIndex(), MultiPhase::m_atomicNumber, MultiPhase::m_atoms, MultiPhase::m_elemAbundances, MultiPhase::m_eloc, MultiPhase::m_enames, MultiPhase::m_init, MultiPhase::m_moleFractions, MultiPhase::m_nel, MultiPhase::m_np, MultiPhase::m_nsp, MultiPhase::m_phase, MultiPhase::m_snames, MultiPhase::m_spphase, MultiPhase::m_spstart, Phase::nAtoms(), Cantera::npos, Phase::nSpecies(), DenseMatrix::resize(), Phase::speciesName(), MultiPhase::updatePhases(), and MultiPhase::uploadMoleFractionsFromPhases().
Referenced by MultiPhase::addPhases(), Cantera::equilibrate(), MultiPhase::equilibrate(), ChemEquil::estimateElementPotentials(), MultiPhase::phase(), ChemEquil::setInitialMoles(), MultiPhase::setMoles(), MultiPhase::setPhaseMoleFractions(), MultiPhase::setState_TP(), MultiPhase::setTemperature(), MultiPhase::speciesIndex(), Cantera::vcs_determine_PhaseStability(), Cantera::vcs_equilibrate(), and Cantera::vcs_equilibrate_1().
std::string phaseName | ( | const index_t | iph | ) | const |
Returns the name of the n'th phase.
iph | phase Index |
Definition at line 1089 of file MultiPhase.cpp.
References Phase::id(), and MultiPhase::m_phase.
Referenced by vcs_MultiPhaseEquil::determine_PhaseStability().
int phaseIndex | ( | const std::string & | pName | ) | const |
Returns the index, given the phase name.
pName | Name of the phase |
Definition at line 1095 of file MultiPhase.cpp.
References Phase::id(), MultiPhase::m_np, and MultiPhase::m_phase.
Referenced by MultiPhase::speciesIndex().
doublereal phaseMoles | ( | const index_t | n | ) | const |
Return the number of moles in phase n.
n | Index of the phase. |
Definition at line 1108 of file MultiPhase.cpp.
References MultiPhase::m_moles.
Referenced by MultiPhaseEquil::computeReactionSteps(), vcs_MultiPhaseEquil::equilibrate_HP(), vcs_MultiPhaseEquil::equilibrate_SP(), and Cantera::operator<<().
void setPhaseMoles | ( | const index_t | n, |
const doublereal | moles | ||
) |
Set the number of moles of phase with index n.
n | Index of the phase |
moles | Number of moles in the phase (kmol) |
Definition at line 1113 of file MultiPhase.cpp.
References MultiPhase::m_moles.
Referenced by vcs_MultiPhaseEquil::equilibrate_TP().
ThermoPhase & phase | ( | index_t | n | ) |
Return a ThermoPhase reference to phase n.
The state of phase n is also updated to match the state stored locally in the mixture object.
n | Phase Index |
Definition at line 260 of file MultiPhase.cpp.
References DATA_PTR, 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(), vcs_MultiPhaseEquil::equilibrate_TP(), Crystal::lattice(), Cantera::operator<<(), and vcs_MultiPhaseEquil::reportCSV().
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 271 of file MultiPhase.cpp.
References MultiPhase::nPhases().
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 278 of file MultiPhase.cpp.
References MultiPhase::nPhases().
doublereal speciesMoles | ( | index_t | kGlob | ) | const |
Returns the moles of global species k
.
Moles of species k
.
Returns the moles of global species k. units = kmol
kGlob | Global species index k |
Definition at line 287 of file MultiPhase.cpp.
References MultiPhase::m_moleFractions, MultiPhase::m_moles, and MultiPhase::m_spphase.
Referenced by MultiPhaseEquil::computeReactionSteps(), and MultiPhaseEquil::MultiPhaseEquil().
Return the global index of the species belonging to phase number p
with local index k
within the phase.
Returns the index of the global species
k | local index of the species within the phase |
p | index of the phase |
Definition at line 259 of file MultiPhase.h.
References MultiPhase::m_spstart.
Referenced by MultiPhase::elementMoles(), MultiPhase::phaseCharge(), and vcs_MultiPhaseEquil::reportCSV().
size_t speciesIndex | ( | std::string | speciesName, |
std::string | phaseName | ||
) |
Return the global index of the species belonging to phase name phaseName
with species name speciesName
.
Returns the index of the global species
speciesName | Species Name |
phaseName | Phase Name |
If the species or phase name is not recognized, this routine throws a CanteraError.
Definition at line 324 of file MultiPhase.cpp.
References MultiPhase::init(), MultiPhase::m_init, MultiPhase::m_phase, MultiPhase::m_spstart, Cantera::npos, and MultiPhase::phaseIndex().
|
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 282 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().
|
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 290 of file MultiPhase.h.
References MultiPhase::m_Tmax.
Referenced by vcs_MultiPhaseEquil::equilibrate(), vcs_MultiPhaseEquil::equilibrate_HP(), vcs_MultiPhaseEquil::equilibrate_SP(), vcs_MultiPhaseEquil::equilibrate_TV(), and MultiPhase::updatePhases().
doublereal charge | ( | ) | const |
Total charge summed over all phases (Coulombs).
Definition at line 314 of file MultiPhase.cpp.
References MultiPhase::m_np, and MultiPhase::phaseCharge().
doublereal phaseCharge | ( | index_t | p | ) | const |
Charge (Coulombs) of phase with index p.
Net charge of one phase (Coulombs).
p | Phase Index |
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.
p | index of the phase for which the charge is desired. |
Definition at line 344 of file MultiPhase.cpp.
References MultiPhase::m_moleFractions, MultiPhase::m_moles, MultiPhase::m_phase, and MultiPhase::speciesIndex().
Referenced by MultiPhase::charge().
doublereal elementMoles | ( | index_t | m | ) | const |
Total moles of global element m, summed over all phases.
m | Index of the global element |
Definition at line 297 of file MultiPhase.cpp.
References MultiPhase::m_atoms, MultiPhase::m_moleFractions, MultiPhase::m_moles, MultiPhase::m_np, MultiPhase::m_phase, and MultiPhase::speciesIndex().
Referenced by MultiPhaseEquil::MultiPhaseEquil().
void getChemPotentials | ( | doublereal * | mu | ) | const |
Returns a vector of Chemical potentials.
Get the chemical potentials of all species in all phases.
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. \).
mu | Chemical potential vector. Length = num global species. Units = J/kmol. |
Definition at line 357 of file MultiPhase.cpp.
References MultiPhase::m_np, MultiPhase::m_phase, and MultiPhase::updatePhases().
Referenced by vcs_MultiPhaseEquil::determine_PhaseStability().
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.
not_mu | Value of the chemical potential to set species in phases, for which the thermo data is not valid |
mu | Vector of chemical potentials length = Global species, units = J kmol-1 |
standard | If 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 394 of file MultiPhase.cpp.
References MultiPhase::m_np, MultiPhase::m_phase, MultiPhase::nSpecies(), MultiPhase::tempOK(), and MultiPhase::updatePhases().
Referenced by MultiPhaseEquil::computeReactionSteps(), MultiPhaseEquil::setInitialMoles(), and MultiPhaseEquil::stepComposition().
|
inline |
Temperature [K].
Definition at line 371 of file MultiPhase.h.
References MultiPhase::m_temp.
Referenced by vcs_MultiPhaseEquil::determine_PhaseStability(), vcs_MultiPhaseEquil::equilibrate(), MultiPhase::equilibrate(), vcs_MultiPhaseEquil::equilibrate_HP(), vcs_MultiPhaseEquil::equilibrate_SP(), vcs_MultiPhaseEquil::equilibrate_TP(), MultiPhaseEquil::MultiPhaseEquil(), and vcs_MultiPhaseEquil::reportCSV().
doublereal equilibrate | ( | int | XY, |
doublereal | err = 1.0e-9 , |
||
int | maxsteps = 1000 , |
||
int | maxiter = 200 , |
||
int | loglevel = -99 |
||
) |
Set the mixture to a state of chemical equilibrium.
XY | Integer flag specifying properties to hold fixed. |
err | Error tolerance for \(\Delta \mu/RT \) for all reactions. Also used as the relative error tolerance for the outer loop. |
maxsteps | Maximum number of steps to take in solving the fixed TP problem. |
maxiter | Maximum number of "outer" iterations for problems holding fixed something other than (T,P). |
loglevel | Level of diagnostic output, written to a file in HTML format. |
Definition at line 679 of file MultiPhase.cpp.
References Cantera::addLogEntry(), Cantera::beginLogGroup(), MultiPhase::cp(), Cantera::endLogGroup(), MultiPhase::enthalpy(), MultiPhase::entropy(), Cantera::fp2str(), MultiPhase::init(), Cantera::int2str(), MultiPhase::m_init, MultiPhase::m_temp, MultiPhase::m_Tmax, MultiPhase::m_Tmin, MultiPhase::pressure(), CanteraError::save(), MultiPhase::setPressure(), MultiPhase::setTemperature(), MultiPhase::temperature(), Cantera::Undef, and MultiPhase::volume().
Referenced by Cantera::equilibrate(), and Cantera::vcs_equilibrate_1().
void setTemperature | ( | const doublereal | T | ) |
Set the temperature [K].
T | value of the temperature (Kelvin) |
Definition at line 1017 of file MultiPhase.cpp.
References MultiPhase::init(), MultiPhase::m_init, MultiPhase::m_temp, and MultiPhase::updatePhases().
Referenced by MultiPhase::equilibrate(), vcs_MultiPhaseEquil::equilibrate_HP(), vcs_MultiPhaseEquil::equilibrate_SP(), and vcs_MultiPhaseEquil::equilibrate_TV().
void setState_TP | ( | const doublereal | T, |
const doublereal | Pres | ||
) |
Set the state of the underlying ThermoPhase objects in one call.
T | Temperature of the system (kelvin) |
Pres | pressure of the system (pascal) (kmol) |
Definition at line 615 of file MultiPhase.cpp.
References MultiPhase::init(), MultiPhase::m_init, MultiPhase::m_press, MultiPhase::m_temp, and MultiPhase::updatePhases().
void setState_TPMoles | ( | const doublereal | T, |
const doublereal | Pres, | ||
const doublereal * | Moles | ||
) |
Set the state of the underlying ThermoPhase objects in one call.
T | Temperature of the system (kelvin) |
Pres | pressure of the system (pascal) |
Moles | Vector of mole numbers of all the species in all the phases (kmol) |
Definition at line 625 of file MultiPhase.cpp.
References MultiPhase::m_press, MultiPhase::m_temp, and MultiPhase::setMoles().
|
inline |
Pressure [Pa].
Definition at line 416 of file MultiPhase.h.
References MultiPhase::m_press.
Referenced by vcs_MultiPhaseEquil::determine_PhaseStability(), MultiPhase::equilibrate(), vcs_MultiPhaseEquil::equilibrate_TP(), vcs_MultiPhaseEquil::equilibrate_TV(), MultiPhaseEquil::MultiPhaseEquil(), and vcs_MultiPhaseEquil::reportCSV().
doublereal volume | ( | ) | const |
Volume [m^3].
The total mixture volume [m^3].
Returns the cummulative sum of the volumes of all the phases in the MultiPhase.
Definition at line 668 of file MultiPhase.cpp.
References MultiPhase::m_moles, MultiPhase::m_np, and MultiPhase::m_phase.
Referenced by MultiPhase::equilibrate(), and vcs_MultiPhaseEquil::equilibrate_TV().
|
inline |
Set the pressure [Pa].
P | Set the pressure in the MultiPhase object (Pa) |
Definition at line 431 of file MultiPhase.h.
References MultiPhase::m_press, and MultiPhase::updatePhases().
Referenced by MultiPhase::equilibrate(), and vcs_MultiPhaseEquil::equilibrate_TV().
doublereal enthalpy | ( | ) | const |
Enthalpy [J].
The enthalpy of the mixture (J).
Definition at line 440 of file MultiPhase.cpp.
References MultiPhase::m_moles, MultiPhase::m_np, MultiPhase::m_phase, and MultiPhase::updatePhases().
Referenced by vcs_MultiPhaseEquil::equilibrate(), MultiPhase::equilibrate(), and vcs_MultiPhaseEquil::equilibrate_HP().
doublereal IntEnergy | ( | ) | const |
Enthalpy [J].
The internal energy of the mixture (J).
Definition at line 454 of file MultiPhase.cpp.
References MultiPhase::m_moles, MultiPhase::m_np, MultiPhase::m_phase, and MultiPhase::updatePhases().
Referenced by vcs_MultiPhaseEquil::equilibrate(), and vcs_MultiPhaseEquil::equilibrate_HP().
doublereal entropy | ( | ) | const |
Entropy [J/K].
The entropy of the mixture (J/K).
Definition at line 468 of file MultiPhase.cpp.
References MultiPhase::m_moles, MultiPhase::m_np, MultiPhase::m_phase, and MultiPhase::updatePhases().
Referenced by vcs_MultiPhaseEquil::equilibrate(), MultiPhase::equilibrate(), and vcs_MultiPhaseEquil::equilibrate_SP().
doublereal gibbs | ( | ) | const |
Gibbs function [J].
The Gibbs free energy of the mixture (J).
Definition at line 426 of file MultiPhase.cpp.
References MultiPhase::m_moles, MultiPhase::m_np, MultiPhase::m_phase, and MultiPhase::updatePhases().
doublereal cp | ( | ) | const |
Heat capacity at constant pressure [J/K].
The specific heat at constant pressure and composition (J/K).
Note that this does not account for changes in composition of the mixture with temperature.
Definition at line 484 of file MultiPhase.cpp.
References MultiPhase::m_moles, MultiPhase::m_np, MultiPhase::m_phase, and MultiPhase::updatePhases().
Referenced by MultiPhase::equilibrate().
|
inline |
Number of phases.
Definition at line 452 of file MultiPhase.h.
References MultiPhase::m_np.
Referenced by MultiPhase::checkPhaseArraySize(), MultiPhase::checkPhaseIndex(), vcs_MultiPhaseEquil::determine_PhaseStability(), vcs_MultiPhaseEquil::equilibrate_TP(), MultiPhaseEquil::MultiPhaseEquil(), Cantera::operator<<(), and vcs_MultiPhaseEquil::vcs_MultiPhaseEquil().
bool solutionSpecies | ( | index_t | kGlob | ) | const |
Return true is species kGlob is a species in a multicomponent solution phase.
True if species k belongs to a solution phase.
kGlob | index of the global species |
Definition at line 416 of file MultiPhase.cpp.
References MultiPhase::m_phase, MultiPhase::m_spphase, and MultiPhase::nSpecies().
Referenced by MultiPhaseEquil::getComponents(), and MultiPhaseEquil::MultiPhaseEquil().
size_t speciesPhaseIndex | ( | const index_t | kGlob | ) | const |
Returns the phase index of the Kth "global" species.
kGlob | Global species index. |
Definition at line 1118 of file MultiPhase.cpp.
References MultiPhase::m_spphase.
Referenced by MultiPhaseEquil::computeReactionSteps(), and MultiPhaseEquil::MultiPhaseEquil().
doublereal moleFraction | ( | const index_t | kGlob | ) | const |
Returns the mole fraction of global species k.
kGlob | Index of the global species. |
Definition at line 1123 of file MultiPhase.cpp.
References MultiPhase::m_moleFractions.
void setPhaseMoleFractions | ( | const index_t | n, |
const doublereal *const | x | ||
) |
Set the Mole fractions of the nth phase.
Set the mole fractions of phase n to the values in array x.
This function sets the mole fractions of the nth phase. Note, the mole number of the phase stays constant
n | ID of the phase |
x | Vector of input mole fractions. |
Definition at line 501 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().
void setMolesByName | ( | compositionMap & | xMap | ) |
Set the number numbers of species in the MultiPhase.
xMap | CompositionMap of the species with nonzero mole numbers units = kmol. |
Definition at line 517 of file MultiPhase.cpp.
References DATA_PTR, MultiPhase::nSpecies(), MultiPhase::setMoles(), and MultiPhase::speciesName().
Referenced by MultiPhase::setMolesByName().
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 by name in the composition map are set to zero.
x | string x in the form of a composition map where values are the moles of the species. |
Definition at line 533 of file MultiPhase.cpp.
References MultiPhase::nSpecies(), Cantera::parseCompString(), MultiPhase::setMolesByName(), and MultiPhase::speciesName().
void getMoles | ( | doublereal * | molNum | ) | const |
Return a vector of global species mole numbers.
Returns a vector of the number of moles of each species in the multiphase object.
molNum | Vector of doubles of length nSpecies containing the global mole numbers (kmol). |
Definition at line 552 of file MultiPhase.cpp.
References MultiPhase::m_moleFractions, MultiPhase::m_moles, MultiPhase::m_np, MultiPhase::m_phase, and Phase::nSpecies().
Referenced by MultiPhase::addSpeciesMoles().
void setMoles | ( | const doublereal * | n | ) |
Sets all of the global species mole numbers.
Set the species moles to the values in array n.
Sets the number of moles of each species in the multiphase object.
n | Vector of doubles of length nSpecies containing the global mole numbers (kmol). |
The state of each phase object is also updated to have the specified composition and the mixture temperature and pressure.
Definition at line 573 of file MultiPhase.cpp.
References DATA_PTR, Phase::getMoleFractions(), MultiPhase::init(), MultiPhase::m_init, MultiPhase::m_moleFractions, MultiPhase::m_moles, MultiPhase::m_np, MultiPhase::m_phase, MultiPhase::m_press, MultiPhase::m_temp, Phase::nSpecies(), and ThermoPhase::setState_TPX().
Referenced by MultiPhase::addSpeciesMoles(), MultiPhaseEquil::finish(), MultiPhase::setMolesByName(), and MultiPhase::setState_TPMoles().
void addSpeciesMoles | ( | const int | indexS, |
const doublereal | addedMoles | ||
) |
Adds moles of a certain species to the mixture.
indexS | Index of the species in the MultiPhase object |
addedMoles | Value of the moles that are added to the species. |
Definition at line 604 of file MultiPhase.cpp.
References DATA_PTR, MultiPhase::getMoles(), MultiPhase::m_nsp, and MultiPhase::setMoles().
void getElemAbundances | ( | doublereal * | elemAbundances | ) | const |
Retrieves a vector of element abundances.
elemAbundances | Vector of element abundances Length = number of elements in the MultiPhase object. Index is the global element index units is in kmol. |
Definition at line 633 of file MultiPhase.cpp.
References MultiPhase::calcElemAbundances(), MultiPhase::m_elemAbundances, and MultiPhase::m_nel.
bool tempOK | ( | index_t | p | ) | const |
Return true if the phase p has valid thermo data for the current temperature.
p | Index of the phase. |
Definition at line 1129 of file MultiPhase.cpp.
References MultiPhase::m_temp_OK.
Referenced by MultiPhase::getValidChemPotentials(), and MultiPhaseEquil::MultiPhaseEquil().
void updateMoleFractions | ( | ) |
Update the locally-stored composition within this object to match the current compositions of the phase objects.
Update the locally-stored species mole fractions.
Definition at line 1135 of file MultiPhase.cpp.
References MultiPhase::uploadMoleFractionsFromPhases().
void uploadMoleFractionsFromPhases | ( | ) |
Update the locally-stored composition within this object to match the current compositions of the phase objects.
Update the locally-stored species mole fractions.
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 1141 of file MultiPhase.cpp.
References MultiPhase::calcElemAbundances(), DATA_PTR, Phase::getMoleFractions(), MultiPhase::m_moleFractions, MultiPhase::m_np, MultiPhase::m_phase, and Phase::nSpecies().
Referenced by vcs_MultiPhaseEquil::determine_PhaseStability(), vcs_MultiPhaseEquil::equilibrate_TP(), MultiPhase::init(), and MultiPhase::updateMoleFractions().
void updatePhases | ( | ) | const |
Set the states of the phase objects to the locally-stored state within this MultiPhase object.
synchronize the phase objects with the mixture state.
Note that if individual phases have T and P different than that stored locally, the phase T and P will be modified.
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.
This method sets each phase to the mixture temperature and pressure, and sets the phase mole fractions based on the mixture mole numbers.
Definition at line 1164 of file MultiPhase.cpp.
References DATA_PTR, MultiPhase::m_moleFractions, MultiPhase::m_np, MultiPhase::m_phase, MultiPhase::m_press, MultiPhase::m_temp, MultiPhase::m_temp_OK, MultiPhase::maxTemp(), and MultiPhase::minTemp().
Referenced by MultiPhase::cp(), MultiPhase::enthalpy(), MultiPhase::entropy(), MultiPhase::getChemPotentials(), MultiPhase::getValidChemPotentials(), MultiPhase::gibbs(), MultiPhase::init(), MultiPhase::IntEnergy(), MultiPhase::setPressure(), MultiPhase::setState_TP(), and MultiPhase::setTemperature().
|
private |
Calculate the element abundance vector.
Definition at line 643 of file MultiPhase.cpp.
References MultiPhase::m_atoms, MultiPhase::m_elemAbundances, MultiPhase::m_moleFractions, MultiPhase::m_moles, MultiPhase::m_nel, MultiPhase::m_np, MultiPhase::m_phase, and Phase::nSpecies().
Referenced by MultiPhase::getElemAbundances(), and MultiPhase::uploadMoleFractionsFromPhases().
|
private |
Vector of the number of moles in each phase.
Length = m_np, number of phases.
Definition at line 612 of file MultiPhase.h.
Referenced by MultiPhase::addPhase(), MultiPhase::addPhases(), MultiPhase::calcElemAbundances(), MultiPhase::cp(), MultiPhase::elementMoles(), MultiPhase::enthalpy(), MultiPhase::entropy(), MultiPhase::getMoles(), MultiPhase::gibbs(), MultiPhase::IntEnergy(), MultiPhase::operator=(), MultiPhase::phaseCharge(), MultiPhase::phaseMoles(), MultiPhase::setMoles(), MultiPhase::setPhaseMoles(), MultiPhase::speciesMoles(), and MultiPhase::volume().
|
private |
Vector of the ThermoPhase Pointers.
Definition at line 617 of file MultiPhase.h.
Referenced by MultiPhase::addPhase(), MultiPhase::addPhases(), MultiPhase::calcElemAbundances(), MultiPhase::cp(), MultiPhase::elementMoles(), MultiPhase::enthalpy(), MultiPhase::entropy(), MultiPhase::getChemPotentials(), MultiPhase::getMoles(), MultiPhase::getValidChemPotentials(), MultiPhase::gibbs(), MultiPhase::init(), MultiPhase::IntEnergy(), MultiPhase::operator=(), MultiPhase::phase(), MultiPhase::phaseCharge(), MultiPhase::phaseIndex(), MultiPhase::phaseName(), MultiPhase::setMoles(), MultiPhase::setPhaseMoleFractions(), MultiPhase::solutionSpecies(), MultiPhase::speciesIndex(), MultiPhase::updatePhases(), MultiPhase::uploadMoleFractionsFromPhases(), and MultiPhase::volume().
|
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 626 of file MultiPhase.h.
Referenced by MultiPhase::calcElemAbundances(), MultiPhase::elementMoles(), MultiPhase::init(), MultiPhase::nAtoms(), and MultiPhase::operator=().
|
private |
Locally stored vector of mole fractions of all species comprising the MultiPhase object.
Definition at line 632 of file MultiPhase.h.
Referenced by MultiPhase::calcElemAbundances(), MultiPhase::elementMoles(), MultiPhase::getMoleFractions(), MultiPhase::getMoles(), MultiPhase::init(), MultiPhase::moleFraction(), MultiPhase::operator=(), MultiPhase::phase(), MultiPhase::phaseCharge(), MultiPhase::setMoles(), MultiPhase::setPhaseMoleFractions(), MultiPhase::speciesMoles(), MultiPhase::updatePhases(), and MultiPhase::uploadMoleFractionsFromPhases().
|
private |
Mapping between the global species number and the phase ID.
m_spphase[kGlobal] = iPhase Length = number of global species
Definition at line 639 of file MultiPhase.h.
Referenced by MultiPhase::init(), MultiPhase::operator=(), MultiPhase::solutionSpecies(), MultiPhase::speciesMoles(), and MultiPhase::speciesPhaseIndex().
|
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 647 of file MultiPhase.h.
Referenced by MultiPhase::init(), MultiPhase::operator=(), MultiPhase::phase(), MultiPhase::setPhaseMoleFractions(), and MultiPhase::speciesIndex().
|
private |
String names of the global elements.
This has a length equal to the number of global elements.
Definition at line 653 of file MultiPhase.h.
Referenced by MultiPhase::addPhase(), MultiPhase::elementIndex(), MultiPhase::elementName(), MultiPhase::init(), and MultiPhase::operator=().
|
private |
Atomic number of each element.
This is the atomic number of each global element.
Definition at line 659 of file MultiPhase.h.
Referenced by MultiPhase::addPhase(), and MultiPhase::init().
|
private |
Vector of species names in the problem.
Vector is over all species defined in the object, the global species index.
Definition at line 666 of file MultiPhase.h.
Referenced by MultiPhase::init(), and MultiPhase::speciesName().
|
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 672 of file MultiPhase.h.
Referenced by MultiPhase::addPhase(), and MultiPhase::operator=().
|
private |
Number of phases in the MultiPhase object.
Definition at line 677 of file MultiPhase.h.
Referenced by MultiPhase::addPhase(), MultiPhase::addPhases(), MultiPhase::calcElemAbundances(), MultiPhase::charge(), MultiPhase::cp(), MultiPhase::elementMoles(), MultiPhase::enthalpy(), MultiPhase::entropy(), MultiPhase::getChemPotentials(), MultiPhase::getMoles(), MultiPhase::getValidChemPotentials(), MultiPhase::gibbs(), MultiPhase::init(), MultiPhase::IntEnergy(), MultiPhase::nPhases(), MultiPhase::operator=(), MultiPhase::phaseIndex(), MultiPhase::setMoles(), MultiPhase::updatePhases(), MultiPhase::uploadMoleFractionsFromPhases(), and MultiPhase::volume().
|
private |
Current value of the temperature (kelvin)
Definition at line 680 of file MultiPhase.h.
Referenced by MultiPhase::addPhase(), MultiPhase::equilibrate(), MultiPhase::operator=(), MultiPhase::phase(), MultiPhase::setMoles(), MultiPhase::setPhaseMoleFractions(), MultiPhase::setState_TP(), MultiPhase::setState_TPMoles(), MultiPhase::setTemperature(), MultiPhase::temperature(), and MultiPhase::updatePhases().
|
private |
Current value of the pressure (Pa)
Definition at line 683 of file MultiPhase.h.
Referenced by MultiPhase::addPhase(), MultiPhase::operator=(), MultiPhase::phase(), MultiPhase::pressure(), MultiPhase::setMoles(), MultiPhase::setPhaseMoleFractions(), MultiPhase::setPressure(), MultiPhase::setState_TP(), MultiPhase::setState_TPMoles(), and MultiPhase::updatePhases().
|
private |
Number of distinct elements in all of the phases.
Definition at line 688 of file MultiPhase.h.
Referenced by MultiPhase::addPhase(), MultiPhase::calcElemAbundances(), MultiPhase::checkElementArraySize(), MultiPhase::checkElementIndex(), MultiPhase::elementIndex(), MultiPhase::getElemAbundances(), MultiPhase::init(), MultiPhase::nElements(), and MultiPhase::operator=().
|
private |
Number of distinct species in all of the phases.
Definition at line 692 of file MultiPhase.h.
Referenced by MultiPhase::addPhase(), MultiPhase::addSpeciesMoles(), MultiPhase::checkSpeciesArraySize(), MultiPhase::checkSpeciesIndex(), MultiPhase::init(), MultiPhase::nSpecies(), and MultiPhase::operator=().
|
private |
True if the init() routine has been called, and the MultiPhase frozen.
Definition at line 695 of file MultiPhase.h.
Referenced by MultiPhase::addPhase(), MultiPhase::equilibrate(), MultiPhase::init(), MultiPhase::operator=(), MultiPhase::phase(), MultiPhase::setMoles(), MultiPhase::setPhaseMoleFractions(), MultiPhase::setState_TP(), MultiPhase::setTemperature(), and MultiPhase::speciesIndex().
|
private |
Global ID of the element corresponding to the electronic charge.
If there is none, then this is equal to -1
Definition at line 701 of file MultiPhase.h.
Referenced by MultiPhase::addPhase(), MultiPhase::init(), and MultiPhase::operator=().
|
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 708 of file MultiPhase.h.
Referenced by MultiPhase::addPhase(), MultiPhase::operator=(), MultiPhase::tempOK(), and MultiPhase::updatePhases().
|
private |
Minimum temperature for which thermo parameterizations are valid.
Stoichiometric phases are ignored in this determination. units Kelvin
Definition at line 715 of file MultiPhase.h.
Referenced by MultiPhase::addPhase(), MultiPhase::equilibrate(), MultiPhase::minTemp(), and MultiPhase::operator=().
|
private |
Minimum temperature for which thermo parameterizations are valid.
Stoichiometric phases are ignored in this determination. units Kelvin
Definition at line 722 of file MultiPhase.h.
Referenced by MultiPhase::addPhase(), MultiPhase::equilibrate(), MultiPhase::maxTemp(), and MultiPhase::operator=().
|
mutableprivate |
Vector of element abundances.
m_elemAbundances[mGlobal] = kmol of element mGlobal summed over all species in all phases.
Definition at line 729 of file MultiPhase.h.
Referenced by MultiPhase::calcElemAbundances(), MultiPhase::getElemAbundances(), MultiPhase::init(), and MultiPhase::operator=().