Cantera 2.6.0
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes | List of all members
Phase Class Reference

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>

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

Public Member Functions

 Phase ()
 Default constructor. More...
 
 Phase (const Phase &)=delete
 
Phaseoperator= (const Phase &)=delete
 
XML_Nodexml () 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_fpmolecularWeights () 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_fpatomicWeights () const
 Return a read-only reference to the vector of atomic weights. More...
 
size_t nElements () const
 Number of elements. More...
 
void checkElementIndex (size_t m) const
 Check that the specified element index is in range. 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< Speciesspecies (const std::string &name) const
 Return the Species object for the named species. More...
 
shared_ptr< Speciesspecies (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_Nodem_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...
 

Detailed Description

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.

Todo:
  • Specify that the input mole, mass, and volume fraction vectors must sum to one on entry to the set state routines. Non-conforming mole/mass fraction vectors are not thermodynamically consistent. Moreover, unless we do this, the calculation of Jacobians will be altered whenever the treatment of non- conforming mole fractions is changed. Add setState functions corresponding to specifying mole numbers, which is actually what is being done (well one of the options, there are many) when non- conforming mole fractions are input. Note, we realize that most numerical Jacobian and some analytical Jacobians use non-conforming calculations. These can easily be changed to the set mole number setState functions.

Definition at line 102 of file Phase.h.

Constructor & Destructor Documentation

◆ Phase()

Phase ( )

Default constructor.

Definition at line 21 of file Phase.cpp.

◆ ~Phase()

~Phase ( )
virtual

Definition at line 37 of file Phase.cpp.

Member Function Documentation

◆ xml()

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.

Deprecated:
The XML input format is deprecated and will be removed in Cantera 3.0.

Definition at line 46 of file Phase.cpp.

References Phase::m_xml.

Referenced by TransportFactory::newTransport().

◆ setXMLdata()

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.

Parameters
xmlPhaseReference to the XML node corresponding to the phase
Deprecated:
The XML input format is deprecated and will be removed in Cantera 3.0.

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().

◆ name()

std::string name ( ) const

◆ setName()

void setName ( const std::string &  nm)

Sets the string name for the phase.

Parameters
nmString 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().

◆ type()

virtual std::string type ( ) const
inlinevirtual

◆ elementName()

string elementName ( size_t  m) const

◆ elementIndex()

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.

Parameters
nameName 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().

◆ elementNames()

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().

◆ atomicWeight()

doublereal atomicWeight ( size_t  m) const

Atomic weight of element m.

Parameters
mElement index

Definition at line 121 of file Phase.cpp.

References Phase::m_atomicWeights.

Referenced by Phase::elementalMassFraction(), ChemEquil::initialize(), and WaterSSTP::initThermo().

◆ entropyElement298()

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.

Parameters
mElement index

Definition at line 126 of file Phase.cpp.

References Phase::checkElementIndex(), and Phase::m_entropy298.

Referenced by PDSS_HKFT::LookupGe().

◆ atomicNumber()

int atomicNumber ( size_t  m) const

Atomic number of element m.

Parameters
mElement index

Definition at line 137 of file Phase.cpp.

References Phase::m_atomicNumbers.

Referenced by MultiPhase::addPhase().

◆ elementType()

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.

Parameters
mElement index
Returns
the element type

Definition at line 142 of file Phase.cpp.

References Phase::m_elem_type.

Referenced by vcs_VolPhase::transferElementsFM().

◆ changeElementType()

int changeElementType ( int  m,
int  elem_type 
)

Change the element type of the mth constraint Reassigns an element type.

Parameters
mElement index
elem_typeNew elem type to be assigned
Returns
the old element type

Definition at line 147 of file Phase.cpp.

References Phase::m_elem_type.

◆ atomicWeights()

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().

◆ nElements()

size_t nElements ( ) const

◆ checkElementIndex()

void checkElementIndex ( size_t  m) const

Check that the specified element index is in range.

Throws an exception if m is greater than nElements()-1

Definition at line 86 of file Phase.cpp.

References Phase::m_mm.

Referenced by Phase::elementalMassFraction(), Phase::elementalMoleFraction(), Phase::elementName(), Phase::entropyElement298(), and Phase::nAtoms().

◆ checkElementArraySize()

void checkElementArraySize ( size_t  mm) const

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

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

Definition at line 93 of file Phase.cpp.

References Phase::m_mm.

◆ nAtoms()

doublereal nAtoms ( size_t  k,
size_t  m 
) const

◆ getAtoms()

void getAtoms ( size_t  k,
double *  atomArray 
) const

Get a vector containing the atomic composition of species k.

Parameters
kspecies index
atomArrayvector 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.

◆ speciesIndex()

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.

Parameters
nameString name of the species. It may also be in the form phaseName:speciesName
Returns
The index of the species. If the name is not found, the value npos is returned.

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().

◆ speciesName()

string speciesName ( size_t  k) const

◆ speciesSPName()

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.

Parameters
kSpecies index within the phase
Returns
The "phaseName:speciesName" string

Definition at line 225 of file Phase.cpp.

References Phase::m_name, and Phase::speciesName().

◆ speciesNames()

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().

◆ nSpecies()

size_t nSpecies ( ) const
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().

◆ checkSpeciesIndex()

void checkSpeciesIndex ( size_t  k) const

◆ checkSpeciesArraySize()

void checkSpeciesArraySize ( size_t  kk) const

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

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

Definition at line 218 of file Phase.cpp.

References Phase::m_kk.

◆ isPure()

virtual bool isPure ( ) const
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().

◆ hasPhaseTransition()

virtual bool hasPhaseTransition ( ) const
inlinevirtual

Return whether phase represents a substance with phase transitions.

Reimplemented in PureFluidPhase.

Definition at line 294 of file Phase.h.

◆ isCompressible()

virtual bool isCompressible ( ) const
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().

◆ nativeState()

std::map< std::string, size_t > nativeState ( ) const
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().

◆ fullStates()

vector< std::string > fullStates ( ) const
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().

◆ partialStates()

vector< std::string > partialStates ( ) const
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().

◆ stateSize()

size_t stateSize ( ) const
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().

◆ saveState() [1/2]

void saveState ( vector_fp state) const

Save the current internal state of the phase.

Write to vector 'state' the current internal state.

Parameters
stateoutput 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().

◆ saveState() [2/2]

void saveState ( size_t  lenstate,
doublereal *  state 
) const
virtual

Write to array 'state' the current internal state.

Parameters
lenstatelength of the state array. Must be >= stateSize()
stateoutput 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().

◆ restoreState() [1/2]

void restoreState ( const vector_fp state)

◆ restoreState() [2/2]

void restoreState ( size_t  lenstate,
const doublereal *  state 
)
virtual

Restore the state of the phase from a previously saved state vector.

Parameters
lenstateLength of the state vector
stateVector 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().

◆ setMoleFractionsByName() [1/2]

void setMoleFractionsByName ( const compositionMap xMap)

Set the species mole fractions by name.

Species not listed by name in xMap are set to zero.

Parameters
xMapmap 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().

◆ setMoleFractionsByName() [2/2]

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.

Parameters
xstring x in the form of a composition map

Definition at line 386 of file Phase.cpp.

References Cantera::parseCompString(), and Phase::setMoleFractionsByName().

◆ setMassFractionsByName() [1/2]

void setMassFractionsByName ( const compositionMap yMap)

Set the species mass fractions by name.

Species not listed by name in yMap are set to zero.

Parameters
yMapmap 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().

◆ setMassFractionsByName() [2/2]

void setMassFractionsByName ( const std::string &  x)

Set the species mass fractions by name.

Species not listed by name in x are set to zero.

Parameters
xString containing a composition map

Definition at line 422 of file Phase.cpp.

References Cantera::parseCompString(), and Phase::setMassFractionsByName().

◆ setState_TRX() [1/2]

void setState_TRX ( doublereal  t,
doublereal  dens,
const doublereal *  x 
)

Set the internally stored temperature (K), density, and mole fractions.

Parameters
tTemperature in kelvin
densDensity (kg/m^3)
xvector of species mole fractions, length m_kk

Definition at line 427 of file Phase.cpp.

References Phase::setDensity(), Phase::setMoleFractions(), and Phase::setTemperature().

◆ setState_TRX() [2/2]

void setState_TRX ( doublereal  t,
doublereal  dens,
const compositionMap x 
)

Set the internally stored temperature (K), density, and mole fractions.

Parameters
tTemperature in kelvin
densDensity (kg/m^3)
xComposition 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().

◆ setState_TRY() [1/2]

void setState_TRY ( doublereal  t,
doublereal  dens,
const doublereal *  y 
)

Set the internally stored temperature (K), density, and mass fractions.

Parameters
tTemperature in kelvin
densDensity (kg/m^3)
yvector of species mass fractions, length m_kk

Definition at line 448 of file Phase.cpp.

References Phase::setDensity(), Phase::setMassFractions(), and Phase::setTemperature().

◆ setState_TRY() [2/2]

void setState_TRY ( doublereal  t,
doublereal  dens,
const compositionMap y 
)

Set the internally stored temperature (K), density, and mass fractions.

Parameters
tTemperature in kelvin
densDensity (kg/m^3)
yComposition 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().

◆ setState_TNX()

void setState_TNX ( doublereal  t,
doublereal  n,
const doublereal *  x 
)

Set the internally stored temperature (K), molar density (kmol/m^3), and mole fractions.

Parameters
tTemperature in kelvin
nmolar density (kmol/m^3)
xvector of species mole fractions, length m_kk

Definition at line 434 of file Phase.cpp.

References Phase::setMolarDensity(), Phase::setMoleFractions(), and Phase::setTemperature().

◆ setState_TR()

void setState_TR ( doublereal  t,
doublereal  rho 
)

Set the internally stored temperature (K) and density (kg/m^3)

Parameters
tTemperature in kelvin
rhoDensity (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().

◆ setState_TX()

void setState_TX ( doublereal  t,
doublereal *  x 
)

Set the internally stored temperature (K) and mole fractions.

Parameters
tTemperature in kelvin
xvector of species mole fractions, length m_kk

Definition at line 468 of file Phase.cpp.

References Phase::setMoleFractions(), and Phase::setTemperature().

◆ setState_TY()

void setState_TY ( doublereal  t,
doublereal *  y 
)

Set the internally stored temperature (K) and mass fractions.

Parameters
tTemperature in kelvin
yvector of species mass fractions, length m_kk

Definition at line 474 of file Phase.cpp.

References Phase::setMassFractions(), and Phase::setTemperature().

◆ setState_RX()

void setState_RX ( doublereal  rho,
doublereal *  x 
)

Set the density (kg/m^3) and mole fractions.

Parameters
rhoDensity (kg/m^3)
xvector of species mole fractions, length m_kk

Definition at line 480 of file Phase.cpp.

References Phase::setDensity(), and Phase::setMoleFractions().

◆ setState_RY()

void setState_RY ( doublereal  rho,
doublereal *  y 
)

Set the density (kg/m^3) and mass fractions.

Parameters
rhoDensity (kg/m^3)
yvector of species mass fractions, length m_kk

Definition at line 486 of file Phase.cpp.

References Phase::setDensity(), and Phase::setMassFractions().

◆ molecularWeight()

doublereal molecularWeight ( size_t  k) const

◆ getMolecularWeights() [1/2]

void getMolecularWeights ( vector_fp weights) const

Copy the vector of molecular weights into vector weights.

Parameters
weightsOutput vector of molecular weights (kg/kmol)

Definition at line 498 of file Phase.cpp.

References Phase::molecularWeights().

◆ getMolecularWeights() [2/2]

void getMolecularWeights ( doublereal *  weights) const

Copy the vector of molecular weights into array weights.

Parameters
weightsOutput array of molecular weights (kg/kmol)

Definition at line 503 of file Phase.cpp.

References Phase::molecularWeights().

◆ molecularWeights()

const vector_fp & molecularWeights ( ) const

◆ getCharges()

void getCharges ( double *  charges) const

Copy the vector of species charges into array charges.

Parameters
chargesOutput array of species charges (elem. charge)

Definition at line 514 of file Phase.cpp.

References Phase::m_speciesCharge.

◆ getMoleFractionsByName()

compositionMap getMoleFractionsByName ( double  threshold = 0.0) const

Get the mole fractions by name.

Parameters
thresholdExclude species with mole fractions less than or equal to this threshold.
Returns
Map of species names to mole fractions

Definition at line 519 of file Phase.cpp.

References Phase::m_kk, Phase::moleFraction(), and Phase::speciesName().

◆ moleFraction() [1/2]

double moleFraction ( size_t  k) const

Return the mole fraction of a single species.

Parameters
kspecies index
Returns
Mole fraction of the species

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().

◆ moleFraction() [2/2]

double moleFraction ( const std::string &  name) const

Return the mole fraction of a single species.

Parameters
nameString name of the species
Returns
Mole fraction of the species

Definition at line 554 of file Phase.cpp.

References Phase::moleFraction(), Cantera::npos, and Phase::speciesIndex().

◆ getMassFractionsByName()

compositionMap getMassFractionsByName ( double  threshold = 0.0) const

Get the mass fractions by name.

Parameters
thresholdExclude species with mass fractions less than or equal to this threshold.
Returns
Map of species names to mass fractions

Definition at line 531 of file Phase.cpp.

References Phase::m_kk, Phase::massFraction(), and Phase::speciesName().

◆ massFraction() [1/2]

double massFraction ( size_t  k) const

Return the mass fraction of a single species.

Parameters
kspecies index
Returns
Mass fraction of the species

Definition at line 569 of file Phase.cpp.

References Phase::checkSpeciesIndex(), and Phase::m_y.

Referenced by Phase::elementalMassFraction(), Phase::getMassFractionsByName(), and ThermoPhase::getParameters().

◆ massFraction() [2/2]

double massFraction ( const std::string &  name) const

Return the mass fraction of a single species.

Parameters
nameString name of the species
Returns
Mass Fraction of the species

Definition at line 575 of file Phase.cpp.

References Phase::massFractions(), Cantera::npos, and Phase::speciesIndex().

◆ getMoleFractions()

void getMoleFractions ( double *const  x) const

◆ setMoleFractions()

void setMoleFractions ( const double *const  x)
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.

Parameters
xArray 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().

◆ setMoleFractions_NoNorm()

void setMoleFractions_NoNorm ( const double *const  x)
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.

Parameters
xInput 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().

◆ getMassFractions()

void getMassFractions ( double *const  y) const

◆ massFractions()

const double * massFractions ( ) const
inline

◆ setMassFractions()

void setMassFractions ( const double *const  y)
virtual

Set the mass fractions to the specified values and normalize them.

Parameters
[in]yArray 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().

◆ setMassFractions_NoNorm()

void setMassFractions_NoNorm ( const double *const  y)
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.

Parameters
yInput 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().

◆ getConcentrations()

void getConcentrations ( double *const  c) const

Get the species concentrations (kmol/m^3).

Parameters
[out]cThe 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().

◆ concentration()

double concentration ( const size_t  k) const

Concentration of species k.

If k is outside the valid range, an exception will be thrown.

Parameters
[in]kIndex of the species within the phase.
Returns
the concentration of species k (kmol m-3).

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().

◆ setConcentrations()

void setConcentrations ( const double *const  conc)
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.

Parameters
[in]concArray 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().

◆ setConcentrationsNoNorm()

void setConcentrationsNoNorm ( const double *const  conc)
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().

◆ elementalMassFraction()

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\).

Parameters
[in]mIndex of the element within the phase. If m is outside the valid range, an exception will be thrown.
Returns
the elemental mass fraction of element m.

Definition at line 642 of file Phase.cpp.

References Phase::atomicWeight(), Phase::checkElementIndex(), Phase::m_kk, Phase::massFraction(), Phase::molecularWeight(), and Phase::nAtoms().

◆ elementalMoleFraction()

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\).

Parameters
[in]mIndex of the element within the phase. If m is outside the valid range, an exception will be thrown.
Returns
the elemental mole fraction of element m.

Definition at line 653 of file Phase.cpp.

References Phase::checkElementIndex(), Phase::m_kk, Phase::moleFraction(), Phase::nAtoms(), and Phase::nElements().

◆ moleFractdivMMW()

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().

◆ charge()

doublereal charge ( size_t  k) const
inline

◆ chargeDensity()

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().

◆ nDim()

size_t nDim ( ) const
inline

◆ setNDim()

void setNDim ( size_t  ndim)
inline

Set the number of spatial dimensions (1, 2, or 3).

The number of spatial dimensions is used for vector involving directions.

Parameters
ndimInput number of dimensions.

Definition at line 645 of file Phase.h.

References Phase::m_ndim.

Referenced by EdgePhase::EdgePhase(), Cantera::importPhase(), and SurfPhase::SurfPhase().

◆ temperature()

doublereal temperature ( ) const
inline

Temperature (K).

Returns
The temperature of the phase

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().

◆ electronTemperature()

virtual double electronTemperature ( ) const
inlinevirtual

Electron Temperature (K)

Returns
The electron temperature of the phase

Reimplemented in PlasmaPhase.

Definition at line 660 of file Phase.h.

References Phase::m_temp.

◆ pressure()

virtual double pressure ( ) const
inlinevirtual

◆ density()

virtual double density ( ) const
inlinevirtual

◆ molarDensity()

double molarDensity ( ) const

◆ molarVolume()

double molarVolume ( ) const

◆ setDensity()

void setDensity ( const double  density_)
virtual

◆ setMolarDensity()

void setMolarDensity ( const double  molarDensity)
virtual

Set the internally stored molar density (kmol/m^3) of the phase.

Parameters
[in]molarDensityInput 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().

◆ setPressure()

virtual void setPressure ( double  p)
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.

Parameters
pinput 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().

◆ setTemperature()

virtual void setTemperature ( double  temp)
inlinevirtual

◆ setElectronTemperature()

virtual void setElectronTemperature ( double  etemp)
inlinevirtual

Set the internally stored electron temperature of the phase (K).

Parameters
etempElectron temperature in Kelvin

Reimplemented in PlasmaPhase.

Definition at line 730 of file Phase.h.

References Phase::type().

◆ mean_X() [1/2]

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.

Parameters
[in]QArray of length m_kk that is to be averaged.
Returns
mole-fraction-weighted mean of Q

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().

◆ mean_X() [2/2]

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.

Parameters
[in]QArray of length m_kk that is to be averaged.
Returns
mole-fraction-weighted mean of Q

Definition at line 722 of file Phase.cpp.

References Phase::m_mmw, and Phase::m_ym.

◆ meanMolecularWeight()

doublereal meanMolecularWeight ( ) const
inline

◆ sum_xlogx()

doublereal sum_xlogx ( ) const

◆ addElement()

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.

Parameters
symbolAtomic symbol std::string.
weightAtomic mass in amu.
atomicNumberAtomic number of the element (unitless)
entropy298Entropy 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_typeSpecifies the type of the element constraint equation. This defaults to CT_ELEM_TYPE_ABSPOS, that is, an element.
Returns
index of the element added

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().

◆ addSpecies()

bool addSpecies ( shared_ptr< Species spec)
virtual

Add a Species to this Phase.

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.

See also
ignoreUndefinedElements addUndefinedElements throwUndefinedElements

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().

◆ modifySpecies()

void modifySpecies ( size_t  k,
shared_ptr< Species spec 
)
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().

◆ addSpeciesAlias()

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.

Parameters
nameoriginal species name std::string.
aliasalternate name std::string.
Returns
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().

◆ findIsomers() [1/2]

vector< std::string > findIsomers ( const compositionMap compMap) const
virtual

Return a vector with isomers names matching a given composition map.

Parameters
compMapcompositionMap of the species.
Returns
A vector of species names for matching species.

Definition at line 932 of file Phase.cpp.

Referenced by Phase::findIsomers().

◆ findIsomers() [2/2]

vector< std::string > findIsomers ( const std::string &  comp) const
virtual

Return a vector with isomers names matching a given composition string.

Parameters
compString containing a composition map
Returns
A vector of species names for matching species.

Definition at line 945 of file Phase.cpp.

References Phase::findIsomers(), and Cantera::parseCompString().

◆ species() [1/2]

shared_ptr< Species > species ( const std::string &  name) const

◆ species() [2/2]

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.

◆ ignoreUndefinedElements()

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().

◆ addUndefinedElements()

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().

◆ throwUndefinedElements()

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().

◆ ready()

bool ready ( ) const
virtual

Returns a bool indicating whether the object is ready for use.

Returns
true if the object is ready for calculation, false otherwise.

Reimplemented in BinarySolutionTabulatedThermo.

Definition at line 979 of file Phase.cpp.

References Phase::m_kk.

Referenced by IdealSolidSolnPhase::addSpecies().

◆ stateMFNumber()

int stateMFNumber ( ) const
inline

◆ invalidateCache()

void invalidateCache ( )
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().

◆ caseSensitiveSpecies()

bool caseSensitiveSpecies ( ) const
inline

Returns true if case sensitive species names are enforced.

Definition at line 867 of file Phase.h.

References Phase::m_caseSensitiveSpecies.

◆ setCaseSensitiveSpecies()

void setCaseSensitiveSpecies ( bool  cflag = true)
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.

◆ setRoot()

void setRoot ( std::shared_ptr< Solution root)
virtual

Set root Solution holding all phase information.

Deprecated:
This function has no effect. To be removed after Cantera 2.6.

Definition at line 1003 of file Phase.cpp.

References Cantera::warn_deprecated().

◆ getCompositionFromMap()

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.

Parameters
[in]compcompositionMap containing the mixture composition
Returns
vector with length m_kk

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().

◆ massFractionsToMoleFractions()

void massFractionsToMoleFractions ( const double *  Y,
double *  X 
) const

Converts a mixture composition from mole fractions to mass fractions.

Parameters
[in]Ymixture composition in mass fractions (length m_kk)
[out]Xmixture composition in mole fractions (length m_kk)

Definition at line 1022 of file Phase.cpp.

References Phase::m_kk, and Phase::m_molwts.

◆ moleFractionsToMassFractions()

void moleFractionsToMassFractions ( const double *  X,
double *  Y 
) const

Converts a mixture composition from mass fractions to mole fractions.

Parameters
[in]Xmixture composition in mole fractions (length m_kk)
[out]Ymixture 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.

◆ assertCompressible()

void assertCompressible ( const std::string &  setter) const
inlineprotected

Ensure that phase is compressible.

An error is raised if the state is incompressible

Parameters
settername 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().

◆ assignDensity()

void assignDensity ( const double  density_)
protected

◆ setMolecularWeight()

void setMolecularWeight ( const int  k,
const double  mw 
)
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.

Parameters
kid of the species
mwMolecular 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().

◆ compositionChanged()

void compositionChanged ( )
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().

◆ findSpeciesLower()

size_t findSpeciesLower ( const std::string &  nameStr) const
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().

Member Data Documentation

◆ m_cache

ValueCache m_cache
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().

◆ m_kk

size_t m_kk
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().

◆ m_ndim

size_t m_ndim
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().

◆ m_speciesComp

vector_fp m_speciesComp
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().

◆ m_speciesCharge

vector_fp m_speciesCharge
protected

◆ m_species

std::map<std::string, shared_ptr<Species> > m_species
protected

Definition at line 956 of file Phase.h.

◆ m_undefinedElementBehavior

UndefElement::behavior m_undefinedElementBehavior
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().

◆ m_caseSensitiveSpecies

bool m_caseSensitiveSpecies
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().

◆ m_xml

XML_Node* m_xml
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().

◆ m_id

std::string m_id
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().

◆ m_name

std::string m_name
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().

◆ m_temp

doublereal m_temp
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().

◆ m_dens

doublereal m_dens
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().

◆ m_mmw

doublereal m_mmw
private

◆ m_ym

vector_fp m_ym
mutableprivate

◆ m_y

vector_fp m_y
mutableprivate

◆ m_molwts

vector_fp m_molwts
private

◆ m_rmolwts

vector_fp m_rmolwts
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().

◆ m_stateNum

int m_stateNum
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().

◆ m_speciesNames

std::vector<std::string> m_speciesNames
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().

◆ m_speciesIndices

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

Map of species names to indices.

Definition at line 1015 of file Phase.h.

Referenced by Phase::addSpecies(), Phase::addSpeciesAlias(), and Phase::speciesIndex().

◆ m_speciesLower

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

Map of lower-case species names to indices.

Definition at line 1018 of file Phase.h.

Referenced by Phase::addSpecies(), and Phase::findSpeciesLower().

◆ m_mm

size_t m_mm
private

◆ m_atomicWeights

vector_fp m_atomicWeights
private

element atomic weights (kg kmol-1)

Definition at line 1021 of file Phase.h.

Referenced by Phase::addElement(), Phase::atomicWeight(), and Phase::atomicWeights().

◆ m_atomicNumbers

vector_int m_atomicNumbers
private

element atomic numbers

Definition at line 1022 of file Phase.h.

Referenced by Phase::addElement(), and Phase::atomicNumber().

◆ m_elementNames

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

element names

Definition at line 1023 of file Phase.h.

Referenced by Phase::addElement(), Phase::elementIndex(), Phase::elementName(), and Phase::elementNames().

◆ m_elem_type

vector_int m_elem_type
private

Vector of element types.

Definition at line 1024 of file Phase.h.

Referenced by Phase::addElement(), Phase::changeElementType(), and Phase::elementType().

◆ m_entropy298

vector_fp m_entropy298
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().


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