Cantera  3.1.0a1
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]

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, and then temperature and pressure. An example of this is the function Phase::setState_TPY(double t, double p, const double* y).

For bulk (3-dimensional) phases, the mass density has units of kg/m^3, and the molar density and concentrations have units of kmol/m^3, and the units listed in the methods of the Phase class assume a bulk phase. However, for surface (2-dimensional) phases have units of kg/m^2 and kmol/m^2, respectively. And for edge (1-dimensional) phases, these units kg/m and kmol/m.

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 94 of file Phase.h.

Public Member Functions

 Phase ()=default
 Default constructor. More...
 
 Phase (const Phase &)=delete
 
Phaseoperator= (const Phase &)=delete
 
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 map< string, size_t > nativeState () const
 Return a map of properties defining the native state of a substance. More...
 
string nativeMode () const
 Return string acronym representing the native state of a Phase. More...
 
virtual vector< string > fullStates () const
 Return a vector containing full states defining a phase. More...
 
virtual vector< 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< double > &state) const
 Save the current internal state of the phase. More...
 
virtual void saveState (size_t lenstate, double *state) const
 Write to array 'state' the current internal state. More...
 
void restoreState (const vector< double > &state)
 Restore a state saved on a previous call to saveState. More...
 
virtual void restoreState (size_t lenstate, const double *state)
 Restore the state of the phase from a previously saved state vector. More...
 
double molecularWeight (size_t k) const
 Molecular weight of species k. More...
 
void getMolecularWeights (double *weights) const
 Copy the vector of molecular weights into array weights. More...
 
const vector< double > & molecularWeights () const
 Return a const reference to the internal vector of molecular weights. More...
 
const vector< double > & inverseMolecularWeights () 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...
 
virtual void setMolesNoTruncate (const double *const N)
 Set the state of the object with moles in [kmol]. More...
 
double elementalMassFraction (const size_t m) const
 Elemental mass fraction of element m. More...
 
double elementalMoleFraction (const size_t m) const
 Elemental mole fraction of element m. More...
 
double 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...
 
double 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...
 
vector< double > getCompositionFromMap (const Composition &comp) const
 Converts a Composition 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.

For phases instantiated from YAML input files, the name is the value of the corresponding key in the phase map.

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.

string name () const
 Return the name of the phase. More...
 
void setName (const string &nm)
 Sets the string name for the phase. More...
 
virtual string type () const
 String indicating the thermodynamic model implemented. More...
 
Element and Species Information
string elementName (size_t m) const
 Name of the element with index m. More...
 
size_t elementIndex (const string &name) const
 Return the index of element named 'name'. More...
 
const vector< string > & elementNames () const
 Return a read-only reference to the vector of element names. More...
 
double atomicWeight (size_t m) const
 Atomic weight of element m. More...
 
double 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< double > & atomicWeights () const
 Return a read-only reference to the vector of atomic weights. More...
 
size_t nElements () const
 Number of elements. More...
 
void checkElementIndex (size_t m) const
 Check that the specified element index is in range. More...
 
void checkElementArraySize (size_t mm) const
 Check that an array size is at least nElements(). More...
 
double nAtoms (size_t k, size_t m) const
 Number of atoms of element m in species k. More...
 
size_t speciesIndex (const string &name) const
 Returns the index of a species named 'name' within the Phase object. More...
 
string speciesName (size_t k) const
 Name of the species with index k. More...
 
const vector< 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 Composition &xMap)
 Set the species mole fractions by name. More...
 
void setMoleFractionsByName (const string &x)
 Set the mole fractions of a group of species by name. More...
 
void setMassFractionsByName (const Composition &yMap)
 Set the species mass fractions by name. More...
 
void setMassFractionsByName (const string &x)
 Set the species mass fractions by name. More...
 
void setState_TD (double t, double rho)
 Set the internally stored temperature (K) and density (kg/m^3) More...
 
Composition
Composition 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 string &name) const
 Return the mole fraction of a single species. More...
 
Composition 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 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...
 
virtual void getConcentrations (double *const c) const
 Get the species concentrations (kmol/m^3). More...
 
virtual 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
double 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...
 
virtual double molarDensity () const
 Molar density (kmol/m^3). More...
 
virtual 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 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
double mean_X (const double *const Q) const
 Evaluate the mole-fraction-weighted mean of an array Q. More...
 
double mean_X (const vector< double > &Q) const
 Evaluate the mole-fraction-weighted mean of an array Q. More...
 
double meanMolecularWeight () const
 The mean molecular weight. Units: (kg/kmol) More...
 
double sum_xlogx () const
 Evaluate \( \sum_k X_k \ln 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 string &symbol, double weight=-12345.0, int atomicNumber=0, double 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 string &name, const string &alias)
 Add a species alias (that is, a user-defined alternative species name). More...
 
virtual vector< string > findIsomers (const Composition &compMap) const
 Return a vector with isomers names matching a given composition map. More...
 
virtual vector< string > findIsomers (const string &comp) const
 Return a vector with isomers names matching a given composition string. More...
 
shared_ptr< Speciesspecies (const 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 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 = 0
 Number of species in the phase. More...
 
size_t m_ndim = 3
 Dimensionality of the phase. More...
 
vector< double > m_speciesComp
 Atomic composition of the species. More...
 
vector< double > m_speciesCharge
 Vector of species charges. length m_kk. More...
 
map< string, shared_ptr< Species > > m_species
 
UndefElement::behavior m_undefinedElementBehavior = UndefElement::add
 Flag determining behavior when adding species with an undefined element. More...
 
bool m_caseSensitiveSpecies = false
 Flag determining whether case sensitive species names are enforced. More...
 

Private Member Functions

size_t findSpeciesLower (const 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

string m_name
 Name of the phase. More...
 
double m_temp = 0.001
 Temperature (K). This is an independent variable. More...
 
double m_dens = 0.001
 Density (kg m-3). More...
 
double m_mmw = 0.0
 mean molecular weight of the mixture (kg kmol-1) More...
 
vector< double > m_ym
 m_ym[k] = mole fraction of species k divided by the mean molecular weight of mixture. More...
 
vector< double > m_y
 Mass fractions of the species. More...
 
vector< double > m_molwts
 species molecular weights (kg kmol-1) More...
 
vector< double > m_rmolwts
 inverse of species molecular weights (kmol kg-1) More...
 
int m_stateNum = -1
 State Change variable. More...
 
vector< string > m_speciesNames
 Vector of the species names. More...
 
map< string, size_t > m_speciesIndices
 Map of species names to indices. More...
 
map< string, size_t > m_speciesLower
 Map of lower-case species names to indices. More...
 
size_t m_mm = 0
 Number of elements. More...
 
vector< double > m_atomicWeights
 element atomic weights (kg kmol-1) More...
 
vector< int > m_atomicNumbers
 element atomic numbers More...
 
vector< string > m_elementNames
 element names More...
 
vector< int > m_elem_type
 Vector of element types. More...
 
vector< double > m_entropy298
 Entropy at 298.15 K and 1 bar of stable state pure elements (J kmol-1) More...
 

Constructor & Destructor Documentation

◆ Phase()

Phase ( )
default

Default constructor.

Member Function Documentation

◆ name()

string name ( ) const

Return the name of the phase.

Names are unique within a Cantera problem.

Definition at line 20 of file Phase.cpp.

◆ setName()

void setName ( const string &  nm)

Sets the string name for the phase.

Parameters
nmString name of the phase

Definition at line 25 of file Phase.cpp.

◆ type()

virtual string type ( ) const
inlinevirtual

String indicating the thermodynamic model implemented.

Usually corresponds to the name of the derived class, less any suffixes such as "Phase", TP", "VPSS", etc.

Since
Starting in Cantera 3.0, the name returned by this method corresponds to the canonical name used in the YAML input format.

Reimplemented in WaterSSTP, ThermoPhase, SurfPhase, StoichSubstance, SingleSpeciesTP, RedlichKwongMFTP, RedlichKisterVPSSTP, PureFluidPhase, PlasmaPhase, PengRobinson, MixtureFugacityTP, MetalPhase, MargulesVPSSTP, LatticeSolidPhase, LatticePhase, IdealSolnGasVPSS, IdealSolidSolnPhase, IdealMolalSoln, IdealGasPhase, HMWSoln, EdgePhase, DebyeHuckel, CoverageDependentSurfPhase, and BinarySolutionTabulatedThermo.

Definition at line 134 of file Phase.h.

◆ elementName()

string elementName ( size_t  m) const

Name of the element with index m.

Parameters
mElement index.

Definition at line 49 of file Phase.cpp.

◆ elementIndex()

size_t elementIndex ( const 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 55 of file Phase.cpp.

◆ elementNames()

const vector< string > & elementNames ( ) const

Return a read-only reference to the vector of element names.

Definition at line 65 of file Phase.cpp.

◆ atomicWeight()

double atomicWeight ( size_t  m) const

Atomic weight of element m.

Parameters
mElement index

Definition at line 70 of file Phase.cpp.

◆ entropyElement298()

double 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 75 of file Phase.cpp.

◆ atomicNumber()

int atomicNumber ( size_t  m) const

Atomic number of element m.

Parameters
mElement index

Definition at line 86 of file Phase.cpp.

◆ 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 91 of file Phase.cpp.

◆ 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 96 of file Phase.cpp.

◆ atomicWeights()

const vector< double > & atomicWeights ( ) const

Return a read-only reference to the vector of atomic weights.

Definition at line 81 of file Phase.cpp.

◆ nElements()

size_t nElements ( ) const

Number of elements.

Definition at line 30 of file Phase.cpp.

◆ 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 35 of file Phase.cpp.

◆ 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 42 of file Phase.cpp.

◆ nAtoms()

double nAtoms ( size_t  k,
size_t  m 
) const

Number of atoms of element m in species k.

Parameters
kspecies index
melement index

Definition at line 103 of file Phase.cpp.

◆ speciesIndex()

size_t speciesIndex ( const 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 129 of file Phase.cpp.

◆ speciesName()

string speciesName ( size_t  k) const

Name of the species with index k.

Parameters
kindex of the species

Definition at line 142 of file Phase.cpp.

◆ speciesNames()

const vector< string > & speciesNames ( ) const

Return a const reference to the vector of species names.

Definition at line 148 of file Phase.cpp.

◆ nSpecies()

size_t nSpecies ( ) const
inline

Returns the number of species in the phase.

Definition at line 231 of file Phase.h.

◆ checkSpeciesIndex()

void checkSpeciesIndex ( size_t  k) const

Check that the specified species index is in range.

Throws an exception if k is greater than nSpecies()-1

Definition at line 153 of file Phase.cpp.

◆ 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 160 of file Phase.cpp.

◆ isPure()

virtual bool isPure ( ) const
inlinevirtual

Return whether phase represents a pure (single species) substance.

Reimplemented in SingleSpeciesTP, and PureFluidPhase.

Definition at line 247 of file Phase.h.

◆ hasPhaseTransition()

virtual bool hasPhaseTransition ( ) const
inlinevirtual

Return whether phase represents a substance with phase transitions.

Reimplemented in PureFluidPhase.

Definition at line 252 of file Phase.h.

◆ isCompressible()

virtual bool isCompressible ( ) const
inlinevirtual

Return whether phase represents a compressible substance.

Reimplemented in VPStandardStateTP, SurfPhase, StoichSubstance, MetalPhase, LatticeSolidPhase, LatticePhase, and IdealSolidSolnPhase.

Definition at line 257 of file Phase.h.

◆ nativeState()

map< 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 LatticeSolidPhase, and LatticePhase.

Definition at line 167 of file Phase.cpp.

◆ nativeMode()

string nativeMode ( ) const

Return string acronym representing the native state of a Phase.

Examples: "TP", "TDY", "TPY".

See also
nativeState
Since
New in Cantera 3.0

Definition at line 184 of file Phase.cpp.

◆ fullStates()

vector< 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 197 of file Phase.cpp.

◆ partialStates()

vector< 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 215 of file Phase.cpp.

◆ stateSize()

size_t stateSize ( ) const
virtual

Return size of vector defining internal state of the phase.

Used by saveState() and restoreState().

Definition at line 228 of file Phase.cpp.

◆ saveState() [1/2]

void saveState ( vector< double > &  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 236 of file Phase.cpp.

◆ saveState() [2/2]

void saveState ( size_t  lenstate,
double *  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 242 of file Phase.cpp.

◆ restoreState() [1/2]

void restoreState ( const vector< double > &  state)

Restore a state saved on a previous call to saveState.

Parameters
stateState vector containing the previously saved state.

Definition at line 260 of file Phase.cpp.

◆ restoreState() [2/2]

void restoreState ( size_t  lenstate,
const double *  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 265 of file Phase.cpp.

◆ setMoleFractionsByName() [1/2]

void setMoleFractionsByName ( const Composition 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 330 of file Phase.cpp.

◆ setMoleFractionsByName() [2/2]

void setMoleFractionsByName ( const 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 336 of file Phase.cpp.

◆ setMassFractionsByName() [1/2]

void setMassFractionsByName ( const Composition 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 366 of file Phase.cpp.

◆ setMassFractionsByName() [2/2]

void setMassFractionsByName ( const 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 372 of file Phase.cpp.

◆ setState_TD()

void setState_TD ( double  t,
double  rho 
)

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

Parameters
tTemperature in kelvin
rhoDensity (kg/m^3)
Since
New in Cantera 3.0.

Definition at line 377 of file Phase.cpp.

◆ molecularWeight()

double molecularWeight ( size_t  k) const

Molecular weight of species k.

Parameters
kindex of species k
Returns
the molecular weight of species k.

Definition at line 383 of file Phase.cpp.

◆ getMolecularWeights()

void getMolecularWeights ( double *  weights) const

Copy the vector of molecular weights into array weights.

Parameters
weightsOutput array of molecular weights (kg/kmol)

Definition at line 389 of file Phase.cpp.

◆ molecularWeights()

const vector< double > & molecularWeights ( ) const

Return a const reference to the internal vector of molecular weights.

units = kg / kmol

Definition at line 395 of file Phase.cpp.

◆ inverseMolecularWeights()

const vector< double > & inverseMolecularWeights ( ) const

Return a const reference to the internal vector of molecular weights.

units = kmol / kg

Definition at line 400 of file Phase.cpp.

◆ 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 405 of file Phase.cpp.

◆ getMoleFractionsByName()

Composition 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 410 of file Phase.cpp.

◆ 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 439 of file Phase.cpp.

◆ moleFraction() [2/2]

double moleFraction ( const 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 445 of file Phase.cpp.

◆ getMassFractionsByName()

Composition 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 422 of file Phase.cpp.

◆ 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 455 of file Phase.cpp.

◆ massFraction() [2/2]

double massFraction ( const 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 461 of file Phase.cpp.

◆ getMoleFractions()

void getMoleFractions ( double *const  x) const

Get the species mole fraction vector.

Parameters
xOn return, x contains the mole fractions. Must have a length greater than or equal to the number of species.

Definition at line 434 of file Phase.cpp.

◆ 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 SingleSpeciesTP, and LatticeSolidPhase.

Definition at line 289 of file Phase.cpp.

◆ 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 321 of file Phase.cpp.

◆ getMassFractions()

void getMassFractions ( double *const  y) const

Get the species mass fractions.

Parameters
[out]yArray of mass fractions, length nSpecies()

Definition at line 471 of file Phase.cpp.

◆ massFractions()

const double* massFractions ( ) const
inline

Return a const pointer to the mass fraction array.

Definition at line 442 of file Phase.h.

◆ 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 SingleSpeciesTP, and LatticeSolidPhase.

Definition at line 341 of file Phase.cpp.

◆ 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 355 of file Phase.cpp.

◆ getConcentrations()

void getConcentrations ( double *const  c) const
virtual

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.

Reimplemented in LatticeSolidPhase.

Definition at line 482 of file Phase.cpp.

◆ concentration()

double concentration ( const size_t  k) const
virtual

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

Reimplemented in LatticeSolidPhase.

Definition at line 476 of file Phase.cpp.

◆ 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 487 of file Phase.cpp.

◆ setConcentrationsNoNorm()

void setConcentrationsNoNorm ( const double *const  conc)
virtual

Set the concentrations without ignoring negative concentrations.

Definition at line 509 of file Phase.cpp.

◆ setMolesNoTruncate()

void setMolesNoTruncate ( const double *const  N)
virtual

Set the state of the object with moles in [kmol].

Definition at line 528 of file Phase.cpp.

◆ elementalMassFraction()

double 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 547 of file Phase.cpp.

◆ elementalMoleFraction()

double 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 558 of file Phase.cpp.

◆ charge()

double charge ( size_t  k) const
inline

Dimensionless electrical charge of a single molecule of species k The charge is normalized by the the magnitude of the electron charge.

Parameters
kspecies index

Definition at line 538 of file Phase.h.

◆ chargeDensity()

double chargeDensity ( ) const

Charge density [C/m^3].

Definition at line 607 of file Phase.cpp.

◆ nDim()

size_t nDim ( ) const
inline

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

Definition at line 546 of file Phase.h.

◆ 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 553 of file Phase.h.

◆ temperature()

double temperature ( ) const
inline

Temperature (K).

Returns
The temperature of the phase

Definition at line 562 of file Phase.h.

◆ electronTemperature()

virtual double electronTemperature ( ) const
inlinevirtual

Electron Temperature (K)

Returns
The electron temperature of the phase

Reimplemented in PlasmaPhase.

Definition at line 568 of file Phase.h.

◆ pressure()

virtual double pressure ( ) const
inlinevirtual

Return the thermodynamic pressure (Pa).

This method must be overloaded in derived classes. Within Cantera, the independent variable is either density or pressure. If the state is defined by temperature, density, and mass fractions, this method should use these values to implement the mechanical equation of state \( P(T, \rho, Y_1, \dots, Y_K) \). Alternatively, it returns a stored value.

Reimplemented in WaterSSTP, VPStandardStateTP, SurfPhase, StoichSubstance, RedlichKwongMFTP, PureFluidPhase, PengRobinson, MetalPhase, LatticeSolidPhase, LatticePhase, IdealSolidSolnPhase, and IdealGasPhase.

Definition at line 580 of file Phase.h.

◆ density()

virtual double density ( ) const
inlinevirtual

Density (kg/m^3).

Returns
The density of the phase

Definition at line 587 of file Phase.h.

◆ molarDensity()

double molarDensity ( ) const
virtual

Molar density (kmol/m^3).

Returns
The molar density of the phase

Definition at line 576 of file Phase.cpp.

◆ molarVolume()

double molarVolume ( ) const
virtual

Molar volume (m^3/kmol).

Returns
The molar volume of the phase

Reimplemented in SurfPhase.

Definition at line 581 of file Phase.cpp.

◆ setDensity()

void setDensity ( const double  density_)
virtual

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

Note the density of a phase is an independent variable.

Parameters
[in]density_density (kg/m^3).

Reimplemented in PureFluidPhase, and WaterSSTP.

Definition at line 586 of file Phase.cpp.

◆ 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 MetalPhase, WaterSSTP, VPStandardStateTP, SurfPhase, StoichSubstance, PureFluidPhase, MixtureFugacityTP, LatticeSolidPhase, LatticePhase, IdealSolnGasVPSS, IdealSolidSolnPhase, and IdealGasPhase.

Definition at line 616 of file Phase.h.

◆ setTemperature()

virtual void setTemperature ( double  temp)
inlinevirtual

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

Parameters
tempTemperature in Kelvin

Reimplemented in WaterSSTP, VPStandardStateTP, MixtureFugacityTP, and PureFluidPhase.

Definition at line 623 of file Phase.h.

◆ 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 634 of file Phase.h.

◆ mean_X() [1/2]

double mean_X ( const double *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 616 of file Phase.cpp.

◆ mean_X() [2/2]

double mean_X ( const vector< double > &  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 621 of file Phase.cpp.

◆ meanMolecularWeight()

double meanMolecularWeight ( ) const
inline

The mean molecular weight. Units: (kg/kmol)

Definition at line 655 of file Phase.h.

◆ sum_xlogx()

double sum_xlogx ( ) const

Evaluate \( \sum_k X_k \ln X_k \).

Returns
The indicated sum. Dimensionless.

Definition at line 626 of file Phase.cpp.

◆ addElement()

size_t addElement ( const string &  symbol,
double  weight = -12345.0,
int  atomicNumber = 0,
double  entropy298 = ENTROPY298_UNKNOWN,
int  elem_type = CT_ELEM_TYPE_ABSPOS 
)

Add an element.

Parameters
symbolAtomic symbol 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 635 of file Phase.cpp.

◆ 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 VPStandardStateTP, ThermoPhase, SurfPhase, SingleSpeciesTP, RedlichKwongMFTP, PlasmaPhase, PengRobinson, MolalityVPSSTP, MixtureFugacityTP, LatticeSolidPhase, LatticePhase, IdealSolnGasVPSS, IdealSolidSolnPhase, IdealMolalSoln, IdealGasPhase, GibbsExcessVPSSTP, DebyeHuckel, CoverageDependentSurfPhase, and BinarySolutionTabulatedThermo.

Definition at line 701 of file Phase.cpp.

◆ 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 805 of file Phase.cpp.

◆ addSpeciesAlias()

void addSpeciesAlias ( const string &  name,
const string &  alias 
)

Add a species alias (that is, a user-defined alternative species name).

Aliases are case-sensitive.

Parameters
nameoriginal species name
aliasalternate name

Definition at line 822 of file Phase.cpp.

◆ findIsomers() [1/2]

vector< string > findIsomers ( const Composition compMap) const
virtual

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

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

Definition at line 838 of file Phase.cpp.

◆ findIsomers() [2/2]

vector< string > findIsomers ( const 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 851 of file Phase.cpp.

◆ species() [1/2]

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

Return the Species object for the named species.

Changes to this object do not affect the ThermoPhase object until the modifySpecies function is called.

Definition at line 856 of file Phase.cpp.

◆ 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 867 of file Phase.cpp.

◆ ignoreUndefinedElements()

void ignoreUndefinedElements ( )

Set behavior when adding a species containing undefined elements to just skip the species.

Definition at line 873 of file Phase.cpp.

◆ 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 877 of file Phase.cpp.

◆ throwUndefinedElements()

void throwUndefinedElements ( )

Set the behavior when adding a species containing undefined elements to throw an exception.

Definition at line 881 of file Phase.cpp.

◆ 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 885 of file Phase.cpp.

◆ stateMFNumber()

int stateMFNumber ( ) const
inline

Return the State Mole Fraction Number.

Definition at line 761 of file Phase.h.

◆ invalidateCache()

void invalidateCache ( )
virtual

Invalidate any cached values which are normally updated only when a change in state is detected.

Reimplemented in VPStandardStateTP, and ThermoPhase.

Definition at line 890 of file Phase.cpp.

◆ caseSensitiveSpecies()

bool caseSensitiveSpecies ( ) const
inline

Returns true if case sensitive species names are enforced.

Definition at line 770 of file Phase.h.

◆ 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 776 of file Phase.h.

◆ getCompositionFromMap()

vector< double > getCompositionFromMap ( const Composition comp) const

Converts a Composition to a vector with entries for each species Species that are not specified are set to zero in the vector.

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

Definition at line 909 of file Phase.cpp.

◆ 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 923 of file Phase.cpp.

◆ 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 938 of file Phase.cpp.

◆ assertCompressible()

void assertCompressible ( const 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 802 of file Phase.h.

◆ assignDensity()

void assignDensity ( const double  density_)
protected

Set the internally stored constant density (kg/m^3) of the phase.

Used for incompressible phases where the density is not an independent variable, that is, density does not affect pressure in state calculations.

Parameters
[in]density_density (kg/m^3).

Definition at line 597 of file Phase.cpp.

◆ 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 895 of file Phase.cpp.

◆ 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 SurfPhase, MixtureFugacityTP, LatticePhase, IdealSolidSolnPhase, GibbsExcessVPSSTP, and BinarySolutionTabulatedThermo.

Definition at line 905 of file Phase.cpp.

◆ findSpeciesLower()

size_t findSpeciesLower ( const 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 110 of file Phase.cpp.

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 822 of file Phase.h.

◆ m_kk

size_t m_kk = 0
protected

Number of species in the phase.

Definition at line 842 of file Phase.h.

◆ m_ndim

size_t m_ndim = 3
protected

Dimensionality of the phase.

Volumetric phases have dimensionality 3 and surface phases have dimensionality 2.

Definition at line 846 of file Phase.h.

◆ m_speciesComp

vector<double> 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 851 of file Phase.h.

◆ m_speciesCharge

vector<double> m_speciesCharge
protected

Vector of species charges. length m_kk.

Definition at line 853 of file Phase.h.

◆ m_undefinedElementBehavior

UndefElement::behavior m_undefinedElementBehavior = UndefElement::add
protected

Flag determining behavior when adding species with an undefined element.

Definition at line 858 of file Phase.h.

◆ m_caseSensitiveSpecies

bool m_caseSensitiveSpecies = false
protected

Flag determining whether case sensitive species names are enforced.

Definition at line 861 of file Phase.h.

◆ m_name

string m_name
private

Name of the phase.

Initially, this is the name specified in the YAML input file. It may be changed to another value during the course of a calculation.

Definition at line 872 of file Phase.h.

◆ m_temp

double m_temp = 0.001
private

Temperature (K). This is an independent variable.

Definition at line 874 of file Phase.h.

◆ m_dens

double m_dens = 0.001
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 880 of file Phase.h.

◆ m_mmw

double m_mmw = 0.0
private

mean molecular weight of the mixture (kg kmol-1)

Definition at line 882 of file Phase.h.

◆ m_ym

vector<double> m_ym
mutableprivate

m_ym[k] = mole fraction of species k divided by the mean molecular weight of mixture.

Definition at line 886 of file Phase.h.

◆ m_y

vector<double> m_y
mutableprivate

Mass fractions of the species.

Note, this vector Length is m_kk

Definition at line 893 of file Phase.h.

◆ m_molwts

vector<double> m_molwts
private

species molecular weights (kg kmol-1)

Definition at line 895 of file Phase.h.

◆ m_rmolwts

vector<double> m_rmolwts
private

inverse of species molecular weights (kmol kg-1)

Definition at line 897 of file Phase.h.

◆ m_stateNum

int m_stateNum = -1
private

State Change variable.

Whenever the mole fraction vector changes, this int is incremented.

Definition at line 901 of file Phase.h.

◆ m_speciesNames

vector<string> m_speciesNames
private

Vector of the species names.

Definition at line 904 of file Phase.h.

◆ m_speciesIndices

map<string, size_t> m_speciesIndices
private

Map of species names to indices.

Definition at line 907 of file Phase.h.

◆ m_speciesLower

map<string, size_t> m_speciesLower
private

Map of lower-case species names to indices.

Definition at line 910 of file Phase.h.

◆ m_mm

size_t m_mm = 0
private

Number of elements.

Definition at line 912 of file Phase.h.

◆ m_atomicWeights

vector<double> m_atomicWeights
private

element atomic weights (kg kmol-1)

Definition at line 913 of file Phase.h.

◆ m_atomicNumbers

vector<int> m_atomicNumbers
private

element atomic numbers

Definition at line 914 of file Phase.h.

◆ m_elementNames

vector<string> m_elementNames
private

element names

Definition at line 915 of file Phase.h.

◆ m_elem_type

vector<int> m_elem_type
private

Vector of element types.

Definition at line 916 of file Phase.h.

◆ m_entropy298

vector<double> m_entropy298
private

Entropy at 298.15 K and 1 bar of stable state pure elements (J kmol-1)

Definition at line 919 of file Phase.h.


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