State Class Reference
[Models of Phases of Matter]

Manages the independent variables of temperature, mass density, and species mass/mole fraction that define the thermodynamic state. More...

#include <State.h>

Inheritance diagram for State:
Inheritance graph
[legend]

List of all members.

Public Member Functions

 State ()
 Constructor.
virtual ~State ()
 Destructor.
 State (const State &right)
 Copy Constructor for the State Class.
Stateoperator= (const State &right)
 Assignment operator for the state class.
bool ready () const
 True if the number species has been set.
void stateMFChangeCalc (bool forceChange=false)
 Every time the mole fractions have changed, this routine will increment the stateMFNumber.
int stateMFNumber () const
 Return the state number.
Species Information

The only thing class State knows about the species is their molecular weights.

const array_fp & molecularWeights () const
 Return a read-only reference to the array of molecular weights.
Composition

void getMoleFractions (doublereal *const x) const
 Get the species mole fraction vector.
doublereal moleFraction (const int k) const
 The mole fraction of species k.
virtual void setMoleFractions (const doublereal *const x)
 Set the mole fractions to the specified values, and then normalize them so that they sum to 1.0.
virtual void setMoleFractions_NoNorm (const doublereal *const x)
 Set the mole fractions to the specified values without normalizing.
void getMassFractions (doublereal *const y) const
 Get the species mass fractions.
doublereal massFraction (const int k) const
 Mass fraction of species k.
virtual void setMassFractions (const doublereal *const y)
 Set the mass fractions to the specified values, and then normalize them so that they sum to 1.0.
virtual void setMassFractions_NoNorm (const doublereal *const y)
 Set the mass fractions to the specified values without normalizing.
void getConcentrations (doublereal *const c) const
 Get the species concentrations (kmol/m^3).
doublereal concentration (const int k) const
 Concentration of species k.
virtual void setConcentrations (const doublereal *const conc)
 Set the concentrations to the specified values within the phase.
const doublereal * massFractions () const
 Returns a read-only pointer to the start of the massFraction array.
const doublereal * moleFractdivMMW () const
 Returns a read-only pointer to the start of the moleFraction/MW array.
Mean Properties

doublereal mean_X (const doublereal *const Q) const
 Evaluate the mole-fraction-weighted mean of Q:

\[ \sum_k X_k Q_k. \]

Array Q should contain pure-species molar property values.

doublereal mean_Y (const doublereal *const Q) const
 Evaluate the mass-fraction-weighted mean of Q:

\[ \sum_k Y_k Q_k \]

.

doublereal meanMolecularWeight () const
 The mean molecular weight.
doublereal sum_xlogx () const
 Evaluate $ \sum_k X_k \log X_k $.
doublereal sum_xlogQ (doublereal *const Q) const
 Evaluate $ \sum_k X_k \log Q_k $.
Thermodynamic Properties

Class State only stores enough thermodynamic data to specify the state.

In addition to composition information, it stores the temperature and mass density.

doublereal temperature () const
 Temperature (K).
virtual doublereal density () const
 Density (kg/m^3).
doublereal molarDensity () const
 Molar density (kmol/m^3).
virtual void setDensity (const doublereal density)
 Set the internally storred density (kg/m^3) of the phase.
virtual void setMolarDensity (const doublereal molarDensity)
 Set the internally storred molar density (kmol/m^3) of the phase.
virtual void setTemperature (const doublereal temp)
 Set the temperature (K).

Protected Member Functions

void init (const array_fp &mw)
void setMolecularWeight (const int k, const double mw)
 Set the molecular weight of a single species to a given value.

Protected Attributes

int m_kk
 m_kk is the number of species in the phase

Private Attributes

doublereal m_temp
 Temperature.
doublereal m_dens
 Density.
doublereal m_mmw
 m_mmw is the mean molecular weight of the mixture (kg kmol-1)
array_fp m_ym
 m_ym[k] = mole fraction of species k divided by the mean molecular weight of mixture.
array_fp m_y
 m_y[k] = mass fraction of species k
array_fp m_molwts
 m_molwts[k] = molecular weight of species k (kg kmol-1)
array_fp m_rmolwts
 m_rmolwts[k] = inverse of the molecular weight of species k units = kmol kg-1.
int m_stateNum
 State Change variable.

Detailed Description

Manages the independent variables of temperature, mass density, and species mass/mole fraction that define the thermodynamic state.

Class State stores just enough information about a multicomponent solution to specify its intensive thermodynamic state. It stores values for the temperature, mass density, and an array of species mass fractions. It also stores an array of species molecular weights, which are used to convert between mole and mass representations of the composition. These are the only properties of the species that class State knows about. For efficiency in mass/mole conversion, the vector of mass fractions divided by molecular weight $ Y_k/M_k $ is also stored.

Class State is not usually used directly in application programs. Its primary use is as a base class for class Phase. Class State has no virtual methods, and none of its methods are meant to be overloaded. However, this is one exception. If the phase is incompressible, then the density must be replaced by the pressure as the independent variable. In this case, functions such as setMassFraction within the class State must actually now calculate the density (at constant T and P) instead of leaving it alone as befits an independent variable. Threfore, these type of functions are virtual functions and need to be overloaded for incompressible phases. Note, for almost incompressible phases (or phases which utilize standard states based on a T and P) this may be advantageous as well, and they need to overload these functions too.

Definition at line 62 of file State.h.


Constructor & Destructor Documentation

State (  ) 

Constructor.

Definition at line 38 of file State.cpp.

~State (  )  [virtual]

Destructor.

Since no memory is allocated by methods of this class, the destructor does nothing.

Definition at line 47 of file State.cpp.

State ( const State right  ) 

Copy Constructor for the State Class.

Parameters:
right Reference to the class to be copied.

Definition at line 51 of file State.cpp.

References State::operator=().


Member Function Documentation

doublereal concentration ( const int  k  )  const

Concentration of species k.

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

Parameters:
k Index of species

Definition at line 134 of file State.cpp.

References State::m_dens, State::m_kk, State::m_rmolwts, and State::m_y.

virtual doublereal density (  )  const [inline, virtual]
void getConcentrations ( doublereal *const   c  )  const

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

Parameters:
c On return, c contains the concentrations for all species. Array c must have a length greater than or equal to the number of species.

Definition at line 219 of file State.cpp.

References State::m_dens, State::m_ym, and Cantera::scale().

Referenced by ConstDensityThermo::getActivityCoefficients(), SurfPhase::getActivityConcentrations(), IdealSolnGasVPSS::getActivityConcentrations(), IdealGasPhase::getActivityConcentrations(), and SurfPhase::getCoverages().

void getMassFractions ( doublereal *const   y  )  const

Get the species mass fractions.

Parameters:
y On return, y contains the mass fractions. Array y must have a length greater than or equal to the number of species.
y Output vector of mass fractions. Length is m_kk.

Definition at line 235 of file State.cpp.

References State::m_y.

Referenced by ThermoPhase::report(), PureFluidPhase::report(), and Phase::saveState().

void getMoleFractions ( doublereal *const   x  )  const
void init ( const array_fp &  mw  )  [protected]

For internal use only.

Initialize. Make a local copy of the vector of molecular weights, and resize the composition arrays to the appropriate size. The only information an instance of State has about the species is their molecular weights.

Parameters:
mw Vector of molecular weights of the species.

Definition at line 244 of file State.cpp.

References Cantera::int2str(), State::m_kk, State::m_mmw, State::m_molwts, State::m_rmolwts, State::m_y, State::m_ym, and Cantera::Tiny.

Referenced by Phase::freezeSpecies().

doublereal massFraction ( const int  k  )  const

Mass fraction of species k.

If k is outside the valid range, an exception will be thrown. Note that it is somewhat more efficent to call getMassFractions if the mass fractions of all species are desired.

Parameters:
k species index

Reimplemented in Phase.

Definition at line 126 of file State.cpp.

References State::m_kk, and State::m_y.

const doublereal* massFractions (  )  const [inline]

Returns a read-only pointer to the start of the massFraction array.

The pointer returned is readonly

Returns:
returns a pointer to a vector of doubles of length m_kk.

Definition at line 242 of file State.h.

References State::m_y.

Referenced by Phase::massFraction().

doublereal mean_X ( const doublereal *const   Q  )  const
doublereal mean_Y ( const doublereal *const   Q  )  const

Evaluate the mass-fraction-weighted mean of Q:

\[ \sum_k Y_k Q_k \]

.

Parameters:
Q Array Q contains a vector of species property values in mass units.
Returns:
Return value containing the mass-fraction-weighted mean of Q.

Definition at line 227 of file State.cpp.

References Cantera::dot(), and State::m_y.

doublereal meanMolecularWeight (  )  const [inline]
doublereal molarDensity (  )  const
const array_fp& molecularWeights (  )  const [inline]

Return a read-only reference to the array of molecular weights.

Reimplemented in Phase.

Definition at line 100 of file State.h.

References State::m_molwts.

const doublereal * moleFractdivMMW (  )  const

Returns a read-only 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 215 of file State.cpp.

References State::m_ym.

Referenced by IdealSolnGasVPSS::calcDensity(), IdealSolidSolnPhase::calcDensity(), and IdealSolidSolnPhase::getActivityConcentrations().

doublereal moleFraction ( const int  k  )  const

The mole fraction of species k.

If k is ouside the valid range, an exception will be thrown. Note that it is somewhat more efficent to call getMoleFractions if the mole fractions of all species are desired.

Parameters:
k species index

Reimplemented in Phase.

Definition at line 91 of file State.cpp.

References State::m_kk, State::m_mmw, and State::m_ym.

State & operator= ( const State right  ) 

Assignment operator for the state class.

Parameters:
right Reference to the class to be copied.

Definition at line 67 of file State.cpp.

References State::m_dens, State::m_kk, State::m_mmw, State::m_molwts, State::m_rmolwts, State::m_stateNum, State::m_temp, State::m_y, and State::m_ym.

Referenced by State::State().

bool ready (  )  const

True if the number species has been set.

Reimplemented in Phase.

Definition at line 279 of file State.cpp.

References State::m_kk.

Referenced by Phase::ready().

void setConcentrations ( const doublereal *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:
conc The input vector to this routine is in dimensional units. For volumetric 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 IdealSolidSolnPhase.

Definition at line 196 of file State.cpp.

References State::m_kk, State::m_mmw, State::m_molwts, State::m_y, State::m_ym, State::setDensity(), and State::stateMFChangeCalc().

Referenced by SurfPhase::setCoverages(), and SurfPhase::setCoveragesNoNorm().

virtual void setDensity ( const doublereal  density  )  [inline, virtual]
void setMassFractions ( const doublereal *const   y  )  [virtual]

Set the mass fractions to the specified values, and then normalize them so that they sum to 1.0.

Parameters:
y Array of unnormalized mass fraction values (input). Must have a length greater than or equal to the number of species.
y Input vector of mass fractions. There is no restriction on the sum of the mass fraction vector. Internally, the State object will normalize this vector before storring its contents. Length is m_kk.

Reimplemented in IdealSolidSolnPhase.

Definition at line 142 of file State.cpp.

References State::m_kk, State::m_mmw, State::m_rmolwts, State::m_y, State::m_ym, Cantera::scale(), and State::stateMFChangeCalc().

Referenced by Phase::setMassFractionsByName(), ThermoPhase::setState_PY(), Phase::setState_RY(), ThermoPhase::setState_TPY(), Phase::setState_TRY(), and Phase::setState_TY().

void setMassFractions_NoNorm ( const doublereal *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:
y Input vector of mass fractions. Length is m_kk.

Reimplemented in IdealSolidSolnPhase.

Definition at line 167 of file State.cpp.

References State::m_kk, State::m_mmw, State::m_rmolwts, State::m_y, State::m_ym, and State::stateMFChangeCalc().

Referenced by Phase::restoreState().

void setMolarDensity ( const doublereal  molarDensity  )  [virtual]

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

Parameters:
molarDensity Input molar density (kmol/m^3).

Reimplemented in DebyeHuckel, HMWSoln, IdealMolalSoln, and IdealSolidSolnPhase.

Definition at line 239 of file State.cpp.

References State::m_dens, and State::meanMolecularWeight().

Referenced by Phase::setState_TNX().

void setMolecularWeight ( const int  k,
const double  mw 
) [inline, protected]

Set the molecular weight of a single species to a given value.

Parameters:
k id of the species
mw Molecular Weight (kg kmol-1)

Definition at line 385 of file State.h.

References State::m_molwts, and State::m_rmolwts.

Referenced by WaterSSTP::initThermoXML().

void setMoleFractions ( const doublereal *const   x  )  [virtual]

Set the mole fractions to the specified values, and then normalize them so that they sum to 1.0.

Parameters:
x Array of unnormalized mole fraction values (input). Must have a length greater than or equal to the number of species.
x Input vector of mole fractions. There is no restriction on the sum of the mole fraction vector. Internally, the State object will normalize this vector before storring its contents. Length is m_kk.

Reimplemented in IdealSolidSolnPhase.

Definition at line 102 of file State.cpp.

References Cantera::dot(), State::m_kk, State::m_mmw, State::m_molwts, State::m_y, State::m_ym, and State::stateMFChangeCalc().

Referenced by SingleSpeciesTP::initThermo(), WaterSSTP::initThermoXML(), MolalityVPSSTP::setMolalities(), MolalityVPSSTP::setMolalitiesByName(), Phase::setMoleFractionsByName(), ThermoPhase::setState_PX(), Phase::setState_RX(), Phase::setState_TNX(), ThermoPhase::setState_TPX(), Phase::setState_TRX(), and Phase::setState_TX().

void setMoleFractions_NoNorm ( const doublereal *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:
x Input vector of mole fractions. Length is m_kk.

Reimplemented in IdealSolidSolnPhase.

Definition at line 115 of file State.cpp.

References Cantera::dot(), State::m_kk, State::m_mmw, State::m_molwts, State::m_y, State::m_ym, and State::stateMFChangeCalc().

virtual void setTemperature ( const doublereal  temp  )  [inline, virtual]
void stateMFChangeCalc ( bool  forceChange = false  )  [inline]

Every time the mole fractions have changed, this routine will increment the stateMFNumber.

Parameters:
forceChange If this is true then the stateMFNumber always changes. This defaults to false.

Definition at line 31 of file State.cpp.

References State::m_stateNum.

Referenced by State::setConcentrations(), State::setMassFractions(), State::setMassFractions_NoNorm(), State::setMoleFractions(), and State::setMoleFractions_NoNorm().

int stateMFNumber (  )  const [inline]

Return the state number.

Return the State Mole Fraction Number.

Definition at line 445 of file State.h.

References State::m_stateNum.

doublereal sum_xlogQ ( doublereal *const  Q  )  const

Evaluate $ \sum_k X_k \log Q_k $.

Parameters:
Q Vector of length m_kk to take the log average of
Returns:
Returns the indicated sum.

Definition at line 188 of file State.cpp.

References State::m_mmw, and State::m_ym.

doublereal sum_xlogx (  )  const

Evaluate $ \sum_k X_k \log X_k $.

Returns:
returns the indicated sum. units are dimensionless.

Definition at line 184 of file State.cpp.

References State::m_mmw, and State::m_ym.

Referenced by IdealSolnGasVPSS::entropy_mole(), IdealSolidSolnPhase::entropy_mole(), IdealGasPhase::entropy_mole(), and IdealSolidSolnPhase::gibbs_mole().

doublereal temperature (  )  const [inline]

Temperature (K).

Definition at line 309 of file State.h.

References State::m_temp.

Referenced by ThermoPhase::_RT(), VPStandardStateTP::_updateStandardStateThermo(), SurfPhase::_updateThermo(), SingleSpeciesTP::_updateThermo(), IdealSolidSolnPhase::_updateThermo(), IdealGasPhase::_updateThermo(), HMWSoln::A_Debye_TP(), DebyeHuckel::A_Debye_TP(), MultiPhase::addPhase(), HMWSoln::ADebye_J(), HMWSoln::ADebye_L(), HMWSoln::ADebye_V(), IdealSolnGasVPSS::calcDensity(), HMWSoln::calcDensity(), LatticePhase::cp_mole(), ConstDensityThermo::cp_mole(), SingleSpeciesTP::cv_mole(), HMWSoln::d2A_DebyedT2_TP(), DebyeHuckel::d2A_DebyedT2_TP(), HMWSoln::dA_DebyedP_TP(), DebyeHuckel::dA_DebyedP_TP(), HMWSoln::dA_DebyedT_TP(), DebyeHuckel::dA_DebyedT_TP(), WaterSSTP::dthermalExpansionCoeffdT(), IdealSolnGasVPSS::enthalpy_mole(), IdealSolidSolnPhase::enthalpy_mole(), IdealGasPhase::enthalpy_mole(), SurfPhase::getChemPotentials(), IdealSolnGasVPSS::getChemPotentials(), IdealSolidSolnPhase::getChemPotentials(), IdealMolalSoln::getChemPotentials(), IdealGasPhase::getChemPotentials(), HMWSoln::getChemPotentials(), DebyeHuckel::getChemPotentials(), SingleSpeciesTP::getChemPotentials_RT(), IdealSolidSolnPhase::getChemPotentials_RT(), WaterSSTP::getCp_R_ref(), ThermoPhase::getElementPotentials(), WaterSSTP::getEnthalpy_RT(), SurfPhase::getEnthalpy_RT(), StoichSubstanceSSTP::getEnthalpy_RT(), MineralEQ3::getEnthalpy_RT(), IdealSolidSolnPhase::getEnthalpy_RT(), WaterSSTP::getEnthalpy_RT_ref(), WaterSSTP::getEntropy_R_ref(), SingleSpeciesTP::getGibbs_ref(), IdealSolidSolnPhase::getGibbs_ref(), WaterSSTP::getGibbs_RT(), SurfPhase::getGibbs_RT(), WaterSSTP::getGibbs_RT_ref(), StoichSubstanceSSTP::getIntEnergy_RT(), MineralEQ3::getIntEnergy_RT(), IdealSolidSolnPhase::getIntEnergy_RT(), StoichSubstanceSSTP::getIntEnergy_RT_ref(), MineralEQ3::getIntEnergy_RT_ref(), MetalSHEelectrons::getIntEnergy_RT_ref(), IdealSolidSolnPhase::getIntEnergy_RT_ref(), HMWSoln::getPartialMolarCp(), DebyeHuckel::getPartialMolarCp(), SurfPhase::getPartialMolarEnthalpies(), SingleSpeciesTP::getPartialMolarEnthalpies(), IdealSolnGasVPSS::getPartialMolarEnthalpies(), IdealSolidSolnPhase::getPartialMolarEnthalpies(), IdealGasPhase::getPartialMolarEnthalpies(), HMWSoln::getPartialMolarEnthalpies(), DebyeHuckel::getPartialMolarEnthalpies(), HMWSoln::getPartialMolarEntropies(), DebyeHuckel::getPartialMolarEntropies(), SingleSpeciesTP::getPartialMolarIntEnergies(), IdealSolnGasVPSS::getPartialMolarIntEnergies(), IdealGasPhase::getPartialMolarIntEnergies(), HMWSoln::getPartialMolarVolumes(), DebyeHuckel::getPartialMolarVolumes(), SingleSpeciesTP::getPureGibbs(), WaterSSTP::getStandardChemPotentials(), StoichSubstanceSSTP::getStandardChemPotentials(), MineralEQ3::getStandardChemPotentials(), MetalSHEelectrons::getStandardChemPotentials(), IdealGasPhase::getStandardChemPotentials(), WaterSSTP::getStandardVolumes_ref(), IdealSolnGasVPSS::gibbs_mole(), IdealSolidSolnPhase::gibbs_mole(), IdealGasPhase::gibbs_mole(), IdealSolidSolnPhase::intEnergy_mole(), IdealGasPhase::intEnergy_mole(), IdealGasPhase::logStandardConc(), IdealGasPhase::pressure(), ThermoPhase::report(), PureFluidPhase::report(), MolalityVPSSTP::report(), HMWSoln::s_updatePitzer_CoeffWRTemp(), HMWSoln::s_updatePitzer_dlnMolalityActCoeff_dP(), HMWSoln::s_updatePitzer_lnMolalityActCoeff(), WaterSSTP::satPressure(), HMWSoln::satPressure(), Phase::saveState(), WaterSSTP::setDensity(), ThermoPhase::setElementPotentials(), WaterSSTP::setPressure(), VPStandardStateTP::setPressure(), IdealMolalSoln::setPressure(), IdealGasPhase::setPressure(), DebyeHuckel::setPressure(), SingleSpeciesTP::setState_HP(), ThermoPhase::setState_HPorUV(), SingleSpeciesTP::setState_SP(), ThermoPhase::setState_SPorSV(), SingleSpeciesTP::setState_SV(), SingleSpeciesTP::setState_UV(), IdealSolnGasVPSS::standardConcentration(), IdealGasPhase::standardConcentration(), MetalSHEelectrons::thermalExpansionCoeff(), IdealGasPhase::thermalExpansionCoeff(), VPStandardStateTP::updateStandardStateThermo(), and WaterSSTP::vaporFraction().


Member Data Documentation

doublereal m_dens [private]

Density.

This is an independent variable except in the incompressible degenerate case. Thus, the pressure is determined from this variable not the other way round. units = kg m-3

Definition at line 405 of file State.h.

Referenced by State::concentration(), State::density(), State::getConcentrations(), State::operator=(), State::setDensity(), and State::setMolarDensity().

int m_kk [protected]
doublereal m_mmw [private]
array_fp m_molwts [private]

m_molwts[k] = molecular weight of species k (kg kmol-1)

Definition at line 427 of file State.h.

Referenced by State::init(), State::molecularWeights(), State::operator=(), State::setConcentrations(), State::setMolecularWeight(), State::setMoleFractions(), and State::setMoleFractions_NoNorm().

array_fp m_rmolwts [private]

m_rmolwts[k] = inverse of the molecular weight of species k units = kmol kg-1.

Definition at line 433 of file State.h.

Referenced by State::concentration(), State::init(), State::operator=(), State::setMassFractions(), State::setMassFractions_NoNorm(), and State::setMolecularWeight().

int m_stateNum [private]

State Change variable.

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

Definition at line 440 of file State.h.

Referenced by State::operator=(), State::stateMFChangeCalc(), and State::stateMFNumber().

doublereal m_temp [private]

Temperature.

This is an independent variable units = Kelvin

Definition at line 396 of file State.h.

Referenced by State::operator=(), State::setTemperature(), and State::temperature().

array_fp m_y [mutable, private]
array_fp m_ym [mutable, private]

The documentation for this class was generated from the following files:
Generated by  doxygen 1.6.3