Cantera 2.6.0
|
Class Phase is the base class for phases of matter, managing the species and elements in a phase, as well as the independent variables of temperature, mass density (compressible substances) or pressure (incompressible substances), species mass/mole fraction, and other generalized forces and intrinsic properties (such as electric potential) that define the thermodynamic state. More...
#include <Phase.h>
Public Member Functions | |
Phase () | |
Default constructor. More... | |
Phase (const Phase &)=delete | |
Phase & | operator= (const Phase &)=delete |
XML_Node & | xml () const |
Returns a const reference to the XML_Node that describes the phase. More... | |
void | setXMLdata (XML_Node &xmlPhase) |
Stores the XML tree information for the current phase. More... | |
virtual bool | isPure () const |
Return whether phase represents a pure (single species) substance. More... | |
virtual bool | hasPhaseTransition () const |
Return whether phase represents a substance with phase transitions. More... | |
virtual bool | isCompressible () const |
Return whether phase represents a compressible substance. More... | |
virtual std::map< std::string, size_t > | nativeState () const |
Return a map of properties defining the native state of a substance. More... | |
virtual std::vector< std::string > | fullStates () const |
Return a vector containing full states defining a phase. More... | |
virtual std::vector< std::string > | partialStates () const |
Return a vector of settable partial property sets within a phase. More... | |
virtual size_t | stateSize () const |
Return size of vector defining internal state of the phase. More... | |
void | saveState (vector_fp &state) const |
Save the current internal state of the phase. More... | |
virtual void | saveState (size_t lenstate, doublereal *state) const |
Write to array 'state' the current internal state. More... | |
void | restoreState (const vector_fp &state) |
Restore a state saved on a previous call to saveState. More... | |
virtual void | restoreState (size_t lenstate, const doublereal *state) |
Restore the state of the phase from a previously saved state vector. More... | |
doublereal | molecularWeight (size_t k) const |
Molecular weight of species k . More... | |
void | getMolecularWeights (vector_fp &weights) const |
Copy the vector of molecular weights into vector weights. More... | |
void | getMolecularWeights (doublereal *weights) const |
Copy the vector of molecular weights into array weights. More... | |
const vector_fp & | molecularWeights () const |
Return a const reference to the internal vector of molecular weights. More... | |
void | getCharges (double *charges) const |
Copy the vector of species charges into array charges. More... | |
doublereal | elementalMassFraction (const size_t m) const |
Elemental mass fraction of element m. More... | |
doublereal | elementalMoleFraction (const size_t m) const |
Elemental mole fraction of element m. More... | |
const double * | moleFractdivMMW () const |
Returns a const pointer to the start of the moleFraction/MW array. More... | |
doublereal | charge (size_t k) const |
Dimensionless electrical charge of a single molecule of species k The charge is normalized by the the magnitude of the electron charge. More... | |
doublereal | chargeDensity () const |
Charge density [C/m^3]. More... | |
size_t | nDim () const |
Returns the number of spatial dimensions (1, 2, or 3) More... | |
void | setNDim (size_t ndim) |
Set the number of spatial dimensions (1, 2, or 3). More... | |
virtual bool | ready () const |
Returns a bool indicating whether the object is ready for use. More... | |
int | stateMFNumber () const |
Return the State Mole Fraction Number. More... | |
virtual void | invalidateCache () |
Invalidate any cached values which are normally updated only when a change in state is detected. More... | |
bool | caseSensitiveSpecies () const |
Returns true if case sensitive species names are enforced. More... | |
void | setCaseSensitiveSpecies (bool cflag=true) |
Set flag that determines whether case sensitive species are enforced in look-up operations, for example speciesIndex. More... | |
virtual void | setRoot (std::shared_ptr< Solution > root) |
Set root Solution holding all phase information. More... | |
vector_fp | getCompositionFromMap (const compositionMap &comp) const |
Converts a compositionMap to a vector with entries for each species Species that are not specified are set to zero in the vector. More... | |
void | massFractionsToMoleFractions (const double *Y, double *X) const |
Converts a mixture composition from mole fractions to mass fractions. More... | |
void | moleFractionsToMassFractions (const double *X, double *Y) const |
Converts a mixture composition from mass fractions to mole fractions. More... | |
Name | |
Class Phase uses the string name to identify a phase. The name is the value of the corresponding key in the phase map (in YAML), name (in CTI), or id (in XML) that is used to initialize a phase when it is read. However, the name field may be changed to another value during the course of a calculation. For example, if duplicates of a phase object are instantiated and used in multiple places (such as a ReactorNet), they will have the same constitutive input, that is, the names of the phases will be the same. Note that this is not a problem for Cantera internally; however, a user may want to rename phase objects in order to clarify. | |
std::string | name () const |
Return the name of the phase. More... | |
void | setName (const std::string &nm) |
Sets the string name for the phase. More... | |
virtual std::string | type () const |
String indicating the thermodynamic model implemented. More... | |
Element and Species Information | |
std::string | elementName (size_t m) const |
Name of the element with index m. More... | |
size_t | elementIndex (const std::string &name) const |
Return the index of element named 'name'. More... | |
const std::vector< std::string > & | elementNames () const |
Return a read-only reference to the vector of element names. More... | |
doublereal | atomicWeight (size_t m) const |
Atomic weight of element m. More... | |
doublereal | entropyElement298 (size_t m) const |
Entropy of the element in its standard state at 298 K and 1 bar. More... | |
int | atomicNumber (size_t m) const |
Atomic number of element m. More... | |
int | elementType (size_t m) const |
Return the element constraint type Possible types include: More... | |
int | changeElementType (int m, int elem_type) |
Change the element type of the mth constraint Reassigns an element type. More... | |
const vector_fp & | atomicWeights () const |
Return a read-only reference to the vector of atomic weights. More... | |
size_t | nElements () const |
Number of elements. More... | |
void | checkElementIndex (size_t m) const |
Check that the specified element index is in range. More... | |
void | checkElementArraySize (size_t mm) const |
Check that an array size is at least nElements(). More... | |
doublereal | nAtoms (size_t k, size_t m) const |
Number of atoms of element m in species k . More... | |
void | getAtoms (size_t k, double *atomArray) const |
Get a vector containing the atomic composition of species k. More... | |
size_t | speciesIndex (const std::string &name) const |
Returns the index of a species named 'name' within the Phase object. More... | |
std::string | speciesName (size_t k) const |
Name of the species with index k. More... | |
std::string | speciesSPName (int k) const |
Returns the expanded species name of a species, including the phase name This is guaranteed to be unique within a Cantera problem. More... | |
const std::vector< std::string > & | speciesNames () const |
Return a const reference to the vector of species names. More... | |
size_t | nSpecies () const |
Returns the number of species in the phase. More... | |
void | checkSpeciesIndex (size_t k) const |
Check that the specified species index is in range. More... | |
void | checkSpeciesArraySize (size_t kk) const |
Check that an array size is at least nSpecies(). More... | |
Set thermodynamic state | |
Set the internal thermodynamic state by setting the internally stored temperature, density and species composition. Note that the composition is always set first. Temperature and density are held constant if not explicitly set. | |
void | setMoleFractionsByName (const compositionMap &xMap) |
Set the species mole fractions by name. More... | |
void | setMoleFractionsByName (const std::string &x) |
Set the mole fractions of a group of species by name. More... | |
void | setMassFractionsByName (const compositionMap &yMap) |
Set the species mass fractions by name. More... | |
void | setMassFractionsByName (const std::string &x) |
Set the species mass fractions by name. More... | |
void | setState_TRX (doublereal t, doublereal dens, const doublereal *x) |
Set the internally stored temperature (K), density, and mole fractions. More... | |
void | setState_TRX (doublereal t, doublereal dens, const compositionMap &x) |
Set the internally stored temperature (K), density, and mole fractions. More... | |
void | setState_TRY (doublereal t, doublereal dens, const doublereal *y) |
Set the internally stored temperature (K), density, and mass fractions. More... | |
void | setState_TRY (doublereal t, doublereal dens, const compositionMap &y) |
Set the internally stored temperature (K), density, and mass fractions. More... | |
void | setState_TNX (doublereal t, doublereal n, const doublereal *x) |
Set the internally stored temperature (K), molar density (kmol/m^3), and mole fractions. More... | |
void | setState_TR (doublereal t, doublereal rho) |
Set the internally stored temperature (K) and density (kg/m^3) More... | |
void | setState_TX (doublereal t, doublereal *x) |
Set the internally stored temperature (K) and mole fractions. More... | |
void | setState_TY (doublereal t, doublereal *y) |
Set the internally stored temperature (K) and mass fractions. More... | |
void | setState_RX (doublereal rho, doublereal *x) |
Set the density (kg/m^3) and mole fractions. More... | |
void | setState_RY (doublereal rho, doublereal *y) |
Set the density (kg/m^3) and mass fractions. More... | |
Composition | |
compositionMap | getMoleFractionsByName (double threshold=0.0) const |
Get the mole fractions by name. More... | |
double | moleFraction (size_t k) const |
Return the mole fraction of a single species. More... | |
double | moleFraction (const std::string &name) const |
Return the mole fraction of a single species. More... | |
compositionMap | getMassFractionsByName (double threshold=0.0) const |
Get the mass fractions by name. More... | |
double | massFraction (size_t k) const |
Return the mass fraction of a single species. More... | |
double | massFraction (const std::string &name) const |
Return the mass fraction of a single species. More... | |
void | getMoleFractions (double *const x) const |
Get the species mole fraction vector. More... | |
virtual void | setMoleFractions (const double *const x) |
Set the mole fractions to the specified values. More... | |
virtual void | setMoleFractions_NoNorm (const double *const x) |
Set the mole fractions to the specified values without normalizing. More... | |
void | getMassFractions (double *const y) const |
Get the species mass fractions. More... | |
const double * | massFractions () const |
Return a const pointer to the mass fraction array. More... | |
virtual void | setMassFractions (const double *const y) |
Set the mass fractions to the specified values and normalize them. More... | |
virtual void | setMassFractions_NoNorm (const double *const y) |
Set the mass fractions to the specified values without normalizing. More... | |
void | getConcentrations (double *const c) const |
Get the species concentrations (kmol/m^3). More... | |
double | concentration (const size_t k) const |
Concentration of species k. More... | |
virtual void | setConcentrations (const double *const conc) |
Set the concentrations to the specified values within the phase. More... | |
virtual void | setConcentrationsNoNorm (const double *const conc) |
Set the concentrations without ignoring negative concentrations. More... | |
Thermodynamic Properties | |
doublereal | temperature () const |
Temperature (K). More... | |
virtual double | electronTemperature () const |
Electron Temperature (K) More... | |
virtual double | pressure () const |
Return the thermodynamic pressure (Pa). More... | |
virtual double | density () const |
Density (kg/m^3). More... | |
double | molarDensity () const |
Molar density (kmol/m^3). More... | |
double | molarVolume () const |
Molar volume (m^3/kmol). More... | |
virtual void | setDensity (const double density_) |
Set the internally stored density (kg/m^3) of the phase. More... | |
virtual void | setMolarDensity (const double molarDensity) |
Set the internally stored molar density (kmol/m^3) of the phase. More... | |
virtual void | setPressure (double p) |
Set the internally stored pressure (Pa) at constant temperature and composition. More... | |
virtual void | setTemperature (double temp) |
Set the internally stored temperature of the phase (K). More... | |
virtual void | setElectronTemperature (double etemp) |
Set the internally stored electron temperature of the phase (K). More... | |
Mean Properties | |
doublereal | mean_X (const doublereal *const Q) const |
Evaluate the mole-fraction-weighted mean of an array Q. More... | |
doublereal | mean_X (const vector_fp &Q) const |
Evaluate the mole-fraction-weighted mean of an array Q. More... | |
doublereal | meanMolecularWeight () const |
The mean molecular weight. Units: (kg/kmol) More... | |
doublereal | sum_xlogx () const |
Evaluate \( \sum_k X_k \log X_k \). More... | |
Adding Elements and Species | |
These methods are used to add new elements or species. These are not usually called by user programs. Since species are checked to insure that they are only composed of declared elements, it is necessary to first add all elements before adding any species. | |
size_t | addElement (const std::string &symbol, doublereal weight=-12345.0, int atomicNumber=0, doublereal entropy298=ENTROPY298_UNKNOWN, int elem_type=CT_ELEM_TYPE_ABSPOS) |
Add an element. More... | |
virtual bool | addSpecies (shared_ptr< Species > spec) |
Add a Species to this Phase. More... | |
virtual void | modifySpecies (size_t k, shared_ptr< Species > spec) |
Modify the thermodynamic data associated with a species. More... | |
void | addSpeciesAlias (const std::string &name, const std::string &alias) |
Add a species alias (that is, a user-defined alternative species name). More... | |
virtual std::vector< std::string > | findIsomers (const compositionMap &compMap) const |
Return a vector with isomers names matching a given composition map. More... | |
virtual std::vector< std::string > | findIsomers (const std::string &comp) const |
Return a vector with isomers names matching a given composition string. More... | |
shared_ptr< Species > | species (const std::string &name) const |
Return the Species object for the named species. More... | |
shared_ptr< Species > | species (size_t k) const |
Return the Species object for species whose index is k. More... | |
void | ignoreUndefinedElements () |
Set behavior when adding a species containing undefined elements to just skip the species. More... | |
void | addUndefinedElements () |
Set behavior when adding a species containing undefined elements to add those elements to the phase. More... | |
void | throwUndefinedElements () |
Set the behavior when adding a species containing undefined elements to throw an exception. More... | |
Protected Member Functions | |
void | assertCompressible (const std::string &setter) const |
Ensure that phase is compressible. More... | |
void | assignDensity (const double density_) |
Set the internally stored constant density (kg/m^3) of the phase. More... | |
void | setMolecularWeight (const int k, const double mw) |
Set the molecular weight of a single species to a given value. More... | |
virtual void | compositionChanged () |
Apply changes to the state which are needed after the composition changes. More... | |
Protected Attributes | |
ValueCache | m_cache |
Cached for saved calculations within each ThermoPhase. More... | |
size_t | m_kk |
Number of species in the phase. More... | |
size_t | m_ndim |
Dimensionality of the phase. More... | |
vector_fp | m_speciesComp |
Atomic composition of the species. More... | |
vector_fp | m_speciesCharge |
Vector of species charges. length m_kk. More... | |
std::map< std::string, shared_ptr< Species > > | m_species |
UndefElement::behavior | m_undefinedElementBehavior |
Flag determining behavior when adding species with an undefined element. More... | |
bool | m_caseSensitiveSpecies |
Flag determining whether case sensitive species names are enforced. More... | |
Private Member Functions | |
size_t | findSpeciesLower (const std::string &nameStr) const |
Find lowercase species name in m_speciesIndices when case sensitive species names are not enforced and a user specifies a non-canonical species name. More... | |
Private Attributes | |
XML_Node * | m_xml |
XML node containing the XML info for this phase. More... | |
std::string | m_id |
ID of the phase. More... | |
std::string | m_name |
Name of the phase. More... | |
doublereal | m_temp |
Temperature (K). This is an independent variable. More... | |
doublereal | m_dens |
Density (kg m-3). More... | |
doublereal | m_mmw |
mean molecular weight of the mixture (kg kmol-1) More... | |
vector_fp | m_ym |
m_ym[k] = mole fraction of species k divided by the mean molecular weight of mixture. More... | |
vector_fp | m_y |
Mass fractions of the species. More... | |
vector_fp | m_molwts |
species molecular weights (kg kmol-1) More... | |
vector_fp | m_rmolwts |
inverse of species molecular weights (kmol kg-1) More... | |
int | m_stateNum |
State Change variable. More... | |
std::vector< std::string > | m_speciesNames |
Vector of the species names. More... | |
std::map< std::string, size_t > | m_speciesIndices |
Map of species names to indices. More... | |
std::map< std::string, size_t > | m_speciesLower |
Map of lower-case species names to indices. More... | |
size_t | m_mm |
Number of elements. More... | |
vector_fp | m_atomicWeights |
element atomic weights (kg kmol-1) More... | |
vector_int | m_atomicNumbers |
element atomic numbers More... | |
std::vector< std::string > | m_elementNames |
element names More... | |
vector_int | m_elem_type |
Vector of element types. More... | |
vector_fp | m_entropy298 |
Entropy at 298.15 K and 1 bar of stable state pure elements (J kmol-1) More... | |
Class Phase is the base class for phases of matter, managing the species and elements in a phase, as well as the independent variables of temperature, mass density (compressible substances) or pressure (incompressible substances), species mass/mole fraction, and other generalized forces and intrinsic properties (such as electric potential) that define the thermodynamic state.
Class Phase provides information about the elements and species in a phase - names, index numbers (location in arrays), atomic or molecular weights, etc. The set of elements must include all those that compose the species, but may include additional elements.
It also stores an array of species molecular weights, which are used to convert between mole and mass representations of the composition. For efficiency in mass/mole conversion, the vector of mass fractions divided by molecular weight \( Y_k/M_k \) is also stored.
Class Phase is not usually used directly. Its primary use is as a base class for class ThermoPhase. It is not generally necessary to overloaded any of class Phase's methods, which handles both compressible and incompressible phases. For incompressible phases, the density is replaced by the pressure as the independent variable, and can no longer be set directly. In this case, the density needs to be calculated from a suitable equation of state, and assigned to the object using the assignDensity() method. This also applies for nearly-incompressible phases or phases which utilize standard states based on a T and P, in which case they need to overload these functions too.
Class Phase contains a number of utility functions that will set the state of the phase in its entirety, by first setting the composition, then the temperature and then density (or pressure for incompressible substances) An example of this is the function Phase::setState_TRY(double t, double dens, const double* y).
Class Phase contains methods for saving and restoring the full internal state of a given phase. These are saveState() and restoreState(). These functions operate on a state vector, which by default uses the first two entries for temperature and density (compressible substances) or temperature and pressure (incompressible substances). If the substance is not pure in a thermodynamic sense (that is, it may contain multiple species), the state also contains nSpecies() entries that specify the composition by corresponding mass fractions. Default definitions can be overloaded by derived classes. For any phase, the native definition of its thermodynamic state is defined the method nativeState(), with the length of the state vector returned by by stateSize(). In addition, methods isPure() and isCompressible() provide information on the implementation of a Phase object.
A species name is referred to via speciesName(), which is unique within a given phase. Note that within multiphase mixtures (MultiPhase()), both a phase name/index as well as species name are required to access information about a species in a particular phase. For surfaces, the species names are unique among the phases.
XML_Node & xml | ( | ) | const |
Returns a const reference to the XML_Node that describes the phase.
The XML_Node for the phase contains all of the input data used to set up the model for the phase during its initialization.
Definition at line 46 of file Phase.cpp.
References Phase::m_xml.
Referenced by TransportFactory::newTransport().
void setXMLdata | ( | XML_Node & | xmlPhase | ) |
Stores the XML tree information for the current phase.
This function now stores the complete XML_Node tree as read into the code via a file. This is needed to move around within the XML tree during construction of transport and kinetics mechanisms after copy construction operations.
xmlPhase | Reference to the XML node corresponding to the phase |
Definition at line 51 of file Phase.cpp.
References XML_Node::copy(), Cantera::findXMLPhase(), XML_Node::id(), Phase::m_xml, and XML_Node::root().
Referenced by Cantera::importPhase().
std::string name | ( | ) | const |
Return the name of the phase.
Names are unique within a Cantera problem.
Definition at line 70 of file Phase.cpp.
References Phase::m_name.
Referenced by MultiPhase::addPhase(), Phase::addSpeciesAlias(), Phase::assertCompressible(), Cantera::checkElectrochemReaction(), Phase::findSpeciesLower(), ThermoPhase::getParameters(), DebyeHuckel::getSpeciesParameters(), IdealSolidSolnPhase::getSpeciesParameters(), LatticePhase::getSpeciesParameters(), LatticeSolidPhase::getSpeciesParameters(), PengRobinson::getSpeciesParameters(), RedlichKwongMFTP::getSpeciesParameters(), StoichSubstance::getSpeciesParameters(), VPStandardStateTP::getSpeciesParameters(), Cantera::importPhase(), BinarySolutionTabulatedThermo::initThermo(), Kinetics::kineticsSpeciesIndex(), Cantera::operator<<(), ThermoPhase::report(), MolalityVPSSTP::report(), PureFluidPhase::report(), vcs_MultiPhaseEquil::reportCSV(), ReactingSurf1D::restore(), ReactingSurf1D::serialize(), StFlow::serialize(), LatticeSolidPhase::setLatticeStoichiometry(), Phase::setName(), PureFluidPhase::setSubstance(), Phase::species(), vcs_VolPhase::transferElementsFM(), and VCS_SOLVE::VCS_SOLVE().
void setName | ( | const std::string & | nm | ) |
Sets the string name for the phase.
nm | String name of the phase |
Definition at line 75 of file Phase.cpp.
References Phase::m_id, Phase::m_name, and Phase::name().
Referenced by Cantera::importPhase(), and Cantera::setupPhase().
|
inlinevirtual |
String indicating the thermodynamic model implemented.
Usually corresponds to the name of the derived class, less any suffixes such as "Phase", TP", "VPSS", etc.
Reimplemented in BinarySolutionTabulatedThermo, DebyeHuckel, EdgePhase, HMWSoln, IdealGasPhase, IdealMolalSoln, IdealSolidSolnPhase, IdealSolnGasVPSS, IonsFromNeutralVPSSTP, LatticePhase, LatticeSolidPhase, MargulesVPSSTP, MaskellSolidSolnPhase, MetalPhase, MixtureFugacityTP, PengRobinson, PlasmaPhase, PureFluidPhase, RedlichKisterVPSSTP, RedlichKwongMFTP, SingleSpeciesTP, StoichSubstance, SurfPhase, ThermoPhase, and WaterSSTP.
Definition at line 164 of file Phase.h.
Referenced by Phase::assertCompressible(), Phase::pressure(), Phase::setElectronTemperature(), and Phase::setPressure().
string elementName | ( | size_t | m | ) | const |
Name of the element with index m.
m | Element index. |
Definition at line 100 of file Phase.cpp.
References Phase::checkElementIndex(), and Phase::m_elementNames.
Referenced by MultiPhase::addPhase(), Reaction::checkBalance(), PDSS_HKFT::convertDGFormation(), Phase::elementIndex(), ChemEquil::equilibrate(), ChemEquil::estimateElementPotentials(), MolalityVPSSTP::findCLMIndex(), ThermoPhase::getParameters(), ChemEquil::initialize(), ChemEquil::setInitialMoles(), and vcs_VolPhase::transferElementsFM().
size_t elementIndex | ( | const std::string & | name | ) | const |
Return the index of element named 'name'.
The index is an integer assigned to each element in the order it was added. Returns npos if the specified element is not found.
name | Name of the element |
Definition at line 106 of file Phase.cpp.
References Phase::elementName(), Phase::m_elementNames, Phase::m_mm, and Cantera::npos.
Referenced by Phase::addSpecies(), MultiPhase::init(), WaterSSTP::initThermo(), PDSS_HKFT::LookupGe(), ThermoPhase::o2Present(), and ThermoPhase::o2Required().
const vector< string > & elementNames | ( | ) | const |
Return a read-only reference to the vector of element names.
Definition at line 116 of file Phase.cpp.
References Phase::m_elementNames.
Referenced by ThermoPhase::getParameters(), and IonsFromNeutralVPSSTP::initThermo().
doublereal atomicWeight | ( | size_t | m | ) | const |
Atomic weight of element m.
m | Element index |
Definition at line 121 of file Phase.cpp.
References Phase::m_atomicWeights.
Referenced by Phase::elementalMassFraction(), ChemEquil::initialize(), and WaterSSTP::initThermo().
doublereal entropyElement298 | ( | size_t | m | ) | const |
Entropy of the element in its standard state at 298 K and 1 bar.
If no entropy value was provided when the phase was constructed, returns the value ENTROPY298_UNKNOWN
.
m | Element index |
Definition at line 126 of file Phase.cpp.
References Phase::checkElementIndex(), and Phase::m_entropy298.
Referenced by PDSS_HKFT::LookupGe().
int atomicNumber | ( | size_t | m | ) | const |
Atomic number of element m.
m | Element index |
Definition at line 137 of file Phase.cpp.
References Phase::m_atomicNumbers.
Referenced by MultiPhase::addPhase().
int elementType | ( | size_t | m | ) | const |
Return the element constraint type Possible types include:
CT_ELEM_TYPE_TURNEDOFF -1
CT_ELEM_TYPE_ABSPOS 0
CT_ELEM_TYPE_ELECTRONCHARGE 1
CT_ELEM_TYPE_CHARGENEUTRALITY 2
CT_ELEM_TYPE_LATTICERATIO 3
CT_ELEM_TYPE_KINETICFROZEN 4
CT_ELEM_TYPE_SURFACECONSTRAINT 5
CT_ELEM_TYPE_OTHERCONSTRAINT 6
The default is CT_ELEM_TYPE_ABSPOS
.
m | Element index |
Definition at line 142 of file Phase.cpp.
References Phase::m_elem_type.
Referenced by vcs_VolPhase::transferElementsFM().
int changeElementType | ( | int | m, |
int | elem_type | ||
) |
Change the element type of the mth constraint Reassigns an element type.
m | Element index |
elem_type | New elem type to be assigned |
Definition at line 147 of file Phase.cpp.
References Phase::m_elem_type.
const vector_fp & atomicWeights | ( | ) | const |
Return a read-only reference to the vector of atomic weights.
Definition at line 132 of file Phase.cpp.
References Phase::m_atomicWeights.
Referenced by Phase::addSpecies().
size_t nElements | ( | ) | const |
Number of elements.
Definition at line 81 of file Phase.cpp.
References Phase::m_mm.
Referenced by MultiPhase::addPhase(), Phase::addSpecies(), Reaction::checkBalance(), PDSS_HKFT::convertDGFormation(), Phase::elementalMoleFraction(), ChemEquil::equilibrate(), MolalityVPSSTP::findCLMIndex(), ThermoPhase::getParameters(), ChemEquil::initialize(), IonsFromNeutralVPSSTP::initThermo(), LatticeSolidPhase::setLatticeStoichiometry(), vcs_VolPhase::transferElementsFM(), and VCS_SOLVE::VCS_SOLVE().
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 86 of file Phase.cpp.
References Phase::m_mm.
Referenced by Phase::elementalMassFraction(), Phase::elementalMoleFraction(), Phase::elementName(), Phase::entropyElement298(), and Phase::nAtoms().
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 93 of file Phase.cpp.
References Phase::m_mm.
doublereal nAtoms | ( | size_t | k, |
size_t | m | ||
) | const |
Number of atoms of element m
in species k
.
k | species index |
m | element index |
Definition at line 154 of file Phase.cpp.
References Phase::checkElementIndex(), Phase::checkSpeciesIndex(), Phase::m_mm, and Phase::m_speciesComp.
Referenced by Reaction::checkBalance(), PDSS_HKFT::convertDGFormation(), Phase::elementalMassFraction(), Phase::elementalMoleFraction(), MolalityVPSSTP::findCLMIndex(), MultiPhase::init(), ChemEquil::initialize(), IonsFromNeutralVPSSTP::initThermo(), ThermoPhase::o2Present(), ThermoPhase::o2Required(), and vcs_VolPhase::transferElementsFM().
void getAtoms | ( | size_t | k, |
double * | atomArray | ||
) | const |
Get a vector containing the atomic composition of species k.
k | species index |
atomArray | vector containing the atomic number in the species. Length: m_mm |
Definition at line 161 of file Phase.cpp.
References Phase::m_mm, and Phase::m_speciesComp.
size_t speciesIndex | ( | const std::string & | name | ) | const |
Returns the index of a species named 'name' within the Phase object.
The first species in the phase will have an index 0, and the last one will have an index of nSpecies() - 1.
name | String name of the species. It may also be in the form phaseName:speciesName |
Definition at line 187 of file Phase.cpp.
References Phase::findSpeciesLower(), Phase::m_caseSensitiveSpecies, Phase::m_speciesIndices, and Cantera::npos.
Referenced by RedlichKisterVPSSTP::addBinaryInteraction(), Phase::addSpeciesAlias(), InterfaceKinetics::buildSurfaceArrhenius(), Reaction::checkBalance(), Cantera::checkElectrochemReaction(), Phase::getCompositionFromMap(), DebyeHuckel::getSpeciesParameters(), IdealSolidSolnPhase::getSpeciesParameters(), LatticePhase::getSpeciesParameters(), PengRobinson::getSpeciesParameters(), RedlichKwongMFTP::getSpeciesParameters(), StoichSubstance::getSpeciesParameters(), VPStandardStateTP::getSpeciesParameters(), BinarySolutionTabulatedThermo::initThermo(), PengRobinson::initThermo(), RedlichKwongMFTP::initThermo(), RedlichKwongMFTP::initThermoXML(), BinarySolutionTabulatedThermo::initThermoXML(), FlowDevice::install(), Kinetics::kineticsSpeciesIndex(), Phase::massFraction(), Phase::moleFraction(), HMWSoln::readXMLBinarySalt(), RedlichKisterVPSSTP::readXMLBinarySpecies(), HMWSoln::readXMLLambdaNeutral(), HMWSoln::readXMLMunnnNeutral(), HMWSoln::readXMLPsi(), RedlichKwongMFTP::readXMLPureFluid(), HMWSoln::readXMLTheta(), HMWSoln::readXMLZetaCation(), OutletRes1D::restore(), ReactingSurf1D::restore(), StFlow::restore(), DebyeHuckel::setBeta(), PengRobinson::setBinaryCoeffs(), RedlichKwongMFTP::setBinaryCoeffs(), StickingCoverage::setContext(), PengRobinson::setSpeciesCoeffs(), RedlichKwongMFTP::setSpeciesCoeffs(), Phase::species(), Kinetics::speciesPhase(), StFlow::StFlow(), and Reaction::usesElectrochemistry().
string speciesName | ( | size_t | k | ) | const |
Name of the species with index k.
k | index of the species |
Definition at line 200 of file Phase.cpp.
References Phase::checkSpeciesIndex(), and Phase::m_speciesNames.
Referenced by ConstPressureReactor::componentName(), ReactingSurf1D::componentName(), StFlow::componentName(), ChemEquil::estimateElementPotentials(), MolalityVPSSTP::findCLMIndex(), Phase::getMassFractionsByName(), Phase::getMoleFractionsByName(), BinarySolutionTabulatedThermo::getParameters(), MargulesVPSSTP::getParameters(), RedlichKisterVPSSTP::getParameters(), ThermoPhase::getParameters(), MultiPhase::init(), ChemEquil::initialize(), PDSS_HKFT::initThermo(), StoichSubstance::initThermo(), RedlichKwongMFTP::initThermoXML(), FlowDevice::install(), Kinetics::kineticsSpeciesName(), Phase::modifySpecies(), ThermoPhase::modifySpecies(), HMWSoln::printCoeffs(), vcs_MultiPhaseEquil::reportCSV(), ThermoPhase::reportCSV(), ReactingSurf1D::restore(), StFlow::restore(), OutletRes1D::save(), ReactingSurf1D::save(), StFlow::save(), OutletRes1D::serialize(), ReactingSurf1D::serialize(), StFlow::serialize(), SurfPhase::setCoveragesByName(), ChemEquil::setInitialMoles(), MolalityVPSSTP::setMolalitiesByName(), Transport::setThermo(), ReactingSurf1D::showSolution(), Phase::species(), Phase::speciesSPName(), and ChemEquil::update().
std::string speciesSPName | ( | int | k | ) | const |
Returns the expanded species name of a species, including the phase name This is guaranteed to be unique within a Cantera problem.
k | Species index within the phase |
Definition at line 225 of file Phase.cpp.
References Phase::m_name, and Phase::speciesName().
const vector< string > & speciesNames | ( | ) | const |
Return a const reference to the vector of species names.
Definition at line 206 of file Phase.cpp.
References Phase::m_speciesNames.
Referenced by ThermoPhase::getParameters(), InterfaceRateBase::setContext(), SurfPhase::setCoveragesByName(), and MolalityVPSSTP::setMolalitiesByName().
|
inline |
Returns the number of species in the phase.
Definition at line 273 of file Phase.h.
References Phase::m_kk.
Referenced by MultiPhase::addPhase(), InterfaceKinetics::applyVoltageKfwdCorrection(), MultiPhase::calcElemAbundances(), ConstPressureReactor::componentName(), MultiPhaseEquil::computeReactionSteps(), MargulesVPSSTP::cp_mole(), MargulesVPSSTP::enthalpy_mole(), MargulesVPSSTP::entropy_mole(), ChemEquil::equilibrate(), ChemEquil::estimateElementPotentials(), ThermoPhase::getActivities(), MetalPhase::getActivityConcentrations(), HighPressureGasTransport::getBinaryDiffCoeffs(), MetalPhase::getChemPotentials(), ImplicitSurfChem::getConcSpecies(), MolalityVPSSTP::getCsvReportData(), ThermoPhase::getCsvReportData(), MetalPhase::getEnthalpy_RT(), MetalPhase::getEntropy_R(), MultiTransport::getMassFluxes(), MultiTransport::getMolarFluxes(), MultiPhase::getMoles(), HighPressureGasTransport::getMultiDiffCoeffs(), MetalPhase::getPartialMolarEnthalpies(), RedlichKwongMFTP::getPartialMolarIntEnergies(), MetalPhase::getStandardChemPotentials(), Cantera::hasChargedSpecies(), MultiPhase::init(), OutletRes1D::init(), IonGasTransport::init(), GasTransport::init(), ChemEquil::initialize(), DustyGasTransport::initialize(), BinarySolutionTabulatedThermo::initThermo(), IdealMolalSoln::initThermo(), BinarySolutionTabulatedThermo::initThermoXML(), FlowDevice::install(), vcs_MultiPhaseEquil::reportCSV(), ThermoPhase::reportCSV(), ThermoPhase::resetHf298(), StFlow::restore(), IonsFromNeutralVPSSTP::s_update_dlnActCoeff_dlnN(), Kinetics::selectPhase(), ImplicitSurfChem::setConcSpecies(), MultiPhase::setMoles(), MultiPhase::setPhaseMoleFractions(), Transport::setThermo(), ReactorBase::setThermoMgr(), Cantera::setupPhase(), solveSP::solveSP(), Phase::stateSize(), StFlow::StFlow(), HighPressureGasTransport::thermalConductivity(), vcs_VolPhase::transferElementsFM(), MixTransport::update_T(), MultiTransport::update_T(), InterfaceKinetics::updateExchangeCurrentQuantities(), InterfaceKinetics::updateMu0(), MultiTransport::updateThermal_T(), MultiPhase::uploadMoleFractionsFromPhases(), VCS_SOLVE::vcs_prob_specifyFully(), VCS_SOLVE::VCS_SOLVE(), and HighPressureGasTransport::viscosity().
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 211 of file Phase.cpp.
References Phase::m_kk.
Referenced by Phase::concentration(), DebyeHuckel::getSpeciesParameters(), PengRobinson::getSpeciesParameters(), RedlichKwongMFTP::getSpeciesParameters(), Phase::massFraction(), Phase::molecularWeight(), Phase::moleFraction(), Phase::nAtoms(), Phase::species(), and Phase::speciesName().
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 218 of file Phase.cpp.
References Phase::m_kk.
|
inlinevirtual |
Return whether phase represents a pure (single species) substance.
Reimplemented in PureFluidPhase, and SingleSpeciesTP.
Definition at line 289 of file Phase.h.
Referenced by Phase::fullStates(), Phase::nativeState(), Phase::partialStates(), and Phase::stateSize().
|
inlinevirtual |
Return whether phase represents a substance with phase transitions.
Reimplemented in PureFluidPhase.
|
inlinevirtual |
Return whether phase represents a compressible substance.
Reimplemented in IdealSolidSolnPhase, LatticePhase, LatticeSolidPhase, MetalPhase, StoichSubstance, and VPStandardStateTP.
Definition at line 299 of file Phase.h.
Referenced by Phase::assertCompressible(), Phase::fullStates(), Phase::nativeState(), Phase::partialStates(), Phase::restoreState(), and Phase::saveState().
|
virtual |
Return a map of properties defining the native state of a substance.
By default, entries include "T", "D", "Y" for a compressible substance and "T", "P", "Y" for an incompressible substance, with offsets 0, 1 and 2, respectively. Mass fractions "Y" are omitted for pure species. In all cases, offsets into the state vector are used by saveState() and restoreState().
Reimplemented in LatticePhase, and LatticeSolidPhase.
Definition at line 230 of file Phase.cpp.
References Phase::isCompressible(), and Phase::isPure().
Referenced by ThermoPhase::getParameters(), Phase::restoreState(), and Phase::saveState().
|
virtual |
Return a vector containing full states defining a phase.
Full states list combinations of properties that allow for the specification of a thermodynamic state based on user input. Properties and states are represented by single letter acronyms, and combinations of letters, respectively (for example, "TDY", "TPX", "SVX"). Supported property acronyms are: "T": temperature "P": pressure "D": density "X": mole fractions "Y": mass fractions "T": temperature "U": specific internal energy "V": specific volume "H": specific enthalpy "S": specific entropy "Q": vapor fraction
Reimplemented in PureFluidPhase.
Definition at line 247 of file Phase.cpp.
References Phase::isCompressible(), and Phase::isPure().
|
virtual |
Return a vector of settable partial property sets within a phase.
Partial states encompass all valid combinations of properties that allow for the specification of a state while ignoring species concentrations (such as "TD", "TP", "SV").
Reimplemented in PureFluidPhase.
Definition at line 265 of file Phase.cpp.
References Phase::isCompressible(), and Phase::isPure().
|
virtual |
Return size of vector defining internal state of the phase.
Used by saveState() and restoreState().
Definition at line 278 of file Phase.cpp.
References Phase::isPure(), and Phase::nSpecies().
Referenced by Phase::restoreState(), and Phase::saveState().
void saveState | ( | vector_fp & | state | ) | const |
Save the current internal state of the phase.
Write to vector 'state' the current internal state.
state | output vector. Will be resized to stateSize(). |
Definition at line 286 of file Phase.cpp.
References Phase::saveState(), and Phase::stateSize().
Referenced by ThermoPhase::equilibrate(), ChemEquil::equilibrate(), ChemEquil::estimateEP_Brinkley(), TransportFactory::newTransport(), and Phase::saveState().
|
virtual |
Write to array 'state' the current internal state.
lenstate | length of the state array. Must be >= stateSize() |
state | output vector. Must be of length stateSizes() or greater. |
Definition at line 292 of file Phase.cpp.
References Phase::density(), Phase::getMassFractions(), Phase::getMoleFractions(), Phase::isCompressible(), Phase::nativeState(), Phase::pressure(), and Phase::temperature().
void restoreState | ( | const vector_fp & | state | ) |
Restore a state saved on a previous call to saveState.
state | State vector containing the previously saved state. |
Definition at line 310 of file Phase.cpp.
References Phase::restoreState().
Referenced by ChemEquil::equilibrate(), ConstPressureReactor::eval(), IdealGasConstPressureReactor::eval(), IdealGasReactor::eval(), MultiTransport::getMassFluxes(), ConstPressureReactor::getState(), IdealGasConstPressureReactor::getState(), IdealGasReactor::getState(), TransportFactory::newTransport(), and Phase::restoreState().
|
virtual |
Restore the state of the phase from a previously saved state vector.
lenstate | Length of the state vector |
state | Vector of state conditions. |
Definition at line 315 of file Phase.cpp.
References Phase::compositionChanged(), Phase::isCompressible(), Phase::nativeState(), Phase::setDensity(), Phase::setMassFractions_NoNorm(), Phase::setMoleFractions_NoNorm(), Phase::setPressure(), Phase::setTemperature(), and Phase::stateSize().
void setMoleFractionsByName | ( | const compositionMap & | xMap | ) |
Set the species mole fractions by name.
Species not listed by name in xMap
are set to zero.
xMap | map from species names to mole fraction values. |
Definition at line 380 of file Phase.cpp.
References Phase::getCompositionFromMap(), and Phase::setMoleFractions().
Referenced by OutletRes1D::setMoleFractions(), Phase::setMoleFractionsByName(), ThermoPhase::setState(), ThermoPhase::setState_RPX(), ThermoPhase::setState_TPX(), Phase::setState_TRX(), MixtureFugacityTP::setStateFromXML(), and ThermoPhase::setStateFromXML().
void setMoleFractionsByName | ( | const std::string & | x | ) |
Set the mole fractions of a group of species by name.
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 |
Definition at line 386 of file Phase.cpp.
References Cantera::parseCompString(), and Phase::setMoleFractionsByName().
void setMassFractionsByName | ( | const compositionMap & | yMap | ) |
Set the species mass fractions by name.
Species not listed by name in yMap
are set to zero.
yMap | map from species names to mass fraction values. |
Definition at line 416 of file Phase.cpp.
References Phase::getCompositionFromMap(), and Phase::setMassFractions().
Referenced by Phase::setMassFractionsByName(), ThermoPhase::setState(), ThermoPhase::setState_RPY(), ThermoPhase::setState_TPY(), Phase::setState_TRY(), MixtureFugacityTP::setStateFromXML(), and ThermoPhase::setStateFromXML().
void setMassFractionsByName | ( | const std::string & | x | ) |
Set the species mass fractions by name.
Species not listed by name in x
are set to zero.
x | String containing a composition map |
Definition at line 422 of file Phase.cpp.
References Cantera::parseCompString(), and Phase::setMassFractionsByName().
void setState_TRX | ( | doublereal | t, |
doublereal | dens, | ||
const doublereal * | x | ||
) |
Set the internally stored temperature (K), density, and mole fractions.
t | Temperature in kelvin |
dens | Density (kg/m^3) |
x | vector of species mole fractions, length m_kk |
Definition at line 427 of file Phase.cpp.
References Phase::setDensity(), Phase::setMoleFractions(), and Phase::setTemperature().
void setState_TRX | ( | doublereal | t, |
doublereal | dens, | ||
const compositionMap & | x | ||
) |
Set the internally stored temperature (K), density, and mole fractions.
t | Temperature in kelvin |
dens | Density (kg/m^3) |
x | Composition Map containing the mole fractions. Species not included in the map are assumed to have a zero mole fraction. |
Definition at line 441 of file Phase.cpp.
References Phase::setDensity(), Phase::setMoleFractionsByName(), and Phase::setTemperature().
void setState_TRY | ( | doublereal | t, |
doublereal | dens, | ||
const doublereal * | y | ||
) |
Set the internally stored temperature (K), density, and mass fractions.
t | Temperature in kelvin |
dens | Density (kg/m^3) |
y | vector of species mass fractions, length m_kk |
Definition at line 448 of file Phase.cpp.
References Phase::setDensity(), Phase::setMassFractions(), and Phase::setTemperature().
void setState_TRY | ( | doublereal | t, |
doublereal | dens, | ||
const compositionMap & | y | ||
) |
Set the internally stored temperature (K), density, and mass fractions.
t | Temperature in kelvin |
dens | Density (kg/m^3) |
y | Composition Map containing the mass fractions. Species not included in the map are assumed to have a zero mass fraction. |
Definition at line 455 of file Phase.cpp.
References Phase::setDensity(), Phase::setMassFractionsByName(), and Phase::setTemperature().
void setState_TNX | ( | doublereal | t, |
doublereal | n, | ||
const doublereal * | x | ||
) |
Set the internally stored temperature (K), molar density (kmol/m^3), and mole fractions.
t | Temperature in kelvin |
n | molar density (kmol/m^3) |
x | vector of species mole fractions, length m_kk |
Definition at line 434 of file Phase.cpp.
References Phase::setMolarDensity(), Phase::setMoleFractions(), and Phase::setTemperature().
void setState_TR | ( | doublereal | t, |
doublereal | rho | ||
) |
Set the internally stored temperature (K) and density (kg/m^3)
t | Temperature in kelvin |
rho | Density (kg/m^3) |
Definition at line 462 of file Phase.cpp.
References Phase::setDensity(), and Phase::setTemperature().
Referenced by ThermoPhase::setState_TP(), MixtureFugacityTP::setStateFromXML(), and IdealGasReactor::updateState().
void setState_TX | ( | doublereal | t, |
doublereal * | x | ||
) |
Set the internally stored temperature (K) and mole fractions.
t | Temperature in kelvin |
x | vector of species mole fractions, length m_kk |
Definition at line 468 of file Phase.cpp.
References Phase::setMoleFractions(), and Phase::setTemperature().
void setState_TY | ( | doublereal | t, |
doublereal * | y | ||
) |
Set the internally stored temperature (K) and mass fractions.
t | Temperature in kelvin |
y | vector of species mass fractions, length m_kk |
Definition at line 474 of file Phase.cpp.
References Phase::setMassFractions(), and Phase::setTemperature().
void setState_RX | ( | doublereal | rho, |
doublereal * | x | ||
) |
Set the density (kg/m^3) and mole fractions.
rho | Density (kg/m^3) |
x | vector of species mole fractions, length m_kk |
Definition at line 480 of file Phase.cpp.
References Phase::setDensity(), and Phase::setMoleFractions().
void setState_RY | ( | doublereal | rho, |
doublereal * | y | ||
) |
Set the density (kg/m^3) and mass fractions.
rho | Density (kg/m^3) |
y | vector of species mass fractions, length m_kk |
Definition at line 486 of file Phase.cpp.
References Phase::setDensity(), and Phase::setMassFractions().
doublereal molecularWeight | ( | size_t | k | ) | const |
Molecular weight of species k
.
k | index of species k |
k
. Definition at line 492 of file Phase.cpp.
References Phase::checkSpeciesIndex(), and Phase::m_molwts.
Referenced by IdealSolidSolnPhase::addSpecies(), LatticePhase::addSpecies(), MolalityVPSSTP::addSpecies(), InterfaceKinetics::buildSurfaceArrhenius(), Phase::elementalMassFraction(), SingleSpeciesTP::getPartialMolarVolumes(), SingleSpeciesTP::getStandardVolumes(), VPStandardStateTP::installPDSS(), and StickingCoverage::setContext().
void getMolecularWeights | ( | vector_fp & | weights | ) | const |
Copy the vector of molecular weights into vector weights.
weights | Output vector of molecular weights (kg/kmol) |
Definition at line 498 of file Phase.cpp.
References Phase::molecularWeights().
void getMolecularWeights | ( | doublereal * | weights | ) | const |
Copy the vector of molecular weights into array weights.
weights | Output array of molecular weights (kg/kmol) |
Definition at line 503 of file Phase.cpp.
References Phase::molecularWeights().
const vector_fp & molecularWeights | ( | ) | const |
Return a const reference to the internal vector of molecular weights.
units = kg / kmol
Definition at line 509 of file Phase.cpp.
References Phase::m_molwts.
Referenced by ConstPressureReactor::eval(), IdealGasConstPressureReactor::eval(), IdealGasReactor::eval(), ReactingSurf1D::eval(), Phase::getMolecularWeights(), MixTransport::getSpeciesFluxes(), DustyGasTransport::initialize(), ThermoPhase::o2Present(), ThermoPhase::o2Required(), and StFlow::StFlow().
void getCharges | ( | double * | charges | ) | const |
Copy the vector of species charges into array charges.
charges | Output array of species charges (elem. charge) |
Definition at line 514 of file Phase.cpp.
References Phase::m_speciesCharge.
compositionMap getMoleFractionsByName | ( | double | threshold = 0.0 | ) | const |
Get the mole fractions by name.
threshold | Exclude species with mole fractions less than or equal to this threshold. |
Definition at line 519 of file Phase.cpp.
References Phase::m_kk, Phase::moleFraction(), and Phase::speciesName().
double moleFraction | ( | size_t | k | ) | const |
Return the mole fraction of a single species.
k | species index |
Definition at line 548 of file Phase.cpp.
References Phase::checkSpeciesIndex(), Phase::m_mmw, and Phase::m_ym.
Referenced by BinarySolutionTabulatedThermo::_updateThermo(), BinarySolutionTabulatedThermo::calcDensity(), Phase::chargeDensity(), Phase::elementalMoleFraction(), SurfPhase::entropy_mole(), ChemEquil::equilibrate(), DebyeHuckel::getActivities(), HMWSoln::getActivities(), IdealMolalSoln::getActivities(), MolalityVPSSTP::getActivityCoefficients(), MixtureFugacityTP::getActivityConcentrations(), IdealSolnGasVPSS::getActivityConcentrations(), PengRobinson::getChemPotentials(), DebyeHuckel::getChemPotentials(), HMWSoln::getChemPotentials(), IdealGasPhase::getChemPotentials(), IdealMolalSoln::getChemPotentials(), IdealSolidSolnPhase::getChemPotentials(), IdealSolnGasVPSS::getChemPotentials(), LatticePhase::getChemPotentials(), RedlichKwongMFTP::getChemPotentials(), IdealSolidSolnPhase::getChemPotentials_RT(), IdealMolalSoln::getMolalityActivityCoefficients(), Phase::getMoleFractionsByName(), ThermoPhase::getParameters(), DebyeHuckel::getPartialMolarEntropies(), HMWSoln::getPartialMolarEntropies(), IdealGasPhase::getPartialMolarEntropies(), IdealMolalSoln::getPartialMolarEntropies(), IdealSolidSolnPhase::getPartialMolarEntropies(), IdealSolnGasVPSS::getPartialMolarEntropies(), LatticePhase::getPartialMolarEntropies(), RedlichKwongMFTP::getPartialMolarEntropies(), BinarySolutionTabulatedThermo::getPartialMolarVolumes(), Phase::moleFraction(), vcs_MultiPhaseEquil::reportCSV(), DebyeHuckel::s_update_d2lnMolalityActCoeff_dT2(), DebyeHuckel::s_update_dlnMolalityActCoeff_dP(), DebyeHuckel::s_update_dlnMolalityActCoeff_dT(), DebyeHuckel::s_update_lnMolalityActCoeff(), HMWSoln::s_updateIMS_lnMolalityActCoeff(), IdealMolalSoln::s_updateIMS_lnMolalityActCoeff(), and ChemEquil::setInitialMoles().
double moleFraction | ( | const std::string & | name | ) | const |
Return the mole fraction of a single species.
name | String name of the species |
Definition at line 554 of file Phase.cpp.
References Phase::moleFraction(), Cantera::npos, and Phase::speciesIndex().
compositionMap getMassFractionsByName | ( | double | threshold = 0.0 | ) | const |
Get the mass fractions by name.
threshold | Exclude species with mass fractions less than or equal to this threshold. |
Definition at line 531 of file Phase.cpp.
References Phase::m_kk, Phase::massFraction(), and Phase::speciesName().
double massFraction | ( | size_t | k | ) | const |
Return the mass fraction of a single species.
k | species index |
Definition at line 569 of file Phase.cpp.
References Phase::checkSpeciesIndex(), and Phase::m_y.
Referenced by Phase::elementalMassFraction(), Phase::getMassFractionsByName(), and ThermoPhase::getParameters().
double massFraction | ( | const std::string & | name | ) | const |
Return the mass fraction of a single species.
name | String name of the species |
Definition at line 575 of file Phase.cpp.
References Phase::massFractions(), Cantera::npos, and Phase::speciesIndex().
void getMoleFractions | ( | double *const | x | ) | const |
Get the species mole fraction vector.
x | On return, x contains the mole fractions. Must have a length greater than or equal to the number of species. |
Definition at line 543 of file Phase.cpp.
References Phase::m_mmw, Phase::m_ym, and Cantera::scale().
Referenced by MolalityVPSSTP::calcMolalities(), HMWSoln::calcMolalitiesCropped(), GibbsExcessVPSSTP::compositionChanged(), MixtureFugacityTP::compositionChanged(), ChemEquil::estimateElementPotentials(), GibbsExcessVPSSTP::getActivities(), LatticePhase::getActivityConcentrations(), HighPressureGasTransport::getBinaryDiffCoeffs(), MolalityVPSSTP::getCsvReportData(), ThermoPhase::getCsvReportData(), MultiTransport::getMassFluxes(), LatticeSolidPhase::getMoleFractions(), HighPressureGasTransport::getMultiDiffCoeffs(), DustyGasTransport::initialize(), HMWSoln::printCoeffs(), HMWSoln::relative_molal_enthalpy(), ThermoPhase::reportCSV(), Phase::saveState(), MolalityVPSSTP::setMolalitiesByName(), MultiPhase::setMoles(), HighPressureGasTransport::thermalConductivity(), ChemEquil::update(), MixTransport::update_C(), MultiTransport::update_C(), solveSP::updateMFKinSpecies(), DustyGasTransport::updateTransport_C(), MultiPhase::uploadMoleFractionsFromPhases(), and HighPressureGasTransport::viscosity().
|
virtual |
Set the mole fractions to the specified values.
There is no restriction on the sum of the mole fraction vector. Internally, the Phase object will normalize this vector before storing its contents.
x | Array of unnormalized mole fraction values (input). Must have a length greater than or equal to the number of species, m_kk. |
Reimplemented in LatticeSolidPhase, and SingleSpeciesTP.
Definition at line 339 of file Phase.cpp.
References Phase::compositionChanged(), Phase::m_kk, Phase::m_mmw, Phase::m_molwts, Phase::m_y, and Phase::m_ym.
Referenced by ChemEquil::calcEmoles(), ChemEquil::equilibrate(), ChemEquil::estimateElementPotentials(), MolalityVPSSTP::setMolalities(), MolalityVPSSTP::setMolalitiesByName(), OutletRes1D::setMoleFractions(), LatticeSolidPhase::setMoleFractions(), Phase::setMoleFractionsByName(), ThermoPhase::setState_PX(), ThermoPhase::setState_RPX(), Phase::setState_RX(), Phase::setState_TNX(), ThermoPhase::setState_TPX(), Phase::setState_TRX(), and Phase::setState_TX().
|
virtual |
Set the mole fractions to the specified values without normalizing.
This is useful when the normalization condition is being handled by some other means, for example by a constraint equation as part of a larger set of equations.
x | Input vector of mole fractions. Length is m_kk. |
Definition at line 371 of file Phase.cpp.
References Phase::compositionChanged(), Cantera::dot(), Phase::m_kk, Phase::m_mmw, Phase::m_molwts, Phase::m_y, Phase::m_ym, and Cantera::scale().
Referenced by Phase::restoreState().
void getMassFractions | ( | double *const | y | ) | const |
Get the species mass fractions.
[out] | y | Array of mass fractions, length nSpecies() |
Definition at line 585 of file Phase.cpp.
References Phase::m_y.
Referenced by StFlow::_getInitialSoln(), ThermoPhase::getCsvReportData(), ConstPressureReactor::getState(), IdealGasConstPressureReactor::getState(), IdealGasReactor::getState(), StFlow::resetBadValues(), Phase::saveState(), and OutletRes1D::setMoleFractions().
|
inline |
Return a const pointer to the mass fraction array.
Definition at line 532 of file Phase.h.
References Phase::m_y.
Referenced by ThermoPhase::equivalenceRatio(), ConstPressureReactor::eval(), IdealGasConstPressureReactor::eval(), IdealGasReactor::eval(), MultiTransport::getMassFluxes(), MixTransport::getSpeciesFluxes(), MultiTransport::getSpeciesFluxes(), and Phase::massFraction().
|
virtual |
Set the mass fractions to the specified values and normalize them.
[in] | y | Array of unnormalized mass fraction values. Length must be greater than or equal to the number of species. The Phase object will normalize this vector before storing its contents. |
Reimplemented in LatticeSolidPhase, and SingleSpeciesTP.
Definition at line 391 of file Phase.cpp.
References Phase::compositionChanged(), Phase::m_kk, Phase::m_mmw, Phase::m_rmolwts, Phase::m_y, Phase::m_ym, and Cantera::scale().
Referenced by StFlow::resetBadValues(), Phase::setMassFractionsByName(), ThermoPhase::setState_PY(), ThermoPhase::setState_RPY(), Phase::setState_RY(), ThermoPhase::setState_TPY(), Phase::setState_TRY(), and Phase::setState_TY().
|
virtual |
Set the mass fractions to the specified values without normalizing.
This is useful when the normalization condition is being handled by some other means, for example by a constraint equation as part of a larger set of equations.
y | Input vector of mass fractions. Length is m_kk. |
Reimplemented in LatticeSolidPhase.
Definition at line 405 of file Phase.cpp.
References Phase::compositionChanged(), Phase::m_kk, Phase::m_mmw, Phase::m_rmolwts, Phase::m_y, and Phase::m_ym.
Referenced by Phase::restoreState(), StFlow::setGas(), StFlow::setGasAtMidpoint(), ConstPressureReactor::updateState(), IdealGasConstPressureReactor::updateState(), and IdealGasReactor::updateState().
void getConcentrations | ( | double *const | c | ) | const |
Get the species concentrations (kmol/m^3).
[out] | c | The vector of species concentrations. Units are kmol/m^3. The length of the vector must be greater than or equal to the number of species within the phase. |
Definition at line 596 of file Phase.cpp.
References Phase::m_dens, Phase::m_ym, and Cantera::scale().
Referenced by InterfaceKinetics::_update_rates_C(), IdealGasPhase::getActivityConcentrations(), SurfPhase::getActivityConcentrations(), ImplicitSurfChem::getConcSpecies(), SurfPhase::getCoverages(), and GasKinetics::update_rates_C().
double concentration | ( | const size_t | k | ) | const |
Concentration of species k.
If k is outside the valid range, an exception will be thrown.
[in] | k | Index of the species within the phase. |
Definition at line 590 of file Phase.cpp.
References Phase::checkSpeciesIndex(), Phase::m_dens, Phase::m_rmolwts, and Phase::m_y.
Referenced by SurfPhase::entropy_mole().
|
virtual |
Set the concentrations to the specified values within the phase.
We set the concentrations here and therefore we set the overall density of the phase. We hold the temperature constant during this operation. Therefore, we have possibly changed the pressure of the phase by calling this routine.
[in] | conc | Array of concentrations in dimensional units. For bulk phases c[k] is the concentration of the kth species in kmol/m3. For surface phases, c[k] is the concentration in kmol/m2. The length of the vector is the number of species in the phase. |
Reimplemented in LatticeSolidPhase.
Definition at line 601 of file Phase.cpp.
References Phase::assertCompressible(), Phase::compositionChanged(), Phase::m_kk, Phase::m_mmw, Phase::m_molwts, Phase::m_y, Phase::m_ym, and Phase::setDensity().
Referenced by ImplicitSurfChem::setConcSpecies(), and SurfPhase::setCoverages().
|
virtual |
Set the concentrations without ignoring negative concentrations.
Definition at line 623 of file Phase.cpp.
References Phase::assertCompressible(), Phase::compositionChanged(), Phase::m_kk, Phase::m_mmw, Phase::m_molwts, Phase::m_y, Phase::m_ym, and Phase::setDensity().
Referenced by SurfPhase::setCoveragesNoNorm().
doublereal elementalMassFraction | ( | const size_t | m | ) | const |
Elemental mass fraction of element m.
The elemental mass fraction \(Z_{\mathrm{mass},m}\) of element \(m\) is defined as
\[ Z_{\mathrm{mass},m} = \sum_k \frac{a_{m,k} M_m}{M_k} Y_k \]
with \(a_{m,k}\) being the number of atoms of element \(m\) in species \(k\), \(M_m\) the atomic weight of element \(m\), \(M_k\) the molecular weight of species \(k\), and \(Y_k\) the mass fraction of species \(k\).
[in] | m | Index of the element within the phase. If m is outside the valid range, an exception will be thrown. |
Definition at line 642 of file Phase.cpp.
References Phase::atomicWeight(), Phase::checkElementIndex(), Phase::m_kk, Phase::massFraction(), Phase::molecularWeight(), and Phase::nAtoms().
doublereal elementalMoleFraction | ( | const size_t | m | ) | const |
Elemental mole fraction of element m.
The elemental mole fraction \(Z_{\mathrm{mole},m}\) of element \(m\) is the number of atoms of element m divided by the total number of atoms. It is defined as:
\[ Z_{\mathrm{mole},m} = \frac{\sum_k a_{m,k} X_k} {\sum_k \sum_j a_{j,k} X_k} \]
with \(a_{m,k}\) being the number of atoms of element \(m\) in species \(k\), \(\sum_j\) being a sum over all elements, and \(X_k\) being the mole fraction of species \(k\).
[in] | m | Index of the element within the phase. If m is outside the valid range, an exception will be thrown. |
Definition at line 653 of file Phase.cpp.
References Phase::checkElementIndex(), Phase::m_kk, Phase::moleFraction(), Phase::nAtoms(), and Phase::nElements().
const double * moleFractdivMMW | ( | ) | const |
Returns a const pointer to the start of the moleFraction/MW array.
This array is the array of mole fractions, each divided by the mean molecular weight.
Definition at line 564 of file Phase.cpp.
References Phase::m_ym.
Referenced by IdealSolidSolnPhase::calcDensity(), IdealSolnGasVPSS::calcDensity(), and IdealSolidSolnPhase::getActivityConcentrations().
|
inline |
Dimensionless electrical charge of a single molecule of species k The charge is normalized by the the magnitude of the electron charge.
k | species index |
Definition at line 630 of file Phase.h.
References Phase::m_speciesCharge.
Referenced by RedlichKisterVPSSTP::addBinaryInteraction(), HMWSoln::applyphScale(), InterfaceKinetics::applyVoltageKfwdCorrection(), HMWSoln::calcMolalitiesCropped(), Phase::chargeDensity(), Cantera::checkElectrochemReaction(), ThermoPhase::getElectrochemPotentials(), DebyeHuckel::getSpeciesParameters(), Cantera::hasChargedSpecies(), PDSS_HKFT::initThermo(), HMWSoln::printCoeffs(), HMWSoln::relative_molal_enthalpy(), HMWSoln::s_updatePitzer_d2lnMolalityActCoeff_dT2(), HMWSoln::s_updatePitzer_dlnMolalityActCoeff_dP(), HMWSoln::s_updatePitzer_dlnMolalityActCoeff_dT(), HMWSoln::s_updatePitzer_lnMolalityActCoeff(), HMWSoln::s_updateScaling_pHScaling(), HMWSoln::s_updateScaling_pHScaling_dP(), HMWSoln::s_updateScaling_pHScaling_dT(), HMWSoln::s_updateScaling_pHScaling_dT2(), InterfaceRateBase::setContext(), MolalityVPSSTP::setMolalitiesByName(), vcs_VolPhase::transferElementsFM(), InterfaceKinetics::updateMu0(), and Reaction::usesElectrochemistry().
doublereal chargeDensity | ( | ) | const |
Charge density [C/m^3].
Definition at line 708 of file Phase.cpp.
References Phase::charge(), Cantera::Faraday, Phase::m_kk, and Phase::moleFraction().
|
inline |
Returns the number of spatial dimensions (1, 2, or 3)
Definition at line 638 of file Phase.h.
References Phase::m_ndim.
Referenced by Kinetics::addPhase(), InterfaceKinetics::buildSurfaceArrhenius(), Reaction::checkBalance(), InterfaceKinetics::init(), Reaction::Reaction(), IdealMolalSoln::standardConcentrationUnits(), IdealSolidSolnPhase::standardConcentrationUnits(), and ThermoPhase::standardConcentrationUnits().
|
inline |
Set the number of spatial dimensions (1, 2, or 3).
The number of spatial dimensions is used for vector involving directions.
ndim | Input number of dimensions. |
Definition at line 645 of file Phase.h.
References Phase::m_ndim.
Referenced by EdgePhase::EdgePhase(), Cantera::importPhase(), and SurfPhase::SurfPhase().
|
inline |
Temperature (K).
Definition at line 654 of file Phase.h.
References Phase::m_temp.
Referenced by StFlow::_getInitialSoln(), MixtureFugacityTP::_updateReferenceStateThermo(), VPStandardStateTP::_updateStandardStateThermo(), BinarySolutionTabulatedThermo::_updateThermo(), IdealSolidSolnPhase::_updateThermo(), LatticePhase::_updateThermo(), LatticeSolidPhase::_updateThermo(), SingleSpeciesTP::_updateThermo(), SurfPhase::_updateThermo(), MultiPhase::addPhase(), HMWSoln::calcDensity(), IonsFromNeutralVPSSTP::calcDensity(), MixtureFugacityTP::calculatePsat(), RedlichKwongMFTP::cp_mole(), RedlichKwongMFTP::cv_mole(), PengRobinson::d2aAlpha_dT2(), ChemEquil::equilibrate(), ChemEquil::estimateElementPotentials(), RedlichKwongMFTP::getActivityCoefficients(), RedlichKwongMFTP::getChemPotentials(), IdealSolidSolnPhase::getChemPotentials_RT(), PureFluidPhase::getEnthalpy_RT_ref(), PureFluidPhase::getEntropy_R_ref(), PureFluidPhase::getGibbs_RT_ref(), ThermoPhase::getParameters(), DebyeHuckel::getPartialMolarCp(), HMWSoln::getPartialMolarCp(), DebyeHuckel::getPartialMolarEnthalpies(), HMWSoln::getPartialMolarEnthalpies(), IonsFromNeutralVPSSTP::getPartialMolarEnthalpies(), MargulesVPSSTP::getPartialMolarEnthalpies(), RedlichKwongMFTP::getPartialMolarEnthalpies(), IonsFromNeutralVPSSTP::getPartialMolarEntropies(), RedlichKwongMFTP::getPartialMolarEntropies(), RedlichKwongMFTP::getPartialMolarVolumes(), WaterSSTP::getStandardChemPotentials(), IdealGasConstPressureReactor::getState(), IdealGasReactor::getState(), ThermoPhase::gibbs_mole(), PureFluidPhase::phaseOfMatter(), MixtureFugacityTP::phaseState(), PlasmaPhase::PlasmaPhase(), IdealGasPhase::pressure(), RedlichKwongMFTP::pressureDerivatives(), ThermoPhase::RT(), HMWSoln::s_update_d2lnMolalityActCoeff_dT2(), MargulesVPSSTP::s_update_dlnActCoeff_dT(), HMWSoln::s_update_dlnMolalityActCoeff_dP(), HMWSoln::s_update_dlnMolalityActCoeff_dT(), HMWSoln::s_updatePitzer_dlnMolalityActCoeff_dP(), HMWSoln::satPressure(), WaterSSTP::satPressure(), Phase::saveState(), WaterSSTP::setDensity(), ChemEquil::setInitialMoles(), MixtureFugacityTP::setPressure(), VPStandardStateTP::setPressure(), vcs_VolPhase::setPtrThermoPhase(), PengRobinson::setSpeciesCoeffs(), SingleSpeciesTP::setState_HP(), ThermoPhase::setState_HPorUV(), SingleSpeciesTP::setState_SP(), ThermoPhase::setState_SPorSV(), SingleSpeciesTP::setState_SV(), ThermoPhase::setState_TP(), SingleSpeciesTP::setState_UV(), MixtureFugacityTP::setStateFromXML(), ImplicitSurfChem::solvePseudoSteadyStateProblem(), HighPressureGasTransport::thermalConductivity(), WaterTransport::thermalConductivity(), IdealGasPhase::thermalExpansionCoeff(), ChemEquil::update(), MixTransport::update_T(), MultiTransport::update_T(), PengRobinson::updateMixingExpressions(), RedlichKwongMFTP::updateMixingExpressions(), VPStandardStateTP::updateStandardStateThermo(), MultiTransport::updateThermal_T(), IdealGasPhase::updateThermo(), PlasmaPhase::updateThermo(), DustyGasTransport::updateTransport_T(), WaterSSTP::vaporFraction(), HighPressureGasTransport::viscosity(), and WaterTransport::viscosity().
|
inlinevirtual |
Electron Temperature (K)
Reimplemented in PlasmaPhase.
Definition at line 660 of file Phase.h.
References Phase::m_temp.
|
inlinevirtual |
Return the thermodynamic pressure (Pa).
This method must be overloaded in derived classes. Within Cantera, the independent variable is either density or pressure. If the state is defined by temperature, density, and mass fractions, this method should use these values to implement the mechanical equation of state \( P(T, \rho, Y_1, \dots, Y_K) \). Alternatively, it returns a stored value.
Reimplemented in IdealGasPhase, IdealSolidSolnPhase, LatticePhase, LatticeSolidPhase, MaskellSolidSolnPhase, MetalPhase, PengRobinson, PureFluidPhase, RedlichKwongMFTP, StoichSubstance, SurfPhase, VPStandardStateTP, and WaterSSTP.
Definition at line 672 of file Phase.h.
References Phase::type().
Referenced by MultiPhase::addPhase(), MixtureFugacityTP::entropy_mole(), ChemEquil::equilibrate(), ChemEquil::estimateElementPotentials(), MixtureFugacityTP::getActivityConcentrations(), HighPressureGasTransport::getBinaryDiffCoeffs(), MixtureFugacityTP::getEntropy_R(), MixtureFugacityTP::getGibbs_RT(), MultiTransport::getMassFluxes(), HighPressureGasTransport::getMultiDiffCoeffs(), ThermoPhase::getParameters(), MixtureFugacityTP::getPureGibbs(), MixtureFugacityTP::getStandardChemPotentials(), MixtureFugacityTP::getStandardVolumes(), ThermoPhase::intEnergy_mole(), Phase::saveState(), ThermoPhase::setEquivalenceRatio(), ChemEquil::setInitialMoles(), ImplicitSurfChem::solvePseudoSteadyStateProblem(), and MixtureFugacityTP::z().
|
inlinevirtual |
Density (kg/m^3).
Definition at line 679 of file Phase.h.
References Phase::m_dens.
Referenced by MixtureFugacityTP::calculatePsat(), WaterSSTP::dthermalExpansionCoeffdT(), ChemEquil::equilibrate(), MultiTransport::getMassFluxes(), UnityLewisTransport::getMixDiffCoeffs(), UnityLewisTransport::getMixDiffCoeffsMass(), MetalPhase::getParameters(), ThermoPhase::getParameters(), StoichSubstance::getParameters(), SingleSpeciesTP::getPartialMolarVolumes(), MultiTransport::getSpeciesFluxes(), SingleSpeciesTP::getStandardVolumes(), ConstPressureReactor::getState(), IdealGasConstPressureReactor::getState(), IdealGasReactor::getState(), RedlichKwongMFTP::hresid(), Phase::molarDensity(), WaterSSTP::satPressure(), Phase::saveState(), MixtureFugacityTP::setPressure(), ThermoPhase::setState_TP(), MixtureFugacityTP::setStateFromXML(), WaterSSTP::setTemperature(), RedlichKwongMFTP::sresid(), WaterTransport::thermalConductivity(), BlowersMaselData::update(), ChemEquil::update(), ConstPressureReactor::updateState(), IdealGasConstPressureReactor::updateState(), StFlow::updateThermo(), StFlow::updateTransport(), WaterSSTP::vaporFraction(), WaterTransport::viscosity(), and MixtureFugacityTP::z().
double molarDensity | ( | ) | const |
Molar density (kmol/m^3).
Definition at line 671 of file Phase.cpp.
References Phase::density(), and Phase::meanMolecularWeight().
Referenced by solveSP::calc_t(), IdealSolidSolnPhase::enthalpy_mole(), LatticePhase::enthalpy_mole(), StoichSubstance::getEnthalpy_RT(), StoichSubstance::getIntEnergy_RT(), StoichSubstance::getIntEnergy_RT_ref(), IdealGasPhase::getPartialMolarVolumes(), PureFluidPhase::getPartialMolarVolumes(), MixTransport::getSpeciesFluxes(), IdealGasPhase::getStandardVolumes(), Phase::molarVolume(), IdealGasPhase::pressure(), GasKinetics::process_ddC(), GasKinetics::process_ddX(), FalloffData::update(), and GasKinetics::update_rates_C().
double molarVolume | ( | ) | const |
Molar volume (m^3/kmol).
Definition at line 682 of file Phase.cpp.
References Phase::molarDensity().
Referenced by RedlichKwongMFTP::cp_mole(), RedlichKwongMFTP::cv_mole(), PengRobinson::getActivityCoefficients(), RedlichKwongMFTP::getActivityCoefficients(), PengRobinson::getChemPotentials(), RedlichKwongMFTP::getChemPotentials(), RedlichKwongMFTP::getPartialMolarEnthalpies(), RedlichKwongMFTP::getPartialMolarEntropies(), PengRobinson::getPartialMolarVolumes(), RedlichKwongMFTP::getPartialMolarVolumes(), PengRobinson::hresid(), ThermoPhase::intEnergy_mole(), RedlichKwongMFTP::pressureDerivatives(), and PengRobinson::sresid().
|
virtual |
Set the internally stored density (kg/m^3) of the phase.
Note the density of a phase is an independent variable.
[in] | density_ | density (kg/m^3). |
Reimplemented in PureFluidPhase, and WaterSSTP.
Definition at line 687 of file Phase.cpp.
References Phase::assertCompressible(), and Phase::m_dens.
Referenced by Phase::restoreState(), Phase::setConcentrations(), Phase::setConcentrationsNoNorm(), PureFluidPhase::setDensity(), WaterSSTP::setDensity(), IdealGasPhase::setPressure(), ThermoPhase::setState_HPorUV(), IdealGasPhase::setState_RP(), Phase::setState_RX(), Phase::setState_RY(), ThermoPhase::setState_SPorSV(), SingleSpeciesTP::setState_SV(), Phase::setState_TR(), Phase::setState_TRX(), Phase::setState_TRY(), SingleSpeciesTP::setState_UV(), and ThermoPhase::setStateFromXML().
|
virtual |
Set the internally stored molar density (kmol/m^3) of the phase.
[in] | molarDensity | Input molar density (kmol/m^3). |
Definition at line 676 of file Phase.cpp.
References Phase::assertCompressible(), Phase::m_dens, and Phase::meanMolecularWeight().
Referenced by Phase::setState_TNX().
|
inlinevirtual |
Set the internally stored pressure (Pa) at constant temperature and composition.
This method must be reimplemented in derived classes, where it may involve the solution of a nonlinear equation. Within Cantera, the independent variable is either density or pressure. Therefore, this function may either solve for the density that will yield the desired input pressure or set an independent variable. The temperature and composition are held constant during this process.
p | input Pressure (Pa) |
Reimplemented in IdealGasPhase, IdealSolidSolnPhase, IdealSolnGasVPSS, LatticePhase, LatticeSolidPhase, MaskellSolidSolnPhase, MixtureFugacityTP, PureFluidPhase, StoichSubstance, SurfPhase, VPStandardStateTP, WaterSSTP, and MetalPhase.
Definition at line 712 of file Phase.h.
References Phase::type().
Referenced by ChemEquil::calcEmoles(), Phase::restoreState(), ThermoPhase::setState_conditional_TP(), SingleSpeciesTP::setState_HP(), ThermoPhase::setState_HPorUV(), ThermoPhase::setState_PX(), ThermoPhase::setState_PY(), SingleSpeciesTP::setState_SP(), ThermoPhase::setState_SPorSV(), ThermoPhase::setState_TP(), ThermoPhase::setStateFromXML(), and ConstPressureReactor::updateState().
|
inlinevirtual |
Set the internally stored temperature of the phase (K).
temp | Temperature in Kelvin |
Reimplemented in PureFluidPhase, MixtureFugacityTP, VPStandardStateTP, and WaterSSTP.
Definition at line 719 of file Phase.h.
References Phase::m_temp.
Referenced by ChemEquil::equilibrate(), ReactingSurf1D::eval(), Phase::restoreState(), StFlow::setGas(), StFlow::setGasAtMidpoint(), ThermoPhase::setState_conditional_TP(), IdealGasPhase::setState_RP(), SingleSpeciesTP::setState_SV(), Phase::setState_TNX(), ThermoPhase::setState_TP(), VPStandardStateTP::setState_TP(), Phase::setState_TR(), Phase::setState_TRX(), Phase::setState_TRY(), Phase::setState_TX(), Phase::setState_TY(), SingleSpeciesTP::setState_UV(), SurfPhase::setStateFromXML(), ThermoPhase::setStateFromXML(), MixtureFugacityTP::setTemperature(), WaterSSTP::setTemperature(), and ChemEquil::setToEquilState().
|
inlinevirtual |
Set the internally stored electron temperature of the phase (K).
etemp | Electron temperature in Kelvin |
Reimplemented in PlasmaPhase.
Definition at line 730 of file Phase.h.
References Phase::type().
doublereal mean_X | ( | const doublereal *const | Q | ) | const |
Evaluate the mole-fraction-weighted mean of an array Q.
\[ \sum_k X_k Q_k. \]
Q should contain pure-species molar property values.
[in] | Q | Array of length m_kk that is to be averaged. |
Definition at line 717 of file Phase.cpp.
References Phase::m_mmw, and Phase::m_ym.
Referenced by DebyeHuckel::calcDensity(), HMWSoln::calcDensity(), IdealMolalSoln::calcDensity(), DebyeHuckel::cp_mole(), HMWSoln::cp_mole(), IdealGasPhase::cp_mole(), IdealMolalSoln::cp_mole(), IdealSolidSolnPhase::cp_mole(), IdealSolnGasVPSS::cp_mole(), IonsFromNeutralVPSSTP::cp_mole(), LatticePhase::cp_mole(), RedlichKwongMFTP::cp_mole(), SurfPhase::cp_mole(), IonsFromNeutralVPSSTP::cv_mole(), RedlichKwongMFTP::cv_mole(), DebyeHuckel::enthalpy_mole(), HMWSoln::enthalpy_mole(), IdealGasPhase::enthalpy_mole(), IdealMolalSoln::enthalpy_mole(), IdealSolidSolnPhase::enthalpy_mole(), IdealSolnGasVPSS::enthalpy_mole(), IonsFromNeutralVPSSTP::enthalpy_mole(), LatticePhase::enthalpy_mole(), MixtureFugacityTP::enthalpy_mole(), SurfPhase::enthalpy_mole(), DebyeHuckel::entropy_mole(), HMWSoln::entropy_mole(), IdealGasPhase::entropy_mole(), IdealMolalSoln::entropy_mole(), IdealSolidSolnPhase::entropy_mole(), IdealSolnGasVPSS::entropy_mole(), IonsFromNeutralVPSSTP::entropy_mole(), LatticePhase::entropy_mole(), MixtureFugacityTP::entropy_mole(), DebyeHuckel::gibbs_mole(), HMWSoln::gibbs_mole(), IdealMolalSoln::gibbs_mole(), IdealSolidSolnPhase::gibbs_mole(), IonsFromNeutralVPSSTP::gibbs_mole(), IdealMolalSoln::intEnergy_mole(), and HMWSoln::relative_enthalpy().
doublereal mean_X | ( | const vector_fp & | Q | ) | const |
Evaluate the mole-fraction-weighted mean of an array Q.
\[ \sum_k X_k Q_k. \]
Q should contain pure-species molar property values.
[in] | Q | Array of length m_kk that is to be averaged. |
Definition at line 722 of file Phase.cpp.
References Phase::m_mmw, and Phase::m_ym.
|
inline |
The mean molecular weight. Units: (kg/kmol)
Definition at line 751 of file Phase.h.
References Phase::m_mmw.
Referenced by BinarySolutionTabulatedThermo::calcDensity(), DebyeHuckel::calcDensity(), GibbsExcessVPSSTP::calcDensity(), HMWSoln::calcDensity(), IdealMolalSoln::calcDensity(), LatticePhase::calcDensity(), MixtureFugacityTP::calculatePsat(), ThermoPhase::cp_mass(), MixtureFugacityTP::critDensity(), ThermoPhase::cv_mass(), RedlichKwongMFTP::densityCalc(), MixtureFugacityTP::densityCalc(), ThermoPhase::enthalpy_mass(), ThermoPhase::entropy_mass(), IdealSolidSolnPhase::getActivityConcentrations(), HighPressureGasTransport::getMultiDiffCoeffs(), MultiTransport::getMultiDiffCoeffs(), ThermoPhase::gibbs_mass(), RedlichKwongMFTP::hresid(), StoichSubstance::initThermo(), ThermoPhase::intEnergy_mass(), Phase::molarDensity(), Phase::setMolarDensity(), IdealGasPhase::setPressure(), IdealGasPhase::setState_RP(), RedlichKwongMFTP::sresid(), StFlow::updateThermo(), StFlow::updateTransport(), HighPressureGasTransport::viscosity(), and MixtureFugacityTP::z().
doublereal sum_xlogx | ( | ) | const |
Evaluate \( \sum_k X_k \log X_k \).
Definition at line 727 of file Phase.cpp.
References Phase::m_kk, Phase::m_mmw, Phase::m_ym, and Cantera::SmallNumber.
Referenced by IdealGasPhase::entropy_mole(), IdealSolidSolnPhase::entropy_mole(), IdealSolnGasVPSS::entropy_mole(), LatticePhase::entropy_mole(), MixtureFugacityTP::entropy_mole(), and IdealSolidSolnPhase::gibbs_mole().
size_t addElement | ( | const std::string & | symbol, |
doublereal | weight = -12345.0 , |
||
int | atomicNumber = 0 , |
||
doublereal | entropy298 = ENTROPY298_UNKNOWN , |
||
int | elem_type = CT_ELEM_TYPE_ABSPOS |
||
) |
Add an element.
symbol | Atomic symbol std::string. |
weight | Atomic mass in amu. |
atomicNumber | Atomic number of the element (unitless) |
entropy298 | Entropy of the element at 298 K and 1 bar in its most stable form. The default is the value ENTROPY298_UNKNOWN, which is interpreted as an unknown, and if used will cause Cantera to throw an error. |
elem_type | Specifies the type of the element constraint equation. This defaults to CT_ELEM_TYPE_ABSPOS, that is, an element. |
Definition at line 736 of file Phase.cpp.
References AnyMap::convert(), CT_ELEM_TYPE_ELECTRONCHARGE, ENTROPY298_UNKNOWN, AnyMap::fromYamlFile(), Cantera::getElementWeight(), Phase::m_atomicNumbers, Phase::m_atomicWeights, Phase::m_elem_type, Phase::m_elementNames, Phase::m_entropy298, Phase::m_kk, Phase::m_mm, and Phase::m_speciesComp.
Referenced by LatticeSolidPhase::addLattice(), Phase::addSpecies(), Cantera::installElements(), and LatticeSolidPhase::setLatticeStoichiometry().
|
virtual |
Returns true
if the species was successfully added, or false
if the species was ignored.
Derived classes which need to size arrays according to the number of species should overload this method. The derived class implementation should call the base class method, and, if this returns true
(indicating that the species has been added), adjust their array sizes accordingly.
Reimplemented in BinarySolutionTabulatedThermo, DebyeHuckel, GibbsExcessVPSSTP, IdealGasPhase, IdealMolalSoln, IdealSolidSolnPhase, IdealSolnGasVPSS, IonsFromNeutralVPSSTP, LatticePhase, LatticeSolidPhase, MixtureFugacityTP, MolalityVPSSTP, PengRobinson, PlasmaPhase, RedlichKwongMFTP, SingleSpeciesTP, SurfPhase, ThermoPhase, VPStandardStateTP, and VPStandardStateTP.
Definition at line 803 of file Phase.cpp.
References Phase::addElement(), Phase::atomicWeights(), CT_ELEM_TYPE_ELECTRONCHARGE, Phase::elementIndex(), Phase::invalidateCache(), Phase::m_kk, Phase::m_mmw, Phase::m_molwts, Phase::m_name, Phase::m_rmolwts, Phase::m_speciesCharge, Phase::m_speciesComp, Phase::m_speciesIndices, Phase::m_speciesLower, Phase::m_speciesNames, Phase::m_undefinedElementBehavior, Phase::m_y, Phase::m_ym, Phase::nElements(), Cantera::npos, Cantera::Tiny, and Cantera::toLowerCopy().
Referenced by ThermoPhase::addSpecies(), and VPStandardStateTP::addSpecies().
|
virtual |
Modify the thermodynamic data associated with a species.
The species name, elemental composition, and type of thermo parameterization must be unchanged. If there are Kinetics objects that depend on this phase, Kinetics::invalidateCache() should be called on those objects after calling this function.
Reimplemented in ThermoPhase.
Definition at line 899 of file Phase.cpp.
References Phase::invalidateCache(), and Phase::speciesName().
Referenced by ThermoPhase::modifySpecies().
void addSpeciesAlias | ( | const std::string & | name, |
const std::string & | alias | ||
) |
Add a species alias (that is, a user-defined alternative species name).
Aliases are case-sensitive.
name | original species name std::string. |
alias | alternate name std::string. |
true
if the alias was successfully added (that is, the original species name is found) Definition at line 916 of file Phase.cpp.
References Phase::m_speciesIndices, Phase::name(), Cantera::npos, and Phase::speciesIndex().
|
virtual |
Return a vector with isomers names matching a given composition map.
compMap | compositionMap of the species. |
Definition at line 932 of file Phase.cpp.
Referenced by Phase::findIsomers().
|
virtual |
Return a vector with isomers names matching a given composition string.
comp | String containing a composition map |
Definition at line 945 of file Phase.cpp.
References Phase::findIsomers(), and Cantera::parseCompString().
shared_ptr< Species > species | ( | const std::string & | name | ) | const |
Return the Species object for the named species.
Changes to this object do not affect the ThermoPhase object until the modifySpecies function is called.
Definition at line 950 of file Phase.cpp.
References Phase::name(), Cantera::npos, Phase::speciesIndex(), and Phase::speciesName().
Referenced by RedlichKwongMFTP::getCoeff(), DebyeHuckel::initThermo(), MargulesVPSSTP::initThermo(), RedlichKisterVPSSTP::initThermo(), StoichSubstance::initThermo(), VPStandardStateTP::installPDSS(), LatticePhase::setSiteDensity(), PengRobinson::setSpeciesCoeffs(), RedlichKwongMFTP::setSpeciesCoeffs(), and Cantera::setupPhase().
shared_ptr< Species > species | ( | size_t | k | ) | const |
Return the Species object for species whose index is k.
Changes to this object do not affect the ThermoPhase object until the modifySpecies function is called.
Definition at line 961 of file Phase.cpp.
References Phase::checkSpeciesIndex(), and Phase::m_speciesNames.
void ignoreUndefinedElements | ( | ) |
Set behavior when adding a species containing undefined elements to just skip the species.
Definition at line 967 of file Phase.cpp.
References Phase::m_undefinedElementBehavior.
Referenced by Cantera::importPhase(), and Cantera::setupPhase().
void addUndefinedElements | ( | ) |
Set behavior when adding a species containing undefined elements to add those elements to the phase.
This is the default behavior.
Definition at line 971 of file Phase.cpp.
References Phase::m_undefinedElementBehavior.
Referenced by Cantera::setupPhase().
void throwUndefinedElements | ( | ) |
Set the behavior when adding a species containing undefined elements to throw an exception.
Definition at line 975 of file Phase.cpp.
References Phase::m_undefinedElementBehavior.
Referenced by Cantera::importPhase(), and Cantera::setupPhase().
|
virtual |
Returns a bool indicating whether the object is ready for use.
Reimplemented in BinarySolutionTabulatedThermo.
Definition at line 979 of file Phase.cpp.
References Phase::m_kk.
Referenced by IdealSolidSolnPhase::addSpecies().
|
inline |
Return the State Mole Fraction Number.
Definition at line 858 of file Phase.h.
References Phase::m_stateNum.
Referenced by HMWSoln::calcDensity(), HMWSoln::s_update_d2lnMolalityActCoeff_dT2(), HMWSoln::s_update_dlnMolalityActCoeff_dP(), HMWSoln::s_update_dlnMolalityActCoeff_dT(), InterfaceData::update(), BlowersMaselData::update(), and FalloffData::update().
|
virtual |
Invalidate any cached values which are normally updated only when a change in state is detected.
Reimplemented in ThermoPhase, and VPStandardStateTP.
Definition at line 984 of file Phase.cpp.
References ValueCache::clear(), Phase::m_cache, and Phase::m_stateNum.
Referenced by Phase::addSpecies(), ThermoPhase::invalidateCache(), and Phase::modifySpecies().
|
inline |
Returns true
if case sensitive species names are enforced.
Definition at line 867 of file Phase.h.
References Phase::m_caseSensitiveSpecies.
|
inline |
Set flag that determines whether case sensitive species are enforced in look-up operations, for example speciesIndex.
Definition at line 873 of file Phase.h.
References Phase::m_caseSensitiveSpecies.
|
virtual |
Set root Solution holding all phase information.
Definition at line 1003 of file Phase.cpp.
References Cantera::warn_deprecated().
vector_fp getCompositionFromMap | ( | const compositionMap & | comp | ) | const |
Converts a compositionMap to a vector with entries for each species Species that are not specified are set to zero in the vector.
[in] | comp | compositionMap containing the mixture composition |
Definition at line 1008 of file Phase.cpp.
References Phase::m_kk, Cantera::npos, and Phase::speciesIndex().
Referenced by ThermoPhase::equivalenceRatio(), ThermoPhase::mixtureFraction(), ThermoPhase::setEquivalenceRatio(), Phase::setMassFractionsByName(), ThermoPhase::setMixtureFraction(), Phase::setMoleFractionsByName(), and ThermoPhase::stoichAirFuelRatio().
void massFractionsToMoleFractions | ( | const double * | Y, |
double * | X | ||
) | const |
Converts a mixture composition from mole fractions to mass fractions.
[in] | Y | mixture composition in mass fractions (length m_kk) |
[out] | X | mixture composition in mole fractions (length m_kk) |
Definition at line 1022 of file Phase.cpp.
References Phase::m_kk, and Phase::m_molwts.
void moleFractionsToMassFractions | ( | const double * | X, |
double * | Y | ||
) | const |
Converts a mixture composition from mass fractions to mole fractions.
[in] | X | mixture composition in mole fractions (length m_kk) |
[out] | Y | mixture composition in mass fractions (length m_kk) |
Definition at line 1037 of file Phase.cpp.
References Cantera::dot(), Phase::m_kk, and Phase::m_molwts.
|
inlineprotected |
Ensure that phase is compressible.
An error is raised if the state is incompressible
setter | name of setter (used for exception handling) |
Definition at line 903 of file Phase.h.
References Phase::isCompressible(), Phase::name(), and Phase::type().
Referenced by Phase::setConcentrations(), Phase::setConcentrationsNoNorm(), Phase::setDensity(), Phase::setMolarDensity(), ThermoPhase::setState_SV(), and ThermoPhase::setState_UV().
|
protected |
Set the internally stored constant density (kg/m^3) of the phase.
Used for incompressible phases where the density is not an independent variable, that is, density does not affect pressure in state calculations.
[in] | density_ | density (kg/m^3). |
Definition at line 698 of file Phase.cpp.
References Phase::m_dens.
Referenced by BinarySolutionTabulatedThermo::calcDensity(), DebyeHuckel::calcDensity(), GibbsExcessVPSSTP::calcDensity(), HMWSoln::calcDensity(), IdealMolalSoln::calcDensity(), IdealSolidSolnPhase::calcDensity(), IdealSolnGasVPSS::calcDensity(), IonsFromNeutralVPSSTP::calcDensity(), LatticePhase::calcDensity(), LatticeSolidPhase::calcDensity(), MetalPhase::initThermo(), StoichSubstance::initThermo(), StoichSubstance::initThermoXML(), StoichSubstance::setParameters(), MetalPhase::setParametersFromXML(), and StoichSubstance::setParametersFromXML().
|
protected |
Set the molecular weight of a single species to a given value.
Used by phases where the equation of state is defined for a specific value of the molecular weight which may not exactly correspond to the value computed from the chemical formula.
k | id of the species |
mw | Molecular Weight (kg kmol-1) |
Definition at line 989 of file Phase.cpp.
References Phase::m_mmw, Phase::m_molwts, Phase::m_rmolwts, Phase::m_y, and Phase::m_ym.
Referenced by WaterSSTP::initThermo().
|
protectedvirtual |
Apply changes to the state which are needed after the composition changes.
This function is called after any call to setMassFractions(), setMoleFractions(), or similar. For phases which need to execute a callback after any change to the composition, it should be done by overriding this function rather than overriding all of the composition- setting functions. Derived class implementations of compositionChanged() should call the parent class method as well.
Reimplemented in BinarySolutionTabulatedThermo, GibbsExcessVPSSTP, IdealSolidSolnPhase, IonsFromNeutralVPSSTP, LatticePhase, and MixtureFugacityTP.
Definition at line 999 of file Phase.cpp.
References Phase::m_stateNum.
Referenced by GibbsExcessVPSSTP::compositionChanged(), IdealSolidSolnPhase::compositionChanged(), LatticePhase::compositionChanged(), MixtureFugacityTP::compositionChanged(), Phase::restoreState(), Phase::setConcentrations(), Phase::setConcentrationsNoNorm(), Phase::setMassFractions(), Phase::setMassFractions_NoNorm(), Phase::setMoleFractions(), and Phase::setMoleFractions_NoNorm().
|
private |
Find lowercase species name in m_speciesIndices when case sensitive species names are not enforced and a user specifies a non-canonical species name.
Raise exception if lowercase name is not unique.
Definition at line 168 of file Phase.cpp.
References Phase::m_speciesLower, Phase::name(), Cantera::npos, and Cantera::toLowerCopy().
Referenced by Phase::speciesIndex().
|
mutableprotected |
Cached for saved calculations within each ThermoPhase.
For more information on how to use this, see examples within the source code and documentation for this within ValueCache class itself.
Definition at line 923 of file Phase.h.
Referenced by HMWSoln::calcDensity(), Phase::invalidateCache(), HMWSoln::s_update_d2lnMolalityActCoeff_dT2(), HMWSoln::s_update_dlnMolalityActCoeff_dP(), HMWSoln::s_update_dlnMolalityActCoeff_dT(), IdealGasPhase::updateThermo(), and PlasmaPhase::updateThermo().
|
protected |
Number of species in the phase.
Definition at line 943 of file Phase.h.
Referenced by DebyeHuckel::_lnactivityWaterHelgesonFixedForm(), MixtureFugacityTP::_updateReferenceStateThermo(), VPStandardStateTP::_updateStandardStateThermo(), BinarySolutionTabulatedThermo::_updateThermo(), IdealSolidSolnPhase::_updateThermo(), LatticePhase::_updateThermo(), SurfPhase::_updateThermo(), Phase::addElement(), BinarySolutionTabulatedThermo::addSpecies(), GibbsExcessVPSSTP::addSpecies(), IdealGasPhase::addSpecies(), IdealSolidSolnPhase::addSpecies(), IonsFromNeutralVPSSTP::addSpecies(), LatticePhase::addSpecies(), MixtureFugacityTP::addSpecies(), MolalityVPSSTP::addSpecies(), PengRobinson::addSpecies(), Phase::addSpecies(), PlasmaPhase::addSpecies(), RedlichKwongMFTP::addSpecies(), SingleSpeciesTP::addSpecies(), SurfPhase::addSpecies(), ThermoPhase::addSpecies(), VPStandardStateTP::addSpecies(), HMWSoln::applyphScale(), GibbsExcessVPSSTP::calcDensity(), IonsFromNeutralVPSSTP::calcIonMoleFractions(), MolalityVPSSTP::calcMolalities(), HMWSoln::calcMolalitiesCropped(), IonsFromNeutralVPSSTP::calcNeutralMoleculeMoleFractions(), PengRobinson::calculateAB(), RedlichKwongMFTP::calculateAB(), Phase::chargeDensity(), GibbsExcessVPSSTP::checkMFSum(), Phase::checkSpeciesArraySize(), Phase::checkSpeciesIndex(), HMWSoln::counterIJ_setup(), RedlichKisterVPSSTP::cp_mole(), PengRobinson::d2aAlpha_dT2(), Phase::elementalMassFraction(), Phase::elementalMoleFraction(), RedlichKisterVPSSTP::enthalpy_mole(), RedlichKisterVPSSTP::entropy_mole(), SurfPhase::entropy_mole(), MolalityVPSSTP::findCLMIndex(), DebyeHuckel::getActivities(), GibbsExcessVPSSTP::getActivities(), HMWSoln::getActivities(), IdealMolalSoln::getActivities(), PengRobinson::getActivityCoefficients(), GibbsExcessVPSSTP::getActivityCoefficients(), IdealGasPhase::getActivityCoefficients(), IdealSolidSolnPhase::getActivityCoefficients(), IdealSolnGasVPSS::getActivityCoefficients(), IonsFromNeutralVPSSTP::getActivityCoefficients(), LatticePhase::getActivityCoefficients(), LatticeSolidPhase::getActivityCoefficients(), MolalityVPSSTP::getActivityCoefficients(), RedlichKwongMFTP::getActivityCoefficients(), ThermoPhase::getActivityCoefficients(), MixtureFugacityTP::getActivityConcentrations(), DebyeHuckel::getActivityConcentrations(), HMWSoln::getActivityConcentrations(), IdealMolalSoln::getActivityConcentrations(), IdealSolidSolnPhase::getActivityConcentrations(), IdealSolnGasVPSS::getActivityConcentrations(), PengRobinson::getChemPotentials(), DebyeHuckel::getChemPotentials(), HMWSoln::getChemPotentials(), IdealGasPhase::getChemPotentials(), IdealMolalSoln::getChemPotentials(), IdealSolidSolnPhase::getChemPotentials(), IdealSolnGasVPSS::getChemPotentials(), LatticePhase::getChemPotentials(), MargulesVPSSTP::getChemPotentials(), RedlichKisterVPSSTP::getChemPotentials(), RedlichKwongMFTP::getChemPotentials(), SurfPhase::getChemPotentials(), IdealSolidSolnPhase::getChemPotentials_RT(), MixtureFugacityTP::getChemPotentials_RT(), RedlichKwongMFTP::getChemPotentials_RT(), VPStandardStateTP::getChemPotentials_RT(), Phase::getCompositionFromMap(), SurfPhase::getCoverages(), IdealSolidSolnPhase::getCp_R_ref(), MargulesVPSSTP::getd2lnActCoeffdT2(), RedlichKisterVPSSTP::getd2lnActCoeffdT2(), IonsFromNeutralVPSSTP::getdlnActCoeffdlnN(), MargulesVPSSTP::getdlnActCoeffdlnN(), RedlichKisterVPSSTP::getdlnActCoeffdlnN(), ThermoPhase::getdlnActCoeffdlnN(), IonsFromNeutralVPSSTP::getdlnActCoeffdlnN_diag(), MargulesVPSSTP::getdlnActCoeffdlnN_diag(), RedlichKisterVPSSTP::getdlnActCoeffdlnN_diag(), IonsFromNeutralVPSSTP::getdlnActCoeffdlnX_diag(), MargulesVPSSTP::getdlnActCoeffdlnX_diag(), RedlichKisterVPSSTP::getdlnActCoeffdlnX_diag(), IonsFromNeutralVPSSTP::getdlnActCoeffds(), RedlichKisterVPSSTP::getdlnActCoeffds(), MargulesVPSSTP::getdlnActCoeffdT(), RedlichKisterVPSSTP::getdlnActCoeffdT(), ThermoPhase::getElectrochemPotentials(), IdealSolidSolnPhase::getEnthalpy_RT(), LatticePhase::getEnthalpy_RT(), IdealSolidSolnPhase::getEnthalpy_RT_ref(), IdealGasPhase::getEntropy_R(), MixtureFugacityTP::getEntropy_R(), IdealSolidSolnPhase::getEntropy_R_ref(), IdealSolidSolnPhase::getGibbs_ref(), LatticePhase::getGibbs_ref(), LatticeSolidPhase::getGibbs_ref(), VPStandardStateTP::getGibbs_ref(), WaterSSTP::getGibbs_ref(), IdealGasPhase::getGibbs_RT(), IdealSolidSolnPhase::getGibbs_RT(), LatticePhase::getGibbs_RT(), MixtureFugacityTP::getGibbs_RT(), IdealSolidSolnPhase::getGibbs_RT_ref(), LatticePhase::getGibbs_RT_ref(), IdealGasPhase::getIntEnergy_RT(), IdealSolidSolnPhase::getIntEnergy_RT(), MixtureFugacityTP::getIntEnergy_RT(), VPStandardStateTP::getIntEnergy_RT(), IdealGasPhase::getIntEnergy_RT_ref(), IdealSolidSolnPhase::getIntEnergy_RT_ref(), MargulesVPSSTP::getLnActivityCoefficients(), RedlichKisterVPSSTP::getLnActivityCoefficients(), ThermoPhase::getLnActivityCoefficients(), Phase::getMassFractionsByName(), MolalityVPSSTP::getMolalities(), DebyeHuckel::getMolalityActivityCoefficients(), IdealMolalSoln::getMolalityActivityCoefficients(), Phase::getMoleFractionsByName(), IonsFromNeutralVPSSTP::getNeutralMoleculeMoleGrads(), ThermoPhase::getParameters(), DebyeHuckel::getPartialMolarCp(), HMWSoln::getPartialMolarCp(), IdealMolalSoln::getPartialMolarCp(), IdealSolidSolnPhase::getPartialMolarCp(), IdealSolnGasVPSS::getPartialMolarCp(), LatticePhase::getPartialMolarCp(), SurfPhase::getPartialMolarCp(), PengRobinson::getPartialMolarEnthalpies(), DebyeHuckel::getPartialMolarEnthalpies(), HMWSoln::getPartialMolarEnthalpies(), IdealMolalSoln::getPartialMolarEnthalpies(), IdealSolidSolnPhase::getPartialMolarEnthalpies(), IdealSolnGasVPSS::getPartialMolarEnthalpies(), IonsFromNeutralVPSSTP::getPartialMolarEnthalpies(), MargulesVPSSTP::getPartialMolarEnthalpies(), RedlichKwongMFTP::getPartialMolarEnthalpies(), SurfPhase::getPartialMolarEnthalpies(), DebyeHuckel::getPartialMolarEntropies(), HMWSoln::getPartialMolarEntropies(), IdealGasPhase::getPartialMolarEntropies(), IdealMolalSoln::getPartialMolarEntropies(), IdealSolidSolnPhase::getPartialMolarEntropies(), IdealSolnGasVPSS::getPartialMolarEntropies(), IonsFromNeutralVPSSTP::getPartialMolarEntropies(), LatticePhase::getPartialMolarEntropies(), RedlichKwongMFTP::getPartialMolarEntropies(), SurfPhase::getPartialMolarEntropies(), PengRobinson::getPartialMolarIntEnergies(), IdealGasPhase::getPartialMolarIntEnergies(), IdealSolnGasVPSS::getPartialMolarIntEnergies(), PengRobinson::getPartialMolarVolumes(), DebyeHuckel::getPartialMolarVolumes(), HMWSoln::getPartialMolarVolumes(), IdealGasPhase::getPartialMolarVolumes(), RedlichKisterVPSSTP::getPartialMolarVolumes(), RedlichKwongMFTP::getPartialMolarVolumes(), IdealGasPhase::getPureGibbs(), IdealSolidSolnPhase::getPureGibbs(), LatticePhase::getPureGibbs(), MixtureFugacityTP::getPureGibbs(), VPStandardStateTP::getPureGibbs(), RedlichKwongMFTP::getSpeciesParameters(), IdealGasPhase::getStandardChemPotentials(), MixtureFugacityTP::getStandardChemPotentials(), VPStandardStateTP::getStandardChemPotentials(), IdealGasPhase::getStandardVolumes(), MixtureFugacityTP::getStandardVolumes(), SurfPhase::getStandardVolumes(), IdealGasPhase::getStandardVolumes_ref(), MixtureFugacityTP::getStandardVolumes_ref(), HMWSoln::getUnscaledMolalityActivityCoefficients(), HMWSoln::initLengths(), MargulesVPSSTP::initLengths(), RedlichKisterVPSSTP::initLengths(), IonsFromNeutralVPSSTP::initThermo(), RedlichKwongMFTP::initThermo(), StoichSubstance::initThermo(), ThermoPhase::initThermo(), VPStandardStateTP::initThermo(), RedlichKwongMFTP::initThermoXML(), Phase::massFractionsToMoleFractions(), Phase::moleFractionsToMassFractions(), Phase::nSpecies(), ThermoPhase::o2Present(), ThermoPhase::o2Required(), MolalityVPSSTP::osmoticCoefficient(), HMWSoln::printCoeffs(), Phase::ready(), HMWSoln::relative_enthalpy(), HMWSoln::relative_molal_enthalpy(), DebyeHuckel::s_update_d2lnMolalityActCoeff_dT2(), HMWSoln::s_update_d2lnMolalityActCoeff_dT2(), IonsFromNeutralVPSSTP::s_update_dlnActCoeff_dlnN(), IonsFromNeutralVPSSTP::s_update_dlnActCoeff_dlnN_diag(), IonsFromNeutralVPSSTP::s_update_dlnActCoeff_dlnX_diag(), MargulesVPSSTP::s_update_dlnActCoeff_dT(), RedlichKisterVPSSTP::s_update_dlnActCoeff_dT(), IonsFromNeutralVPSSTP::s_update_dlnActCoeffdT(), DebyeHuckel::s_update_dlnMolalityActCoeff_dP(), HMWSoln::s_update_dlnMolalityActCoeff_dP(), DebyeHuckel::s_update_dlnMolalityActCoeff_dT(), HMWSoln::s_update_dlnMolalityActCoeff_dT(), IonsFromNeutralVPSSTP::s_update_lnActCoeff(), DebyeHuckel::s_update_lnMolalityActCoeff(), HMWSoln::s_updateIMS_lnMolalityActCoeff(), IdealMolalSoln::s_updateIMS_lnMolalityActCoeff(), HMWSoln::s_updatePitzer_d2lnMolalityActCoeff_dT2(), HMWSoln::s_updatePitzer_dlnMolalityActCoeff_dP(), HMWSoln::s_updatePitzer_dlnMolalityActCoeff_dT(), HMWSoln::s_updatePitzer_lnMolalityActCoeff(), HMWSoln::s_updateScaling_pHScaling(), HMWSoln::s_updateScaling_pHScaling_dP(), HMWSoln::s_updateScaling_pHScaling_dT(), HMWSoln::s_updateScaling_pHScaling_dT2(), RedlichKwongMFTP::setBinaryCoeffs(), Phase::setConcentrations(), Phase::setConcentrationsNoNorm(), SurfPhase::setCoverages(), SurfPhase::setCoveragesByName(), SurfPhase::setCoveragesNoNorm(), DebyeHuckel::setDefaultIonicRadius(), Phase::setMassFractions(), Phase::setMassFractions_NoNorm(), MolalityVPSSTP::setMolalities(), MolalityVPSSTP::setMolalitiesByName(), Phase::setMoleFractions(), Phase::setMoleFractions_NoNorm(), LatticePhase::setSiteDensity(), PengRobinson::setSpeciesCoeffs(), RedlichKwongMFTP::setSpeciesCoeffs(), IdealSolnGasVPSS::setToEquilState(), IdealGasPhase::setToEquilState(), IdealSolidSolnPhase::setToEquilState(), ThermoPhase::speciesData(), IdealSolidSolnPhase::standardConcentration(), Phase::sum_xlogx(), PengRobinson::updateMixingExpressions(), RedlichKwongMFTP::updateMixingExpressions(), and IdealGasPhase::updateThermo().
|
protected |
Dimensionality of the phase.
Volumetric phases have dimensionality 3 and surface phases have dimensionality 2.
Definition at line 947 of file Phase.h.
Referenced by SurfPhase::getParameters(), SurfPhase::initThermo(), Phase::nDim(), and Phase::setNDim().
|
protected |
Atomic composition of the species.
The number of atoms of element i in species k is equal to m_speciesComp[k * m_mm + i] The length of this vector is equal to m_kk * m_mm
Definition at line 952 of file Phase.h.
Referenced by Phase::addElement(), Phase::addSpecies(), Phase::getAtoms(), Phase::nAtoms(), and LatticeSolidPhase::setLatticeStoichiometry().
|
protected |
Vector of species charges. length m_kk.
Definition at line 954 of file Phase.h.
Referenced by Phase::addSpecies(), Phase::charge(), Phase::getCharges(), IonsFromNeutralVPSSTP::getDissociationCoeffs(), and DebyeHuckel::s_update_lnMolalityActCoeff().
|
protected |
|
protected |
Flag determining behavior when adding species with an undefined element.
Definition at line 959 of file Phase.h.
Referenced by Phase::addSpecies(), Phase::addUndefinedElements(), Phase::ignoreUndefinedElements(), and Phase::throwUndefinedElements().
|
protected |
Flag determining whether case sensitive species names are enforced.
Definition at line 962 of file Phase.h.
Referenced by Phase::caseSensitiveSpecies(), Phase::setCaseSensitiveSpecies(), and Phase::speciesIndex().
|
private |
XML node containing the XML info for this phase.
Definition at line 970 of file Phase.h.
Referenced by Phase::setXMLdata(), and Phase::xml().
|
private |
ID of the phase.
This is the value of the ID attribute of the XML phase node. The field will stay that way even if the name is changed.
Definition at line 974 of file Phase.h.
Referenced by Phase::setName().
|
private |
Name of the phase.
Initially, this is the name specified in the YAML or CTI input file, or the value of the ID attribute of the XML phase node. It may be changed to another value during the course of a calculation.
Definition at line 980 of file Phase.h.
Referenced by Phase::addSpecies(), Phase::name(), Phase::setName(), and Phase::speciesSPName().
|
private |
Temperature (K). This is an independent variable.
Definition at line 982 of file Phase.h.
Referenced by Phase::electronTemperature(), Phase::setTemperature(), and Phase::temperature().
|
private |
Density (kg m-3).
This is an independent variable except in the case of incompressible phases, where it has to be changed using the assignDensity() method. For compressible substances, the pressure is determined from this variable rather than other way round.
Definition at line 988 of file Phase.h.
Referenced by Phase::assignDensity(), Phase::concentration(), Phase::density(), Phase::getConcentrations(), Phase::setDensity(), and Phase::setMolarDensity().
|
private |
mean molecular weight of the mixture (kg kmol-1)
Definition at line 990 of file Phase.h.
Referenced by Phase::addSpecies(), Phase::getMoleFractions(), Phase::mean_X(), Phase::meanMolecularWeight(), Phase::moleFraction(), Phase::setConcentrations(), Phase::setConcentrationsNoNorm(), Phase::setMassFractions(), Phase::setMassFractions_NoNorm(), Phase::setMolecularWeight(), Phase::setMoleFractions(), Phase::setMoleFractions_NoNorm(), and Phase::sum_xlogx().
|
mutableprivate |
m_ym[k] = mole fraction of species k divided by the mean molecular weight of mixture.
Definition at line 994 of file Phase.h.
Referenced by Phase::addSpecies(), Phase::getConcentrations(), Phase::getMoleFractions(), Phase::mean_X(), Phase::moleFractdivMMW(), Phase::moleFraction(), Phase::setConcentrations(), Phase::setConcentrationsNoNorm(), Phase::setMassFractions(), Phase::setMassFractions_NoNorm(), Phase::setMolecularWeight(), Phase::setMoleFractions(), Phase::setMoleFractions_NoNorm(), and Phase::sum_xlogx().
|
mutableprivate |
Mass fractions of the species.
Note, this vector Length is m_kk
Definition at line 1001 of file Phase.h.
Referenced by Phase::addSpecies(), Phase::concentration(), Phase::getMassFractions(), Phase::massFraction(), Phase::massFractions(), Phase::setConcentrations(), Phase::setConcentrationsNoNorm(), Phase::setMassFractions(), Phase::setMassFractions_NoNorm(), Phase::setMolecularWeight(), Phase::setMoleFractions(), and Phase::setMoleFractions_NoNorm().
|
private |
species molecular weights (kg kmol-1)
Definition at line 1003 of file Phase.h.
Referenced by Phase::addSpecies(), Phase::massFractionsToMoleFractions(), Phase::molecularWeight(), Phase::molecularWeights(), Phase::moleFractionsToMassFractions(), Phase::setConcentrations(), Phase::setConcentrationsNoNorm(), Phase::setMolecularWeight(), Phase::setMoleFractions(), and Phase::setMoleFractions_NoNorm().
|
private |
inverse of species molecular weights (kmol kg-1)
Definition at line 1005 of file Phase.h.
Referenced by Phase::addSpecies(), Phase::concentration(), Phase::setMassFractions(), Phase::setMassFractions_NoNorm(), and Phase::setMolecularWeight().
|
private |
State Change variable.
Whenever the mole fraction vector changes, this int is incremented.
Definition at line 1009 of file Phase.h.
Referenced by Phase::compositionChanged(), Phase::invalidateCache(), and Phase::stateMFNumber().
|
private |
Vector of the species names.
Definition at line 1012 of file Phase.h.
Referenced by Phase::addSpecies(), Phase::species(), Phase::speciesName(), and Phase::speciesNames().
|
private |
Map of species names to indices.
Definition at line 1015 of file Phase.h.
Referenced by Phase::addSpecies(), Phase::addSpeciesAlias(), and Phase::speciesIndex().
|
private |
Map of lower-case species names to indices.
Definition at line 1018 of file Phase.h.
Referenced by Phase::addSpecies(), and Phase::findSpeciesLower().
|
private |
Number of elements.
Definition at line 1020 of file Phase.h.
Referenced by Phase::addElement(), Phase::checkElementArraySize(), Phase::checkElementIndex(), Phase::elementIndex(), Phase::getAtoms(), Phase::nAtoms(), and Phase::nElements().
|
private |
element atomic weights (kg kmol-1)
Definition at line 1021 of file Phase.h.
Referenced by Phase::addElement(), Phase::atomicWeight(), and Phase::atomicWeights().
|
private |
element atomic numbers
Definition at line 1022 of file Phase.h.
Referenced by Phase::addElement(), and Phase::atomicNumber().
|
private |
element names
Definition at line 1023 of file Phase.h.
Referenced by Phase::addElement(), Phase::elementIndex(), Phase::elementName(), and Phase::elementNames().
|
private |
Vector of element types.
Definition at line 1024 of file Phase.h.
Referenced by Phase::addElement(), Phase::changeElementType(), and Phase::elementType().
|
private |
Entropy at 298.15 K and 1 bar of stable state pure elements (J kmol-1)
Definition at line 1027 of file Phase.h.
Referenced by Phase::addElement(), and Phase::entropyElement298().