Cantera
3.0.0
|
Class MaskellSolidSolnPhase represents a condensed phase non-ideal solution with 2 species following the thermodynamic model described in Maskell, Shaw, and Tye, Manganese Dioxide Electrode – IX, Electrochimica Acta 28(2) pp 231-235, 1983. More...
#include <MaskellSolidSolnPhase.h>
Class MaskellSolidSolnPhase represents a condensed phase non-ideal solution with 2 species following the thermodynamic model described in Maskell, Shaw, and Tye, Manganese Dioxide Electrode – IX, Electrochimica Acta 28(2) pp 231-235, 1983.
Definition at line 30 of file MaskellSolidSolnPhase.h.
Public Member Functions | |
string | type () const override |
String indicating the thermodynamic model implemented. | |
Units | standardConcentrationUnits () const override |
Returns the units of the "standard concentration" for this phase. | |
void | getActivityConcentrations (double *c) const override |
This method returns an array of generalized concentrations. | |
double | standardConcentration (size_t k=0) const override |
Return the standard concentration for the kth species. | |
double | logStandardConc (size_t k=0) const override |
Natural logarithm of the standard concentration of the kth species. | |
Molar Thermodynamic Properties of the Solution | |
double | enthalpy_mole () const override |
Molar enthalpy. Units: J/kmol. | |
double | entropy_mole () const override |
Molar entropy. Units: J/kmol/K. | |
Mechanical Equation of State Properties | |
In this equation of state implementation, the density is a function only of the mole fractions. Therefore, it can't be an independent variable. Instead, the pressure is used as the independent variable. Functions which try to set the thermodynamic state by calling setDensity() will cause an exception to be thrown. | |
double | pressure () const override |
Pressure. | |
void | setPressure (double p) override |
Set the pressure at constant temperature. | |
void | calcDensity () override |
Calculate the density of the mixture using the partial molar volumes and mole fractions as input. | |
Chemical Potentials and Activities | |
void | getActivityCoefficients (double *ac) const override |
Get the array of non-dimensional molar-based activity coefficients at the current solution temperature, pressure, and solution concentration. | |
void | getChemPotentials (double *mu) const override |
Get the species chemical potentials. Units: J/kmol. | |
void | getChemPotentials_RT (double *mu) const override |
Partial Molar Properties of the Solution | |
void | getPartialMolarEnthalpies (double *hbar) const override |
Returns an array of partial molar enthalpies for the species in the mixture. | |
void | getPartialMolarEntropies (double *sbar) const override |
Returns an array of partial molar entropies of the species in the solution. | |
void | getPartialMolarCp (double *cpbar) const override |
Return an array of partial molar heat capacities for the species in the mixture. | |
void | getPartialMolarVolumes (double *vbar) const override |
Return an array of partial molar volumes for the species in the mixture. | |
void | getPureGibbs (double *gpure) const override |
Get the Gibbs functions for the standard state of the species at the current T and P of the solution. | |
void | getStandardChemPotentials (double *mu) const override |
Get the array of chemical potentials at unit activity for the species at their standard states at the current T and P of the solution. | |
Utility Functions | |
void | initThermo () override |
Initialize the ThermoPhase object after all species have been set up. | |
void | getParameters (AnyMap &phaseNode) const override |
Store the parameters of a ThermoPhase object such that an identical one could be reconstructed using the newThermo(AnyMap&) function. | |
void | set_h_mix (const double hmix) |
void | setProductSpecies (const string &name) |
Set the product Species. Must be called after species have been added. | |
Public Member Functions inherited from VPStandardStateTP | |
void | setTemperature (const double temp) override |
Set the temperature of the phase. | |
void | setPressure (double p) override |
Set the internally stored pressure (Pa) at constant temperature and composition. | |
void | setState_TP (double T, double pres) override |
Set the temperature and pressure at the same time. | |
double | pressure () const override |
Returns the current pressure of the phase. | |
virtual void | updateStandardStateThermo () const |
Updates the standard state thermodynamic functions at the current T and P of the solution. | |
double | minTemp (size_t k=npos) const override |
Minimum temperature for which the thermodynamic data for the species or phase are valid. | |
double | maxTemp (size_t k=npos) const override |
Maximum temperature for which the thermodynamic data for the species are valid. | |
PDSS * | providePDSS (size_t k) |
const PDSS * | providePDSS (size_t k) const |
VPStandardStateTP () | |
Constructor. | |
bool | isCompressible () const override |
Return whether phase represents a compressible substance. | |
int | standardStateConvention () const override |
This method returns the convention used in specification of the standard state, of which there are currently two, temperature based, and variable pressure based. | |
void | getChemPotentials_RT (double *mu) const override |
Get the array of non-dimensional species chemical potentials. | |
void | getStandardChemPotentials (double *mu) const override |
Get the array of chemical potentials at unit activity for the species at their standard states at the current T and P of the solution. | |
void | getEnthalpy_RT (double *hrt) const override |
Get the nondimensional Enthalpy functions for the species at their standard states at the current T and P of the solution. | |
void | getEntropy_R (double *sr) const override |
Get the array of nondimensional Entropy functions for the standard state species at the current T and P of the solution. | |
void | getGibbs_RT (double *grt) const override |
Get the nondimensional Gibbs functions for the species in their standard states at the current T and P of the solution. | |
void | getPureGibbs (double *gpure) const override |
Get the Gibbs functions for the standard state of the species at the current T and P of the solution. | |
void | getIntEnergy_RT (double *urt) const override |
Returns the vector of nondimensional Internal Energies of the standard state species at the current T and P of the solution. | |
void | getCp_R (double *cpr) const override |
Get the nondimensional Heat Capacities at constant pressure for the species standard states at the current T and P of the solution. | |
void | getStandardVolumes (double *vol) const override |
Get the molar volumes of the species standard states at the current T and P of the solution. | |
virtual const vector< double > & | getStandardVolumes () const |
void | initThermo () override |
Initialize the ThermoPhase object after all species have been set up. | |
void | getSpeciesParameters (const string &name, AnyMap &speciesNode) const override |
Get phase-specific parameters of a Species object such that an identical one could be reconstructed and added to this phase. | |
bool | addSpecies (shared_ptr< Species > spec) override |
Add a Species to this Phase. | |
void | installPDSS (size_t k, unique_ptr< PDSS > &&pdss) |
Install a PDSS object for species k | |
virtual bool | addSpecies (shared_ptr< Species > spec) |
Add a Species to this Phase. | |
void | getEnthalpy_RT_ref (double *hrt) const override |
Returns the vector of nondimensional enthalpies of the reference state at the current temperature of the solution and the reference pressure for the species. | |
void | getGibbs_RT_ref (double *grt) const override |
Returns the vector of nondimensional Gibbs Free Energies of the reference state at the current temperature of the solution and the reference pressure for the species. | |
void | getGibbs_ref (double *g) const override |
Returns the vector of the Gibbs function of the reference state at the current temperature of the solution and the reference pressure for the species. | |
void | getEntropy_R_ref (double *er) const override |
Returns the vector of nondimensional entropies of the reference state at the current temperature of the solution and the reference pressure for each species. | |
void | getCp_R_ref (double *cprt) const override |
Returns the vector of nondimensional constant pressure heat capacities of the reference state at the current temperature of the solution and reference pressure for each species. | |
void | getStandardVolumes_ref (double *vol) const override |
Get the molar volumes of the species reference states at the current T and P_ref of the solution. | |
Public Member Functions inherited from ThermoPhase | |
ThermoPhase ()=default | |
Constructor. | |
double | RT () const |
Return the Gas Constant multiplied by the current temperature. | |
double | equivalenceRatio () const |
Compute the equivalence ratio for the current mixture from available oxygen and required oxygen. | |
string | type () const override |
String indicating the thermodynamic model implemented. | |
virtual bool | isIdeal () const |
Boolean indicating whether phase is ideal. | |
virtual string | phaseOfMatter () const |
String indicating the mechanical phase of the matter in this Phase. | |
virtual double | refPressure () const |
Returns the reference pressure in Pa. | |
double | Hf298SS (const size_t k) const |
Report the 298 K Heat of Formation of the standard state of one species (J kmol-1) | |
virtual void | modifyOneHf298SS (const size_t k, const double Hf298New) |
Modify the value of the 298 K Heat of Formation of one species in the phase (J kmol-1) | |
virtual void | resetHf298 (const size_t k=npos) |
Restore the original heat of formation of one or more species. | |
bool | chargeNeutralityNecessary () const |
Returns the chargeNeutralityNecessity boolean. | |
virtual double | intEnergy_mole () const |
Molar internal energy. Units: J/kmol. | |
virtual double | gibbs_mole () const |
Molar Gibbs function. Units: J/kmol. | |
virtual double | cp_mole () const |
Molar heat capacity at constant pressure. Units: J/kmol/K. | |
virtual double | cv_mole () const |
Molar heat capacity at constant volume. Units: J/kmol/K. | |
virtual double | isothermalCompressibility () const |
Returns the isothermal compressibility. Units: 1/Pa. | |
virtual double | thermalExpansionCoeff () const |
Return the volumetric thermal expansion coefficient. Units: 1/K. | |
virtual double | soundSpeed () const |
Return the speed of sound. Units: m/s. | |
void | setElectricPotential (double v) |
Set the electric potential of this phase (V). | |
double | electricPotential () const |
Returns the electric potential of this phase (V). | |
virtual int | activityConvention () const |
This method returns the convention used in specification of the activities, of which there are currently two, molar- and molality-based conventions. | |
virtual void | getActivities (double *a) const |
Get the array of non-dimensional activities at the current solution temperature, pressure, and solution concentration. | |
virtual void | getLnActivityCoefficients (double *lnac) const |
Get the array of non-dimensional molar-based ln activity coefficients at the current solution temperature, pressure, and solution concentration. | |
void | getElectrochemPotentials (double *mu) const |
Get the species electrochemical potentials. | |
virtual void | getPartialMolarIntEnergies (double *ubar) const |
Return an array of partial molar internal energies for the species in the mixture. | |
virtual void | getIntEnergy_RT_ref (double *urt) const |
Returns the vector of nondimensional internal Energies of the reference state at the current temperature of the solution and the reference pressure for each species. | |
double | enthalpy_mass () const |
Specific enthalpy. Units: J/kg. | |
double | intEnergy_mass () const |
Specific internal energy. Units: J/kg. | |
double | entropy_mass () const |
Specific entropy. Units: J/kg/K. | |
double | gibbs_mass () const |
Specific Gibbs function. Units: J/kg. | |
double | cp_mass () const |
Specific heat at constant pressure. Units: J/kg/K. | |
double | cv_mass () const |
Specific heat at constant volume. Units: J/kg/K. | |
virtual void | setState_TPX (double t, double p, const double *x) |
Set the temperature (K), pressure (Pa), and mole fractions. | |
virtual void | setState_TPX (double t, double p, const Composition &x) |
Set the temperature (K), pressure (Pa), and mole fractions. | |
virtual void | setState_TPX (double t, double p, const string &x) |
Set the temperature (K), pressure (Pa), and mole fractions. | |
virtual void | setState_TPY (double t, double p, const double *y) |
Set the internally stored temperature (K), pressure (Pa), and mass fractions of the phase. | |
virtual void | setState_TPY (double t, double p, const Composition &y) |
Set the internally stored temperature (K), pressure (Pa), and mass fractions of the phase. | |
virtual void | setState_TPY (double t, double p, const string &y) |
Set the internally stored temperature (K), pressure (Pa), and mass fractions of the phase. | |
virtual void | setState_PX (double p, double *x) |
Set the pressure (Pa) and mole fractions. | |
virtual void | setState_PY (double p, double *y) |
Set the internally stored pressure (Pa) and mass fractions. | |
virtual void | setState_HP (double h, double p, double tol=1e-9) |
Set the internally stored specific enthalpy (J/kg) and pressure (Pa) of the phase. | |
virtual void | setState_UV (double u, double v, double tol=1e-9) |
Set the specific internal energy (J/kg) and specific volume (m^3/kg). | |
virtual void | setState_SP (double s, double p, double tol=1e-9) |
Set the specific entropy (J/kg/K) and pressure (Pa). | |
virtual void | setState_SV (double s, double v, double tol=1e-9) |
Set the specific entropy (J/kg/K) and specific volume (m^3/kg). | |
virtual void | setState_ST (double s, double t, double tol=1e-9) |
Set the specific entropy (J/kg/K) and temperature (K). | |
virtual void | setState_TV (double t, double v, double tol=1e-9) |
Set the temperature (K) and specific volume (m^3/kg). | |
virtual void | setState_PV (double p, double v, double tol=1e-9) |
Set the pressure (Pa) and specific volume (m^3/kg). | |
virtual void | setState_UP (double u, double p, double tol=1e-9) |
Set the specific internal energy (J/kg) and pressure (Pa). | |
virtual void | setState_VH (double v, double h, double tol=1e-9) |
Set the specific volume (m^3/kg) and the specific enthalpy (J/kg) | |
virtual void | setState_TH (double t, double h, double tol=1e-9) |
Set the temperature (K) and the specific enthalpy (J/kg) | |
virtual void | setState_SH (double s, double h, double tol=1e-9) |
Set the specific entropy (J/kg/K) and the specific enthalpy (J/kg) | |
void | setState_RP (double rho, double p) |
Set the density (kg/m**3) and pressure (Pa) at constant composition. | |
virtual void | setState_DP (double rho, double p) |
Set the density (kg/m**3) and pressure (Pa) at constant composition. | |
virtual void | setState_RPX (double rho, double p, const double *x) |
Set the density (kg/m**3), pressure (Pa) and mole fractions. | |
virtual void | setState_RPX (double rho, double p, const Composition &x) |
Set the density (kg/m**3), pressure (Pa) and mole fractions. | |
virtual void | setState_RPX (double rho, double p, const string &x) |
Set the density (kg/m**3), pressure (Pa) and mole fractions. | |
virtual void | setState_RPY (double rho, double p, const double *y) |
Set the density (kg/m**3), pressure (Pa) and mass fractions. | |
virtual void | setState_RPY (double rho, double p, const Composition &y) |
Set the density (kg/m**3), pressure (Pa) and mass fractions. | |
virtual void | setState_RPY (double rho, double p, const string &y) |
Set the density (kg/m**3), pressure (Pa) and mass fractions. | |
virtual void | setState (const AnyMap &state) |
Set the state using an AnyMap containing any combination of properties supported by the thermodynamic model. | |
void | setMixtureFraction (double mixFrac, const double *fuelComp, const double *oxComp, ThermoBasis basis=ThermoBasis::molar) |
Set the mixture composition according to the mixture fraction = kg fuel / (kg oxidizer + kg fuel) | |
void | setMixtureFraction (double mixFrac, const string &fuelComp, const string &oxComp, ThermoBasis basis=ThermoBasis::molar) |
Set the mixture composition according to the mixture fraction = kg fuel / (kg oxidizer + kg fuel) | |
void | setMixtureFraction (double mixFrac, const Composition &fuelComp, const Composition &oxComp, ThermoBasis basis=ThermoBasis::molar) |
Set the mixture composition according to the mixture fraction = kg fuel / (kg oxidizer + kg fuel) | |
double | mixtureFraction (const double *fuelComp, const double *oxComp, ThermoBasis basis=ThermoBasis::molar, const string &element="Bilger") const |
Compute the mixture fraction = kg fuel / (kg oxidizer + kg fuel) for the current mixture given fuel and oxidizer compositions. | |
double | mixtureFraction (const string &fuelComp, const string &oxComp, ThermoBasis basis=ThermoBasis::molar, const string &element="Bilger") const |
Compute the mixture fraction = kg fuel / (kg oxidizer + kg fuel) for the current mixture given fuel and oxidizer compositions. | |
double | mixtureFraction (const Composition &fuelComp, const Composition &oxComp, ThermoBasis basis=ThermoBasis::molar, const string &element="Bilger") const |
Compute the mixture fraction = kg fuel / (kg oxidizer + kg fuel) for the current mixture given fuel and oxidizer compositions. | |
void | setEquivalenceRatio (double phi, const double *fuelComp, const double *oxComp, ThermoBasis basis=ThermoBasis::molar) |
Set the mixture composition according to the equivalence ratio. | |
void | setEquivalenceRatio (double phi, const string &fuelComp, const string &oxComp, ThermoBasis basis=ThermoBasis::molar) |
Set the mixture composition according to the equivalence ratio. | |
void | setEquivalenceRatio (double phi, const Composition &fuelComp, const Composition &oxComp, ThermoBasis basis=ThermoBasis::molar) |
Set the mixture composition according to the equivalence ratio. | |
double | equivalenceRatio (const double *fuelComp, const double *oxComp, ThermoBasis basis=ThermoBasis::molar) const |
Compute the equivalence ratio for the current mixture given the compositions of fuel and oxidizer. | |
double | equivalenceRatio (const string &fuelComp, const string &oxComp, ThermoBasis basis=ThermoBasis::molar) const |
Compute the equivalence ratio for the current mixture given the compositions of fuel and oxidizer. | |
double | equivalenceRatio (const Composition &fuelComp, const Composition &oxComp, ThermoBasis basis=ThermoBasis::molar) const |
Compute the equivalence ratio for the current mixture given the compositions of fuel and oxidizer. | |
double | stoichAirFuelRatio (const double *fuelComp, const double *oxComp, ThermoBasis basis=ThermoBasis::molar) const |
Compute the stoichiometric air to fuel ratio (kg oxidizer / kg fuel) given fuel and oxidizer compositions. | |
double | stoichAirFuelRatio (const string &fuelComp, const string &oxComp, ThermoBasis basis=ThermoBasis::molar) const |
Compute the stoichiometric air to fuel ratio (kg oxidizer / kg fuel) given fuel and oxidizer compositions. | |
double | stoichAirFuelRatio (const Composition &fuelComp, const Composition &oxComp, ThermoBasis basis=ThermoBasis::molar) const |
Compute the stoichiometric air to fuel ratio (kg oxidizer / kg fuel) given fuel and oxidizer compositions. | |
void | equilibrate (const string &XY, const string &solver="auto", double rtol=1e-9, int max_steps=50000, int max_iter=100, int estimate_equil=0, int log_level=0) |
Equilibrate a ThermoPhase object. | |
virtual void | setToEquilState (const double *mu_RT) |
This method is used by the ChemEquil equilibrium solver. | |
virtual bool | compatibleWithMultiPhase () const |
Indicates whether this phase type can be used with class MultiPhase for equilibrium calculations. | |
virtual double | critTemperature () const |
Critical temperature (K). | |
virtual double | critPressure () const |
Critical pressure (Pa). | |
virtual double | critVolume () const |
Critical volume (m3/kmol). | |
virtual double | critCompressibility () const |
Critical compressibility (unitless). | |
virtual double | critDensity () const |
Critical density (kg/m3). | |
virtual double | satTemperature (double p) const |
Return the saturation temperature given the pressure. | |
virtual double | satPressure (double t) |
Return the saturation pressure given the temperature. | |
virtual double | vaporFraction () const |
Return the fraction of vapor at the current conditions. | |
virtual void | setState_Tsat (double t, double x) |
Set the state to a saturated system at a particular temperature. | |
virtual void | setState_Psat (double p, double x) |
Set the state to a saturated system at a particular pressure. | |
void | setState_TPQ (double T, double P, double Q) |
Set the temperature, pressure, and vapor fraction (quality). | |
void | modifySpecies (size_t k, shared_ptr< Species > spec) override |
Modify the thermodynamic data associated with a species. | |
virtual MultiSpeciesThermo & | speciesThermo (int k=-1) |
Return a changeable reference to the calculation manager for species reference-state thermodynamic properties. | |
virtual const MultiSpeciesThermo & | speciesThermo (int k=-1) const |
void | initThermoFile (const string &inputFile, const string &id) |
Initialize a ThermoPhase object using an input file. | |
virtual void | setParameters (const AnyMap &phaseNode, const AnyMap &rootNode=AnyMap()) |
Set equation of state parameters from an AnyMap phase description. | |
AnyMap | parameters (bool withInput=true) const |
Returns the parameters of a ThermoPhase object such that an identical one could be reconstructed using the newThermo(AnyMap&) function. | |
const AnyMap & | input () const |
Access input data associated with the phase description. | |
AnyMap & | input () |
virtual void | getdlnActCoeffds (const double dTds, const double *const dXds, double *dlnActCoeffds) const |
Get the change in activity coefficients wrt changes in state (temp, mole fraction, etc) along a line in parameter space or along a line in physical space. | |
virtual void | getdlnActCoeffdlnX_diag (double *dlnActCoeffdlnX_diag) const |
Get the array of ln mole fraction derivatives of the log activity coefficients - diagonal component only. | |
virtual void | getdlnActCoeffdlnN_diag (double *dlnActCoeffdlnN_diag) const |
Get the array of log species mole number derivatives of the log activity coefficients. | |
virtual void | getdlnActCoeffdlnN (const size_t ld, double *const dlnActCoeffdlnN) |
Get the array of derivatives of the log activity coefficients with respect to the log of the species mole numbers. | |
virtual void | getdlnActCoeffdlnN_numderiv (const size_t ld, double *const dlnActCoeffdlnN) |
virtual string | report (bool show_thermo=true, double threshold=-1e-14) const |
returns a summary of the state of the phase as a string | |
virtual void | reportCSV (std::ofstream &csvFile) const |
returns a summary of the state of the phase to a comma separated file. | |
Public Member Functions inherited from Phase | |
Phase ()=default | |
Default constructor. | |
Phase (const Phase &)=delete | |
Phase & | operator= (const Phase &)=delete |
virtual bool | isPure () const |
Return whether phase represents a pure (single species) substance. | |
virtual bool | hasPhaseTransition () const |
Return whether phase represents a substance with phase transitions. | |
virtual bool | isCompressible () const |
Return whether phase represents a compressible substance. | |
virtual map< string, size_t > | nativeState () const |
Return a map of properties defining the native state of a substance. | |
string | nativeMode () const |
Return string acronym representing the native state of a Phase. | |
virtual vector< string > | fullStates () const |
Return a vector containing full states defining a phase. | |
virtual vector< string > | partialStates () const |
Return a vector of settable partial property sets within a phase. | |
virtual size_t | stateSize () const |
Return size of vector defining internal state of the phase. | |
void | saveState (vector< double > &state) const |
Save the current internal state of the phase. | |
virtual void | saveState (size_t lenstate, double *state) const |
Write to array 'state' the current internal state. | |
void | restoreState (const vector< double > &state) |
Restore a state saved on a previous call to saveState. | |
virtual void | restoreState (size_t lenstate, const double *state) |
Restore the state of the phase from a previously saved state vector. | |
double | molecularWeight (size_t k) const |
Molecular weight of species k . | |
void | getMolecularWeights (vector< double > &weights) const |
Copy the vector of molecular weights into vector weights. | |
void | getMolecularWeights (double *weights) const |
Copy the vector of molecular weights into array weights. | |
const vector< double > & | molecularWeights () const |
Return a const reference to the internal vector of molecular weights. | |
const vector< double > & | inverseMolecularWeights () const |
Return a const reference to the internal vector of molecular weights. | |
void | getCharges (double *charges) const |
Copy the vector of species charges into array charges. | |
virtual void | setMolesNoTruncate (const double *const N) |
Set the state of the object with moles in [kmol]. | |
double | elementalMassFraction (const size_t m) const |
Elemental mass fraction of element m. | |
double | elementalMoleFraction (const size_t m) const |
Elemental mole fraction of element m. | |
const double * | moleFractdivMMW () const |
Returns a const pointer to the start of the moleFraction/MW array. | |
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. | |
double | chargeDensity () const |
Charge density [C/m^3]. | |
size_t | nDim () const |
Returns the number of spatial dimensions (1, 2, or 3) | |
void | setNDim (size_t ndim) |
Set the number of spatial dimensions (1, 2, or 3). | |
virtual bool | ready () const |
Returns a bool indicating whether the object is ready for use. | |
int | stateMFNumber () const |
Return the State Mole Fraction Number. | |
virtual void | invalidateCache () |
Invalidate any cached values which are normally updated only when a change in state is detected. | |
bool | caseSensitiveSpecies () const |
Returns true if case sensitive species names are enforced. | |
void | setCaseSensitiveSpecies (bool cflag=true) |
Set flag that determines whether case sensitive species are enforced in look-up operations, for example speciesIndex. | |
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. | |
void | massFractionsToMoleFractions (const double *Y, double *X) const |
Converts a mixture composition from mole fractions to mass fractions. | |
void | moleFractionsToMassFractions (const double *X, double *Y) const |
Converts a mixture composition from mass fractions to mole fractions. | |
string | name () const |
Return the name of the phase. | |
void | setName (const string &nm) |
Sets the string name for the phase. | |
string | elementName (size_t m) const |
Name of the element with index m. | |
size_t | elementIndex (const string &name) const |
Return the index of element named 'name'. | |
const vector< string > & | elementNames () const |
Return a read-only reference to the vector of element names. | |
double | atomicWeight (size_t m) const |
Atomic weight of element m. | |
double | entropyElement298 (size_t m) const |
Entropy of the element in its standard state at 298 K and 1 bar. | |
int | atomicNumber (size_t m) const |
Atomic number of element m. | |
int | elementType (size_t m) const |
Return the element constraint type Possible types include: | |
int | changeElementType (int m, int elem_type) |
Change the element type of the mth constraint Reassigns an element type. | |
const vector< double > & | atomicWeights () const |
Return a read-only reference to the vector of atomic weights. | |
size_t | nElements () const |
Number of elements. | |
void | checkElementIndex (size_t m) const |
Check that the specified element index is in range. | |
void | checkElementArraySize (size_t mm) const |
Check that an array size is at least nElements(). | |
double | nAtoms (size_t k, size_t m) const |
Number of atoms of element m in species k . | |
void | getAtoms (size_t k, double *atomArray) const |
Get a vector containing the atomic composition of species k. | |
size_t | speciesIndex (const string &name) const |
Returns the index of a species named 'name' within the Phase object. | |
string | speciesName (size_t k) const |
Name of the species with index k. | |
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. | |
const vector< string > & | speciesNames () const |
Return a const reference to the vector of species names. | |
size_t | nSpecies () const |
Returns the number of species in the phase. | |
void | checkSpeciesIndex (size_t k) const |
Check that the specified species index is in range. | |
void | checkSpeciesArraySize (size_t kk) const |
Check that an array size is at least nSpecies(). | |
void | setMoleFractionsByName (const Composition &xMap) |
Set the species mole fractions by name. | |
void | setMoleFractionsByName (const string &x) |
Set the mole fractions of a group of species by name. | |
void | setMassFractionsByName (const Composition &yMap) |
Set the species mass fractions by name. | |
void | setMassFractionsByName (const string &x) |
Set the species mass fractions by name. | |
void | setState_TRX (double t, double dens, const double *x) |
Set the internally stored temperature (K), density, and mole fractions. | |
void | setState_TRX (double t, double dens, const Composition &x) |
Set the internally stored temperature (K), density, and mole fractions. | |
void | setState_TRY (double t, double dens, const double *y) |
Set the internally stored temperature (K), density, and mass fractions. | |
void | setState_TRY (double t, double dens, const Composition &y) |
Set the internally stored temperature (K), density, and mass fractions. | |
void | setState_TNX (double t, double n, const double *x) |
Set the internally stored temperature (K), molar density (kmol/m^3), and mole fractions. | |
void | setState_TR (double t, double rho) |
Set the internally stored temperature (K) and density (kg/m^3) | |
void | setState_TD (double t, double rho) |
Set the internally stored temperature (K) and density (kg/m^3) | |
void | setState_TX (double t, double *x) |
Set the internally stored temperature (K) and mole fractions. | |
void | setState_TY (double t, double *y) |
Set the internally stored temperature (K) and mass fractions. | |
void | setState_RX (double rho, double *x) |
Set the density (kg/m^3) and mole fractions. | |
void | setState_RY (double rho, double *y) |
Set the density (kg/m^3) and mass fractions. | |
Composition | getMoleFractionsByName (double threshold=0.0) const |
Get the mole fractions by name. | |
double | moleFraction (size_t k) const |
Return the mole fraction of a single species. | |
double | moleFraction (const string &name) const |
Return the mole fraction of a single species. | |
Composition | getMassFractionsByName (double threshold=0.0) const |
Get the mass fractions by name. | |
double | massFraction (size_t k) const |
Return the mass fraction of a single species. | |
double | massFraction (const string &name) const |
Return the mass fraction of a single species. | |
void | getMoleFractions (double *const x) const |
Get the species mole fraction vector. | |
virtual void | setMoleFractions (const double *const x) |
Set the mole fractions to the specified values. | |
virtual void | setMoleFractions_NoNorm (const double *const x) |
Set the mole fractions to the specified values without normalizing. | |
void | getMassFractions (double *const y) const |
Get the species mass fractions. | |
const double * | massFractions () const |
Return a const pointer to the mass fraction array. | |
virtual void | setMassFractions (const double *const y) |
Set the mass fractions to the specified values and normalize them. | |
virtual void | setMassFractions_NoNorm (const double *const y) |
Set the mass fractions to the specified values without normalizing. | |
virtual void | getConcentrations (double *const c) const |
Get the species concentrations (kmol/m^3). | |
virtual double | concentration (const size_t k) const |
Concentration of species k. | |
virtual void | setConcentrations (const double *const conc) |
Set the concentrations to the specified values within the phase. | |
virtual void | setConcentrationsNoNorm (const double *const conc) |
Set the concentrations without ignoring negative concentrations. | |
double | temperature () const |
Temperature (K). | |
virtual double | electronTemperature () const |
Electron Temperature (K) | |
virtual double | density () const |
Density (kg/m^3). | |
virtual double | molarDensity () const |
Molar density (kmol/m^3). | |
virtual double | molarVolume () const |
Molar volume (m^3/kmol). | |
virtual void | setDensity (const double density_) |
Set the internally stored density (kg/m^3) of the phase. | |
virtual void | setMolarDensity (const double molarDensity) |
Set the internally stored molar density (kmol/m^3) of the phase. | |
virtual void | setElectronTemperature (double etemp) |
Set the internally stored electron temperature of the phase (K). | |
double | mean_X (const double *const Q) const |
Evaluate the mole-fraction-weighted mean of an array Q. | |
double | mean_X (const vector< double > &Q) const |
Evaluate the mole-fraction-weighted mean of an array Q. | |
double | meanMolecularWeight () const |
The mean molecular weight. Units: (kg/kmol) | |
double | sum_xlogx () const |
Evaluate \( \sum_k X_k \ln X_k \). | |
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. | |
void | addSpeciesAlias (const string &name, const string &alias) |
Add a species alias (that is, a user-defined alternative species name). | |
virtual vector< string > | findIsomers (const Composition &compMap) const |
Return a vector with isomers names matching a given composition map. | |
virtual vector< string > | findIsomers (const string &comp) const |
Return a vector with isomers names matching a given composition string. | |
shared_ptr< Species > | species (const string &name) const |
Return the Species object for the named species. | |
shared_ptr< Species > | species (size_t k) const |
Return the Species object for species whose index is k. | |
void | ignoreUndefinedElements () |
Set behavior when adding a species containing undefined elements to just skip the species. | |
void | addUndefinedElements () |
Set behavior when adding a species containing undefined elements to add those elements to the phase. | |
void | throwUndefinedElements () |
Set the behavior when adding a species containing undefined elements to throw an exception. | |
Private Member Functions | |
double | s () const |
double | fm (const double r) const |
double | p (const double r) const |
Private Attributes | |
double | m_Pcurrent = OneAtm |
m_Pcurrent = The current pressure. | |
double | h_mixing = 0.0 |
Value of the enthalpy change on mixing due to protons changing from type B to type A configurations. | |
int | product_species_index = -1 |
Index of the species whose mole fraction defines the extent of reduction r. | |
int | reactant_species_index = -1 |
Additional Inherited Members | |
Protected Member Functions inherited from VPStandardStateTP | |
virtual void | calcDensity () |
Calculate the density of the mixture using the partial molar volumes and mole fractions as input. | |
virtual void | _updateStandardStateThermo () const |
Updates the standard state thermodynamic functions at the current T and P of the solution. | |
void | invalidateCache () override |
Invalidate any cached values which are normally updated only when a change in state is detected. | |
const vector< double > & | Gibbs_RT_ref () const |
Protected Member Functions inherited from ThermoPhase | |
virtual void | getParameters (AnyMap &phaseNode) const |
Store the parameters of a ThermoPhase object such that an identical one could be reconstructed using the newThermo(AnyMap&) function. | |
virtual void | getCsvReportData (vector< string > &names, vector< vector< double > > &data) const |
Fills names and data with the column names and species thermo properties to be included in the output of the reportCSV method. | |
Protected Member Functions inherited from Phase | |
void | assertCompressible (const string &setter) const |
Ensure that phase is compressible. | |
void | assignDensity (const double density_) |
Set the internally stored constant density (kg/m^3) of the phase. | |
void | setMolecularWeight (const int k, const double mw) |
Set the molecular weight of a single species to a given value. | |
virtual void | compositionChanged () |
Apply changes to the state which are needed after the composition changes. | |
Protected Attributes inherited from VPStandardStateTP | |
double | m_Pcurrent = OneAtm |
Current value of the pressure - state variable. | |
double | m_minTemp = 0.0 |
The minimum temperature at which data for all species is valid. | |
double | m_maxTemp = BigNumber |
The maximum temperature at which data for all species is valid. | |
double | m_Tlast_ss = -1.0 |
The last temperature at which the standard state thermodynamic properties were calculated at. | |
double | m_Plast_ss = -1.0 |
The last pressure at which the Standard State thermodynamic properties were calculated at. | |
vector< unique_ptr< PDSS > > | m_PDSS_storage |
Storage for the PDSS objects for the species. | |
vector< double > | m_h0_RT |
Vector containing the species reference enthalpies at T = m_tlast and P = p_ref. | |
vector< double > | m_cp0_R |
Vector containing the species reference constant pressure heat capacities at T = m_tlast and P = p_ref. | |
vector< double > | m_g0_RT |
Vector containing the species reference Gibbs functions at T = m_tlast and P = p_ref. | |
vector< double > | m_s0_R |
Vector containing the species reference entropies at T = m_tlast and P = p_ref. | |
vector< double > | m_V0 |
Vector containing the species reference molar volumes. | |
vector< double > | m_hss_RT |
Vector containing the species Standard State enthalpies at T = m_tlast and P = m_plast. | |
vector< double > | m_cpss_R |
Vector containing the species Standard State constant pressure heat capacities at T = m_tlast and P = m_plast. | |
vector< double > | m_gss_RT |
Vector containing the species Standard State Gibbs functions at T = m_tlast and P = m_plast. | |
vector< double > | m_sss_R |
Vector containing the species Standard State entropies at T = m_tlast and P = m_plast. | |
vector< double > | m_Vss |
Vector containing the species standard state volumes at T = m_tlast and P = m_plast. | |
Protected Attributes inherited from ThermoPhase | |
MultiSpeciesThermo | m_spthermo |
Pointer to the calculation manager for species reference-state thermodynamic properties. | |
AnyMap | m_input |
Data supplied via setParameters. | |
double | m_phi = 0.0 |
Stored value of the electric potential for this phase. Units are Volts. | |
bool | m_chargeNeutralityNecessary = false |
Boolean indicating whether a charge neutrality condition is a necessity. | |
int | m_ssConvention = cSS_CONVENTION_TEMPERATURE |
Contains the standard state convention. | |
double | m_tlast = 0.0 |
last value of the temperature processed by reference state | |
Protected Attributes inherited from Phase | |
ValueCache | m_cache |
Cached for saved calculations within each ThermoPhase. | |
size_t | m_kk = 0 |
Number of species in the phase. | |
size_t | m_ndim = 3 |
Dimensionality of the phase. | |
vector< double > | m_speciesComp |
Atomic composition of the species. | |
vector< double > | m_speciesCharge |
Vector of species charges. length m_kk. | |
map< string, shared_ptr< Species > > | m_species |
UndefElement::behavior | m_undefinedElementBehavior = UndefElement::add |
Flag determining behavior when adding species with an undefined element. | |
bool | m_caseSensitiveSpecies = false |
Flag determining whether case sensitive species names are enforced. | |
Definition at line 20 of file MaskellSolidSolnPhase.cpp.
|
inlineoverridevirtual |
String indicating the thermodynamic model implemented.
Usually corresponds to the name of the derived class, less any suffixes such as "Phase", TP", "VPSS", etc.
Reimplemented from Phase.
Definition at line 35 of file MaskellSolidSolnPhase.h.
|
inlineoverridevirtual |
Returns the units of the "standard concentration" for this phase.
These are the units of the values returned by the functions getActivityConcentrations() and standardConcentration(), which can vary between different ThermoPhase-derived classes, or change within a single class depending on input options. See the documentation for standardConcentration() for the derived class for specific details.
Reimplemented from ThermoPhase.
Definition at line 39 of file MaskellSolidSolnPhase.h.
|
overridevirtual |
This method returns an array of generalized concentrations.
\( C^a_k \) are defined such that \( a_k = C^a_k / C^0_k, \) where \( C^0_k \) is a standard concentration defined below and \( a_k \) are activities used in the thermodynamic functions. These activity (or generalized) concentrations are used by kinetics manager classes to compute the forward and reverse rates of elementary reactions. Note that they may or may not have units of concentration — they might be partial pressures, mole fractions, or surface coverages, for example.
c | Output array of generalized concentrations. The units depend upon the implementation of the reaction rate expressions within the phase. |
Reimplemented from ThermoPhase.
Definition at line 25 of file MaskellSolidSolnPhase.cpp.
|
inlineoverridevirtual |
Return the standard concentration for the kth species.
The standard concentration \( C^0_k \) used to normalize the activity (that is, generalized) concentration. In many cases, this quantity will be the same for all species in a phase - for example, for an ideal gas \( C^0_k = P/\hat R T \). For this reason, this method returns a single value, instead of an array. However, for phases in which the standard concentration is species-specific (such as surface species of different sizes), this method may be called with an optional parameter indicating the species.
k | Optional parameter indicating the species. The default is to assume this refers to species 0. |
Reimplemented from ThermoPhase.
Definition at line 41 of file MaskellSolidSolnPhase.h.
|
inlineoverridevirtual |
Natural logarithm of the standard concentration of the kth species.
k | index of the species (defaults to zero) |
Reimplemented from ThermoPhase.
Definition at line 42 of file MaskellSolidSolnPhase.h.
|
overridevirtual |
Molar enthalpy. Units: J/kmol.
Reimplemented from ThermoPhase.
Definition at line 35 of file MaskellSolidSolnPhase.cpp.
|
overridevirtual |
Molar entropy. Units: J/kmol/K.
Reimplemented from ThermoPhase.
Definition at line 48 of file MaskellSolidSolnPhase.cpp.
|
inlineoverridevirtual |
Pressure.
Units: Pa. For this incompressible system, we return the internally stored independent value of the pressure.
Reimplemented from Phase.
Definition at line 65 of file MaskellSolidSolnPhase.h.
|
overridevirtual |
Set the pressure at constant temperature.
Units: Pa. This method sets a constant within the object. The mass density is not a function of pressure.
p | Input Pressure (Pa) |
Reimplemented from Phase.
Definition at line 72 of file MaskellSolidSolnPhase.cpp.
|
overridevirtual |
Calculate the density of the mixture using the partial molar volumes and mole fractions as input.
The formula for this is
\[ \rho = \frac{\sum_k{X_k W_k}}{\sum_k{X_k V_k}} \]
where \( X_k \) are the mole fractions, \( W_k \) are the molecular weights, and \( V_k \) are the pure species molar volumes.
Note, the basis behind this formula is that in an ideal solution the partial molar volumes are equal to the pure species molar volumes. We have additionally specified in this class that the pure species molar volumes are independent of temperature and pressure.
NOTE: This function is not a member of the ThermoPhase base class.
Reimplemented from VPStandardStateTP.
Definition at line 59 of file MaskellSolidSolnPhase.cpp.
|
overridevirtual |
Get the array of non-dimensional molar-based activity coefficients at the current solution temperature, pressure, and solution concentration.
ac | Output vector of activity coefficients. Length: m_kk. |
Reimplemented from ThermoPhase.
Definition at line 79 of file MaskellSolidSolnPhase.cpp.
|
overridevirtual |
Get the species chemical potentials. Units: J/kmol.
This function returns a vector of chemical potentials of the species in solution at the current temperature, pressure and mole fraction of the solution.
mu | Output vector of species chemical potentials. Length: m_kk. Units: J/kmol |
Reimplemented from ThermoPhase.
Definition at line 98 of file MaskellSolidSolnPhase.cpp.
|
overridevirtual |
Reimplemented from ThermoPhase.
Definition at line 111 of file MaskellSolidSolnPhase.cpp.
|
overridevirtual |
Returns an array of partial molar enthalpies for the species in the mixture.
Units (J/kmol)
hbar | Output vector of species partial molar enthalpies. Length: m_kk. units are J/kmol. |
Reimplemented from ThermoPhase.
Definition at line 123 of file MaskellSolidSolnPhase.cpp.
|
overridevirtual |
Returns an array of partial molar entropies of the species in the solution.
Units: J/kmol/K.
sbar | Output vector of species partial molar entropies. Length = m_kk. units are J/kmol/K. |
Reimplemented from ThermoPhase.
Definition at line 128 of file MaskellSolidSolnPhase.cpp.
|
overridevirtual |
Return an array of partial molar heat capacities for the species in the mixture.
Units: J/kmol/K
cpbar | Output vector of species partial molar heat capacities at constant pressure. Length = m_kk. units are J/kmol/K. |
Reimplemented from ThermoPhase.
Definition at line 133 of file MaskellSolidSolnPhase.cpp.
|
overridevirtual |
Return an array of partial molar volumes for the species in the mixture.
Units: m^3/kmol.
vbar | Output vector of species partial molar volumes. Length = m_kk. units are m^3/kmol. |
Reimplemented from ThermoPhase.
Definition at line 138 of file MaskellSolidSolnPhase.cpp.
|
overridevirtual |
Get the Gibbs functions for the standard state of the species at the current T and P of the solution.
Units are Joules/kmol
gpure | Output vector of standard state Gibbs free energies. Length: m_kk. |
Reimplemented from ThermoPhase.
Definition at line 143 of file MaskellSolidSolnPhase.cpp.
|
overridevirtual |
Get the array of chemical potentials at unit activity for the species at their standard states at the current T and P of the solution.
These are the standard state chemical potentials \( \mu^0_k(T,P) \). The values are evaluated at the current temperature and pressure of the solution
mu | Output vector of chemical potentials. Length: m_kk. |
Reimplemented from ThermoPhase.
Definition at line 150 of file MaskellSolidSolnPhase.cpp.
|
overridevirtual |
Initialize the ThermoPhase object after all species have been set up.
This method is provided to allow subclasses to perform any initialization required after all species have been added. For example, it might be used to resize internal work arrays that must have an entry for each species. The base class implementation does nothing, and subclasses that do not require initialization do not need to overload this method. Derived classes which do override this function should call their parent class's implementation of this function as their last action.
When importing from an AnyMap phase description (or from a YAML file), setupPhase() adds all the species, stores the input data in m_input, and then calls this method to set model parameters from the data stored in m_input.
Reimplemented from ThermoPhase.
Definition at line 159 of file MaskellSolidSolnPhase.cpp.
|
overridevirtual |
Store the parameters of a ThermoPhase object such that an identical one could be reconstructed using the newThermo(AnyMap&) function.
This does not include user-defined fields available in input().
Reimplemented from ThermoPhase.
Definition at line 168 of file MaskellSolidSolnPhase.cpp.
|
inline |
Definition at line 107 of file MaskellSolidSolnPhase.h.
void setProductSpecies | ( | const string & | name | ) |
Set the product Species. Must be called after species have been added.
Definition at line 175 of file MaskellSolidSolnPhase.cpp.
|
private |
Definition at line 185 of file MaskellSolidSolnPhase.cpp.
|
private |
Definition at line 190 of file MaskellSolidSolnPhase.cpp.
|
private |
Definition at line 195 of file MaskellSolidSolnPhase.cpp.
|
private |
m_Pcurrent = The current pressure.
Since the density isn't a function of pressure, but only of the mole fractions, we need to independently specify the pressure.
Definition at line 119 of file MaskellSolidSolnPhase.h.
|
private |
Value of the enthalpy change on mixing due to protons changing from type B to type A configurations.
Definition at line 123 of file MaskellSolidSolnPhase.h.
|
private |
Index of the species whose mole fraction defines the extent of reduction r.
Definition at line 126 of file MaskellSolidSolnPhase.h.
|
private |
Definition at line 127 of file MaskellSolidSolnPhase.h.