Generated CLib API for Cantera's ThermoPhase class. More...
Generated CLib API for Cantera's ThermoPhase class.
Functions | |
int32_t | thermo_name (int32_t handle, int32_t bufLen, char *buf) |
Return the name of the phase. | |
int32_t | thermo_setName (int32_t handle, const char *nm) |
Sets the string name for the phase. | |
int32_t | thermo_type (int32_t handle, int32_t bufLen, char *buf) |
String indicating the thermodynamic model implemented. | |
int32_t | thermo_nElements (int32_t handle) |
Number of elements. | |
int32_t | thermo_nSpecies (int32_t handle) |
Returns the number of species in the phase. | |
double | thermo_temperature (int32_t handle) |
Temperature (K). | |
int32_t | thermo_setTemperature (int32_t handle, double temp) |
Set the internally stored temperature of the phase (K). | |
double | thermo_pressure (int32_t handle) |
Return the thermodynamic pressure (Pa). | |
int32_t | thermo_setPressure (int32_t handle, double p) |
Set the internally stored pressure (Pa) at constant temperature and composition. | |
double | thermo_density (int32_t handle) |
Density (kg/m^3). | |
int32_t | thermo_setDensity (int32_t handle, const double density_) |
Set the internally stored density (kg/m^3) of the phase. | |
double | thermo_molarDensity (int32_t handle) |
Molar density (kmol/m^3). | |
double | thermo_meanMolecularWeight (int32_t handle) |
The mean molecular weight. | |
double | thermo_moleFraction (int32_t handle, int32_t k) |
Return the mole fraction of a single species. | |
double | thermo_massFraction (int32_t handle, int32_t k) |
Return the mass fraction of a single species. | |
int32_t | thermo_getMoleFractions (int32_t handle, int32_t xLen, double *x) |
Get the species mole fraction vector. | |
int32_t | thermo_getMassFractions (int32_t handle, int32_t yLen, double *y) |
Get the species mass fractions. | |
int32_t | thermo_setMoleFractions (int32_t handle, int32_t xLen, const double *x) |
Set the mole fractions to the specified values. | |
int32_t | thermo_setMassFractions (int32_t handle, int32_t yLen, const double *y) |
Set the mass fractions to the specified values and normalize them. | |
int32_t | thermo_setMoleFractionsByName (int32_t handle, const char *x) |
Set the mole fractions of a group of species by name. | |
int32_t | thermo_setMassFractionsByName (int32_t handle, const char *x) |
Set the species mass fractions by name. | |
int32_t | thermo_atomicWeights (int32_t handle, int32_t bufLen, double *buf) |
Return a read-only reference to the vector of atomic weights. | |
int32_t | thermo_getMolecularWeights (int32_t handle, int32_t weightsLen, double *weights) |
Copy the vector of molecular weights into array weights. | |
int32_t | thermo_getCharges (int32_t handle, int32_t chargesLen, double *charges) |
Copy the vector of species charges into array charges. | |
int32_t | thermo_elementName (int32_t handle, int32_t m, int32_t bufLen, char *buf) |
Name of the element with index m. | |
int32_t | thermo_speciesName (int32_t handle, int32_t k, int32_t bufLen, char *buf) |
Name of the species with index k. | |
int32_t | thermo_elementIndex (int32_t handle, const char *name) |
Return the index of element named 'name'. | |
int32_t | thermo_speciesIndex (int32_t handle, const char *name) |
Returns the index of a species named 'name' within the Phase object. | |
double | thermo_nAtoms (int32_t handle, int32_t k, int32_t m) |
Number of atoms of element. | |
int32_t | thermo_addElement (int32_t handle, const char *symbol, double weight, int32_t atomicNumber, double entropy298, int32_t elem_type) |
Add an element. | |
double | thermo_refPressure (int32_t handle) |
Returns the reference pressure in Pa. | |
double | thermo_minTemp (int32_t handle, int32_t k) |
Minimum temperature for which the thermodynamic data for the species or phase are valid. | |
double | thermo_maxTemp (int32_t handle, int32_t k) |
Maximum temperature for which the thermodynamic data for the species are valid. | |
double | thermo_enthalpy_mole (int32_t handle) |
Molar enthalpy. | |
double | thermo_enthalpy_mass (int32_t handle) |
Specific enthalpy. | |
double | thermo_entropy_mole (int32_t handle) |
Molar entropy. | |
double | thermo_entropy_mass (int32_t handle) |
Specific entropy. | |
double | thermo_intEnergy_mole (int32_t handle) |
Molar internal energy. | |
double | thermo_intEnergy_mass (int32_t handle) |
Specific internal energy. | |
double | thermo_gibbs_mole (int32_t handle) |
Molar Gibbs function. | |
double | thermo_gibbs_mass (int32_t handle) |
Specific Gibbs function. | |
double | thermo_cp_mole (int32_t handle) |
Molar heat capacity at constant pressure. | |
double | thermo_cp_mass (int32_t handle) |
Specific heat at constant pressure. | |
double | thermo_cv_mole (int32_t handle) |
Molar heat capacity at constant volume. | |
double | thermo_cv_mass (int32_t handle) |
Specific heat at constant volume. | |
int32_t | thermo_getChemPotentials (int32_t handle, int32_t muLen, double *mu) |
Get the species chemical potentials. | |
int32_t | thermo_getElectrochemPotentials (int32_t handle, int32_t muLen, double *mu) |
Get the species electrochemical potentials. | |
double | thermo_electricPotential (int32_t handle) |
Returns the electric potential of this phase (V). | |
int32_t | thermo_setElectricPotential (int32_t handle, double v) |
Set the electric potential of this phase (V). | |
double | thermo_thermalExpansionCoeff (int32_t handle) |
Return the volumetric thermal expansion coefficient. | |
double | thermo_isothermalCompressibility (int32_t handle) |
Returns the isothermal compressibility. | |
int32_t | thermo_getPartialMolarEnthalpies (int32_t handle, int32_t hbarLen, double *hbar) |
Returns an array of partial molar enthalpies for the species in the mixture. | |
int32_t | thermo_getPartialMolarEntropies (int32_t handle, int32_t sbarLen, double *sbar) |
Returns an array of partial molar entropies of the species in the solution. | |
int32_t | thermo_getPartialMolarIntEnergies (int32_t handle, int32_t ubarLen, double *ubar) |
Return an array of partial molar internal energies for the species in the mixture. | |
int32_t | thermo_getPartialMolarCp (int32_t handle, int32_t cpbarLen, double *cpbar) |
Return an array of partial molar heat capacities for the species in the mixture. | |
int32_t | thermo_getPartialMolarVolumes (int32_t handle, int32_t vbarLen, double *vbar) |
Return an array of partial molar volumes for the species in the mixture. | |
int32_t | thermo_setState_TPX (int32_t handle, double t, double p, int32_t xLen, const double *x) |
Set the temperature (K), pressure (Pa), and mole fractions. | |
int32_t | thermo_setState_TPY (int32_t handle, double t, double p, int32_t yLen, const double *y) |
Set the internally stored temperature (K), pressure (Pa), and mass fractions of the phase. | |
int32_t | thermo_setState_TP (int32_t handle, double t, double p) |
Set the temperature (K) and pressure (Pa) | |
int32_t | thermo_setState_TD (int32_t handle, double t, double rho) |
Set the internally stored temperature (K) and density (kg/m^3) | |
int32_t | thermo_setState_DP (int32_t handle, double rho, double p) |
Set the density (kg/m**3) and pressure (Pa) at constant composition. | |
int32_t | thermo_setState_HP (int32_t handle, double h, double p) |
Set the internally stored specific enthalpy (J/kg) and pressure (Pa) of the phase. | |
int32_t | thermo_setState_UV (int32_t handle, double u, double v) |
Set the specific internal energy (J/kg) and specific volume (m^3/kg). | |
int32_t | thermo_setState_SV (int32_t handle, double s, double v) |
Set the specific entropy (J/kg/K) and specific volume (m^3/kg). | |
int32_t | thermo_setState_SP (int32_t handle, double s, double p) |
Set the specific entropy (J/kg/K) and pressure (Pa). | |
int32_t | thermo_setState_ST (int32_t handle, double s, double t) |
Set the specific entropy (J/kg/K) and temperature (K). | |
int32_t | thermo_setState_TV (int32_t handle, double t, double v) |
Set the temperature (K) and specific volume (m^3/kg). | |
int32_t | thermo_setState_PV (int32_t handle, double p, double v) |
Set the pressure (Pa) and specific volume (m^3/kg). | |
int32_t | thermo_setState_UP (int32_t handle, double u, double p) |
Set the specific internal energy (J/kg) and pressure (Pa). | |
int32_t | thermo_setState_VH (int32_t handle, double v, double h) |
Set the specific volume (m^3/kg) and the specific enthalpy (J/kg) | |
int32_t | thermo_setState_TH (int32_t handle, double t, double h) |
Set the temperature (K) and the specific enthalpy (J/kg) | |
int32_t | thermo_setState_SH (int32_t handle, double s, double h) |
Set the specific entropy (J/kg/K) and the specific enthalpy (J/kg) | |
int32_t | thermo_equilibrate (int32_t handle, const char *XY, const char *solver, double rtol, int32_t max_steps, int32_t max_iter, int32_t estimate_equil) |
Equilibrate a ThermoPhase object. | |
double | thermo_critTemperature (int32_t handle) |
Critical temperature (K). | |
double | thermo_critPressure (int32_t handle) |
Critical pressure (Pa). | |
double | thermo_critDensity (int32_t handle) |
Critical density (kg/m3). | |
double | thermo_vaporFraction (int32_t handle) |
Return the fraction of vapor at the current conditions. | |
double | thermo_satTemperature (int32_t handle, double p) |
Return the saturation temperature given the pressure. | |
double | thermo_satPressure (int32_t handle, double t) |
Return the saturation pressure given the temperature. | |
int32_t | thermo_setState_Psat (int32_t handle, double p, double x) |
Set the state to a saturated system at a particular pressure. | |
int32_t | thermo_setState_Tsat (int32_t handle, double t, double x) |
Set the state to a saturated system at a particular temperature. | |
int32_t | surf_getCoverages (int32_t handle, int32_t thetaLen, double *theta) |
Return a vector of surface coverages. | |
int32_t | surf_setCoverages (int32_t handle, int32_t thetaLen, const double *theta) |
Set the surface site fractions to a specified state. | |
int32_t | thermo_getConcentrations (int32_t handle, int32_t cLen, double *c) |
Get the species concentrations (kmol/m^3). | |
int32_t | thermo_setConcentrations (int32_t handle, int32_t concLen, const double *conc) |
Set the concentrations to the specified values within the phase. | |
double | surf_siteDensity (int32_t handle) |
Returns the site density. | |
int32_t | surf_setSiteDensity (int32_t handle, double n0) |
Set the site density of the surface phase (kmol m-2) | |
int32_t | surf_setCoveragesByName (int32_t handle, const char *cov) |
Set the coverages from a string of colon-separated name:value pairs. | |
int32_t | thermo_report (int32_t handle, int32_t show_thermo, double threshold, int32_t bufLen, char *buf) |
returns a summary of the state of the phase as a string | |
int32_t | thermo_print (int32_t handle, int32_t showThermo, double threshold) |
Print a summary of the state of the phase to the logger. | |
int32_t | thermo_del (int32_t handle) |
Destructor; required by some APIs although object is managed by Solution. | |
int32_t | thermo_cabinetSize () |
Return size of ThermoPhase storage. | |
int32_t | thermo_parentHandle (int32_t handle) |
Return handle to parent of ThermoPhase object. | |
int32_t thermo_name | ( | int32_t | handle, |
int32_t | bufLen, | ||
char * | buf | ||
) |
Return the name of the phase.
Wraps C++ getter: string Phase::name()
handle | Handle to queried Phase object. | |
[in] | bufLen | Length of reserved array. |
[out] | buf | Returned string value. |
Definition at line 38 of file ctthermo.cpp.
int32_t thermo_setName | ( | int32_t | handle, |
const char * | nm | ||
) |
Sets the string name for the phase.
Wraps C++ setter: void Phase::setName(const string&)
handle | Handle to queried Phase object. |
nm | String name of the phase |
Definition at line 50 of file ctthermo.cpp.
int32_t thermo_type | ( | int32_t | handle, |
int32_t | bufLen, | ||
char * | buf | ||
) |
String indicating the thermodynamic model implemented.
Wraps C++ getter: string ThermoPhase::type()
handle | Handle to queried ThermoPhase object. | |
[in] | bufLen | Length of reserved array. |
[out] | buf | Returned string value. |
Definition at line 61 of file ctthermo.cpp.
int32_t thermo_nElements | ( | int32_t | handle | ) |
Number of elements.
Wraps C++ getter: size_t Phase::nElements()
handle | Handle to queried Phase object. |
Definition at line 73 of file ctthermo.cpp.
int32_t thermo_nSpecies | ( | int32_t | handle | ) |
Returns the number of species in the phase.
Wraps C++ getter: size_t Phase::nSpecies()
handle | Handle to queried Phase object. |
Definition at line 83 of file ctthermo.cpp.
double thermo_temperature | ( | int32_t | handle | ) |
Temperature (K).
Wraps C++ getter: double Phase::temperature()
handle | Handle to queried Phase object. |
Definition at line 93 of file ctthermo.cpp.
int32_t thermo_setTemperature | ( | int32_t | handle, |
double | temp | ||
) |
Set the internally stored temperature of the phase (K).
Wraps C++ setter: virtual void Phase::setTemperature(double)
handle | Handle to queried Phase object. |
temp | Temperature in Kelvin |
Definition at line 103 of file ctthermo.cpp.
double thermo_pressure | ( | int32_t | handle | ) |
Return the thermodynamic pressure (Pa).
Wraps C++ getter: virtual double Phase::pressure()
handle | Handle to queried Phase object. |
Definition at line 114 of file ctthermo.cpp.
int32_t thermo_setPressure | ( | int32_t | handle, |
double | p | ||
) |
Set the internally stored pressure (Pa) at constant temperature and composition.
Wraps C++ setter: virtual void Phase::setPressure(double)
handle | Handle to queried Phase object. |
p | input Pressure (Pa) |
Definition at line 124 of file ctthermo.cpp.
double thermo_density | ( | int32_t | handle | ) |
Density (kg/m^3).
Wraps C++ getter: virtual double Phase::density()
handle | Handle to queried Phase object. |
Definition at line 135 of file ctthermo.cpp.
int32_t thermo_setDensity | ( | int32_t | handle, |
const double | density_ | ||
) |
Set the internally stored density (kg/m^3) of the phase.
Wraps C++ setter: virtual void Phase::setDensity(const double)
handle | Handle to queried Phase object. | |
[in] | density_ | density (kg/m^3). |
Definition at line 145 of file ctthermo.cpp.
double thermo_molarDensity | ( | int32_t | handle | ) |
Molar density (kmol/m^3).
Wraps C++ getter: virtual double Phase::molarDensity()
handle | Handle to queried Phase object. |
Definition at line 156 of file ctthermo.cpp.
double thermo_meanMolecularWeight | ( | int32_t | handle | ) |
The mean molecular weight.
Units: (kg/kmol)
Wraps C++ getter: double Phase::meanMolecularWeight()
handle | Handle to queried Phase object. |
Definition at line 166 of file ctthermo.cpp.
double thermo_moleFraction | ( | int32_t | handle, |
int32_t | k | ||
) |
Return the mole fraction of a single species.
Wraps C++ method: double Phase::moleFraction(size_t)
handle | Handle to queried Phase object. |
k | species index |
Definition at line 176 of file ctthermo.cpp.
double thermo_massFraction | ( | int32_t | handle, |
int32_t | k | ||
) |
Return the mass fraction of a single species.
Wraps C++ method: double Phase::massFraction(size_t)
handle | Handle to queried Phase object. |
k | species index |
Definition at line 186 of file ctthermo.cpp.
int32_t thermo_getMoleFractions | ( | int32_t | handle, |
int32_t | xLen, | ||
double * | x | ||
) |
Get the species mole fraction vector.
Wraps C++ getter: void Phase::getMoleFractions(double* const)
Uses:
size_t Phase::nSpecies()
handle | Handle to queried Phase object. | |
[in] | xLen | Length of array reserved for x. |
x | On return, x contains the mole fractions. Must have a length greater than or equal to the number of species. |
Definition at line 196 of file ctthermo.cpp.
int32_t thermo_getMassFractions | ( | int32_t | handle, |
int32_t | yLen, | ||
double * | y | ||
) |
Get the species mass fractions.
Wraps C++ getter: void Phase::getMassFractions(double* const)
Uses:
size_t Phase::nSpecies()
handle | Handle to queried Phase object. | |
[in] | yLen | Length of array reserved for y. |
[out] | y | Array of mass fractions, length nSpecies() |
Definition at line 211 of file ctthermo.cpp.
int32_t thermo_setMoleFractions | ( | int32_t | handle, |
int32_t | xLen, | ||
const double * | x | ||
) |
Set the mole fractions to the specified values.
Wraps C++ setter: virtual void Phase::setMoleFractions(const double* const)
Uses:
size_t Phase::nSpecies()
handle | Handle to queried Phase object. | |
[in] | xLen | Length of array reserved for x. |
x | Array of unnormalized mole fraction values (input). Must have a length greater than or equal to the number of species, m_kk. |
Definition at line 226 of file ctthermo.cpp.
int32_t thermo_setMassFractions | ( | int32_t | handle, |
int32_t | yLen, | ||
const double * | y | ||
) |
Set the mass fractions to the specified values and normalize them.
Wraps C++ setter: virtual void Phase::setMassFractions(const double* const)
Uses:
size_t Phase::nSpecies()
handle | Handle to queried Phase object. | |
[in] | yLen | Length of array reserved for y. |
[in] | y | Array of unnormalized mass fraction values. Length must be greater than or equal to the number of species. The Phase object will normalize this vector before storing its contents. |
Definition at line 241 of file ctthermo.cpp.
int32_t thermo_setMoleFractionsByName | ( | int32_t | handle, |
const char * | x | ||
) |
Set the mole fractions of a group of species by name.
Wraps C++ setter: void Phase::setMoleFractionsByName(const string&)
handle | Handle to queried Phase object. |
x | string x in the form of a composition map |
Definition at line 256 of file ctthermo.cpp.
int32_t thermo_setMassFractionsByName | ( | int32_t | handle, |
const char * | x | ||
) |
Set the species mass fractions by name.
Wraps C++ setter: void Phase::setMassFractionsByName(const string&)
handle | Handle to queried Phase object. |
x | String containing a composition map |
Definition at line 267 of file ctthermo.cpp.
int32_t thermo_atomicWeights | ( | int32_t | handle, |
int32_t | bufLen, | ||
double * | buf | ||
) |
Return a read-only reference to the vector of atomic weights.
Wraps C++ getter: const vector<double>& Phase::atomicWeights()
Uses:
size_t Phase::nElements()
handle | Handle to queried Phase object. | |
[in] | bufLen | Length of reserved array. |
[out] | buf | Returned array value. |
Definition at line 278 of file ctthermo.cpp.
int32_t thermo_getMolecularWeights | ( | int32_t | handle, |
int32_t | weightsLen, | ||
double * | weights | ||
) |
Copy the vector of molecular weights into array weights.
Wraps C++ getter: void Phase::getMolecularWeights(double*)
Uses:
size_t Phase::nSpecies()
handle | Handle to queried Phase object. | |
[in] | weightsLen | Length of array reserved for weights. |
weights | Output array of molecular weights (kg/kmol) |
Definition at line 294 of file ctthermo.cpp.
int32_t thermo_getCharges | ( | int32_t | handle, |
int32_t | chargesLen, | ||
double * | charges | ||
) |
Copy the vector of species charges into array charges.
Wraps C++ getter: void Phase::getCharges(double*)
Uses:
size_t Phase::nElements()
handle | Handle to queried Phase object. | |
[in] | chargesLen | Length of array reserved for charges. |
charges | Output array of species charges (elem. charge) |
Definition at line 309 of file ctthermo.cpp.
int32_t thermo_elementName | ( | int32_t | handle, |
int32_t | m, | ||
int32_t | bufLen, | ||
char * | buf | ||
) |
Name of the element with index m.
Wraps C++ method: string Phase::elementName(size_t)
handle | Handle to queried Phase object. | |
m | Element index. | |
[in] | bufLen | Length of reserved array. |
[out] | buf | Returned string value. |
Definition at line 324 of file ctthermo.cpp.
int32_t thermo_speciesName | ( | int32_t | handle, |
int32_t | k, | ||
int32_t | bufLen, | ||
char * | buf | ||
) |
Name of the species with index k.
Wraps C++ method: string Phase::speciesName(size_t)
handle | Handle to queried Phase object. | |
k | index of the species | |
[in] | bufLen | Length of reserved array. |
[out] | buf | Returned string value. |
Definition at line 336 of file ctthermo.cpp.
int32_t thermo_elementIndex | ( | int32_t | handle, |
const char * | name | ||
) |
Return the index of element named 'name'.
Wraps C++ method: size_t Phase::elementIndex(const string&)
handle | Handle to queried Phase object. |
name | Name of the element |
Definition at line 348 of file ctthermo.cpp.
int32_t thermo_speciesIndex | ( | int32_t | handle, |
const char * | name | ||
) |
Returns the index of a species named 'name' within the Phase object.
Wraps C++ method: size_t Phase::speciesIndex(const string&)
handle | Handle to queried Phase object. |
name | String name of the species. It may also be in the form phaseName:speciesName |
Definition at line 358 of file ctthermo.cpp.
double thermo_nAtoms | ( | int32_t | handle, |
int32_t | k, | ||
int32_t | m | ||
) |
Number of atoms of element.
Wraps C++ method: double Phase::nAtoms(size_t, size_t)
handle | Handle to queried Phase object. |
k | species index |
m | element index |
Definition at line 368 of file ctthermo.cpp.
int32_t thermo_addElement | ( | int32_t | handle, |
const char * | symbol, | ||
double | weight, | ||
int32_t | atomicNumber, | ||
double | entropy298, | ||
int32_t | elem_type | ||
) |
Add an element.
Wraps C++ method: size_t Phase::addElement(const string&, double, int, double, int)
handle | Handle to queried Phase object. |
symbol | Atomic symbol string. |
weight | Atomic mass in amu. |
atomicNumber | Atomic number of the element (unitless) |
entropy298 | Entropy of the element at 298 K and 1 bar in its most stable form. The default is the value ENTROPY298_UNKNOWN, which is interpreted as an unknown, and if used will cause Cantera to throw an error. |
elem_type | Specifies the type of the element constraint equation. This defaults to CT_ELEM_TYPE_ABSPOS, that is, an element. |
Definition at line 378 of file ctthermo.cpp.
double thermo_refPressure | ( | int32_t | handle | ) |
Returns the reference pressure in Pa.
Wraps C++ getter: virtual double ThermoPhase::refPressure()
handle | Handle to queried ThermoPhase object. |
Definition at line 388 of file ctthermo.cpp.
double thermo_minTemp | ( | int32_t | handle, |
int32_t | k | ||
) |
Minimum temperature for which the thermodynamic data for the species or phase are valid.
Wraps C++ method: virtual double ThermoPhase::minTemp(size_t)
handle | Handle to queried ThermoPhase object. |
k | index of the species. Default is -1, which will return the max of the min value over all species. |
Definition at line 398 of file ctthermo.cpp.
double thermo_maxTemp | ( | int32_t | handle, |
int32_t | k | ||
) |
Maximum temperature for which the thermodynamic data for the species are valid.
Wraps C++ method: virtual double ThermoPhase::maxTemp(size_t)
handle | Handle to queried ThermoPhase object. |
k | index of the species. Default is -1, which will return the min of the max value over all species. |
Definition at line 408 of file ctthermo.cpp.
double thermo_enthalpy_mole | ( | int32_t | handle | ) |
Molar enthalpy.
Units: J/kmol.
Wraps C++ getter: virtual double ThermoPhase::enthalpy_mole()
handle | Handle to queried ThermoPhase object. |
Definition at line 418 of file ctthermo.cpp.
double thermo_enthalpy_mass | ( | int32_t | handle | ) |
Specific enthalpy.
Units: J/kg.
Wraps C++ getter: double ThermoPhase::enthalpy_mass()
handle | Handle to queried ThermoPhase object. |
Definition at line 428 of file ctthermo.cpp.
double thermo_entropy_mole | ( | int32_t | handle | ) |
Molar entropy.
Units: J/kmol/K.
Wraps C++ getter: virtual double ThermoPhase::entropy_mole()
handle | Handle to queried ThermoPhase object. |
Definition at line 438 of file ctthermo.cpp.
double thermo_entropy_mass | ( | int32_t | handle | ) |
Specific entropy.
Units: J/kg/K.
Wraps C++ getter: double ThermoPhase::entropy_mass()
handle | Handle to queried ThermoPhase object. |
Definition at line 448 of file ctthermo.cpp.
double thermo_intEnergy_mole | ( | int32_t | handle | ) |
Molar internal energy.
Units: J/kmol.
Wraps C++ getter: virtual double ThermoPhase::intEnergy_mole()
handle | Handle to queried ThermoPhase object. |
Definition at line 458 of file ctthermo.cpp.
double thermo_intEnergy_mass | ( | int32_t | handle | ) |
Specific internal energy.
Units: J/kg.
Wraps C++ getter: double ThermoPhase::intEnergy_mass()
handle | Handle to queried ThermoPhase object. |
Definition at line 468 of file ctthermo.cpp.
double thermo_gibbs_mole | ( | int32_t | handle | ) |
Molar Gibbs function.
Units: J/kmol.
Wraps C++ getter: virtual double ThermoPhase::gibbs_mole()
handle | Handle to queried ThermoPhase object. |
Definition at line 478 of file ctthermo.cpp.
double thermo_gibbs_mass | ( | int32_t | handle | ) |
Specific Gibbs function.
Units: J/kg.
Wraps C++ getter: double ThermoPhase::gibbs_mass()
handle | Handle to queried ThermoPhase object. |
Definition at line 488 of file ctthermo.cpp.
double thermo_cp_mole | ( | int32_t | handle | ) |
Molar heat capacity at constant pressure.
Units: J/kmol/K.
Wraps C++ getter: virtual double ThermoPhase::cp_mole()
handle | Handle to queried ThermoPhase object. |
Definition at line 498 of file ctthermo.cpp.
double thermo_cp_mass | ( | int32_t | handle | ) |
Specific heat at constant pressure.
Units: J/kg/K.
Wraps C++ getter: double ThermoPhase::cp_mass()
handle | Handle to queried ThermoPhase object. |
Definition at line 508 of file ctthermo.cpp.
double thermo_cv_mole | ( | int32_t | handle | ) |
Molar heat capacity at constant volume.
Units: J/kmol/K.
Wraps C++ getter: virtual double ThermoPhase::cv_mole()
handle | Handle to queried ThermoPhase object. |
Definition at line 518 of file ctthermo.cpp.
double thermo_cv_mass | ( | int32_t | handle | ) |
Specific heat at constant volume.
Units: J/kg/K.
Wraps C++ getter: double ThermoPhase::cv_mass()
handle | Handle to queried ThermoPhase object. |
Definition at line 528 of file ctthermo.cpp.
int32_t thermo_getChemPotentials | ( | int32_t | handle, |
int32_t | muLen, | ||
double * | mu | ||
) |
Get the species chemical potentials.
Units: J/kmol.
Wraps C++ getter: virtual void ThermoPhase::getChemPotentials(double*)
handle | Handle to queried ThermoPhase object. | |
[in] | muLen | Length of array reserved for mu. |
mu | Output vector of species chemical potentials. Length: m_kk. Units: J/kmol |
Definition at line 538 of file ctthermo.cpp.
int32_t thermo_getElectrochemPotentials | ( | int32_t | handle, |
int32_t | muLen, | ||
double * | mu | ||
) |
Get the species electrochemical potentials.
Wraps C++ getter: void ThermoPhase::getElectrochemPotentials(double*)
handle | Handle to queried ThermoPhase object. | |
[in] | muLen | Length of array reserved for mu. |
mu | Output vector of species electrochemical potentials. Length: m_kk. Units: J/kmol |
Definition at line 551 of file ctthermo.cpp.
double thermo_electricPotential | ( | int32_t | handle | ) |
Returns the electric potential of this phase (V).
Wraps C++ getter: double ThermoPhase::electricPotential()
handle | Handle to queried ThermoPhase object. |
Definition at line 564 of file ctthermo.cpp.
int32_t thermo_setElectricPotential | ( | int32_t | handle, |
double | v | ||
) |
Set the electric potential of this phase (V).
Wraps C++ setter: void ThermoPhase::setElectricPotential(double)
handle | Handle to queried ThermoPhase object. |
v | Input value of the electric potential in Volts |
Definition at line 574 of file ctthermo.cpp.
double thermo_thermalExpansionCoeff | ( | int32_t | handle | ) |
Return the volumetric thermal expansion coefficient.
Units: 1/K.
Wraps C++ getter: virtual double ThermoPhase::thermalExpansionCoeff()
handle | Handle to queried ThermoPhase object. |
Definition at line 585 of file ctthermo.cpp.
double thermo_isothermalCompressibility | ( | int32_t | handle | ) |
Returns the isothermal compressibility.
Units: 1/Pa.
Wraps C++ getter: virtual double ThermoPhase::isothermalCompressibility()
handle | Handle to queried ThermoPhase object. |
Definition at line 595 of file ctthermo.cpp.
int32_t thermo_getPartialMolarEnthalpies | ( | int32_t | handle, |
int32_t | hbarLen, | ||
double * | hbar | ||
) |
Returns an array of partial molar enthalpies for the species in the mixture.
Wraps C++ getter: virtual void ThermoPhase::getPartialMolarEnthalpies(double*)
handle | Handle to queried ThermoPhase object. | |
[in] | hbarLen | Length of array reserved for hbar. |
hbar | Output vector of species partial molar enthalpies. Length: m_kk. units are J/kmol. |
Definition at line 605 of file ctthermo.cpp.
int32_t thermo_getPartialMolarEntropies | ( | int32_t | handle, |
int32_t | sbarLen, | ||
double * | sbar | ||
) |
Returns an array of partial molar entropies of the species in the solution.
Wraps C++ getter: virtual void ThermoPhase::getPartialMolarEntropies(double*)
handle | Handle to queried ThermoPhase object. | |
[in] | sbarLen | Length of array reserved for sbar. |
sbar | Output vector of species partial molar entropies. Length = m_kk. units are J/kmol/K. |
Definition at line 618 of file ctthermo.cpp.
int32_t thermo_getPartialMolarIntEnergies | ( | int32_t | handle, |
int32_t | ubarLen, | ||
double * | ubar | ||
) |
Return an array of partial molar internal energies for the species in the mixture.
Wraps C++ getter: virtual void ThermoPhase::getPartialMolarIntEnergies(double*)
handle | Handle to queried ThermoPhase object. | |
[in] | ubarLen | Length of array reserved for ubar. |
ubar | Output vector of species partial molar internal energies. Length = m_kk. units are J/kmol. |
Definition at line 631 of file ctthermo.cpp.
int32_t thermo_getPartialMolarCp | ( | int32_t | handle, |
int32_t | cpbarLen, | ||
double * | cpbar | ||
) |
Return an array of partial molar heat capacities for the species in the mixture.
Wraps C++ getter: virtual void ThermoPhase::getPartialMolarCp(double*)
handle | Handle to queried ThermoPhase object. | |
[in] | cpbarLen | Length of array reserved for cpbar. |
cpbar | Output vector of species partial molar heat capacities at constant pressure. Length = m_kk. units are J/kmol/K. |
Definition at line 644 of file ctthermo.cpp.
int32_t thermo_getPartialMolarVolumes | ( | int32_t | handle, |
int32_t | vbarLen, | ||
double * | vbar | ||
) |
Return an array of partial molar volumes for the species in the mixture.
Wraps C++ getter: virtual void ThermoPhase::getPartialMolarVolumes(double*)
handle | Handle to queried ThermoPhase object. | |
[in] | vbarLen | Length of array reserved for vbar. |
vbar | Output vector of species partial molar volumes. Length = m_kk. units are m^3/kmol. |
Definition at line 657 of file ctthermo.cpp.
int32_t thermo_setState_TPX | ( | int32_t | handle, |
double | t, | ||
double | p, | ||
int32_t | xLen, | ||
const double * | x | ||
) |
Set the temperature (K), pressure (Pa), and mole fractions.
Wraps C++ method: virtual void ThermoPhase::setState_TPX(double, double, const double*)
handle | Handle to queried ThermoPhase object. | |
t | Temperature (K) | |
p | Pressure (Pa) | |
[in] | xLen | Length of array reserved for x. |
x | Vector of mole fractions. Length is equal to m_kk. |
Definition at line 670 of file ctthermo.cpp.
int32_t thermo_setState_TPY | ( | int32_t | handle, |
double | t, | ||
double | p, | ||
int32_t | yLen, | ||
const double * | y | ||
) |
Set the internally stored temperature (K), pressure (Pa), and mass fractions of the phase.
Wraps C++ method: virtual void ThermoPhase::setState_TPY(double, double, const double*)
handle | Handle to queried ThermoPhase object. | |
t | Temperature (K) | |
p | Pressure (Pa) | |
[in] | yLen | Length of array reserved for y. |
y | Vector of mass fractions. Length is equal to m_kk. |
Definition at line 681 of file ctthermo.cpp.
int32_t thermo_setState_TP | ( | int32_t | handle, |
double | t, | ||
double | p | ||
) |
Set the temperature (K) and pressure (Pa)
Wraps C++ method: virtual void ThermoPhase::setState_TP(double, double)
handle | Handle to queried ThermoPhase object. |
t | Temperature (K) |
p | Pressure (Pa) |
Definition at line 692 of file ctthermo.cpp.
int32_t thermo_setState_TD | ( | int32_t | handle, |
double | t, | ||
double | rho | ||
) |
Set the internally stored temperature (K) and density (kg/m^3)
Wraps C++ method: void Phase::setState_TD(double, double)
handle | Handle to queried Phase object. |
t | Temperature in kelvin |
rho | Density (kg/m^3) |
Definition at line 703 of file ctthermo.cpp.
int32_t thermo_setState_DP | ( | int32_t | handle, |
double | rho, | ||
double | p | ||
) |
Set the density (kg/m**3) and pressure (Pa) at constant composition.
Wraps C++ method: virtual void ThermoPhase::setState_DP(double, double)
handle | Handle to queried ThermoPhase object. |
rho | Density (kg/m^3) |
p | Pressure (Pa) |
Definition at line 714 of file ctthermo.cpp.
int32_t thermo_setState_HP | ( | int32_t | handle, |
double | h, | ||
double | p | ||
) |
Set the internally stored specific enthalpy (J/kg) and pressure (Pa) of the phase.
Wraps C++ method: virtual void ThermoPhase::setState_HP(double, double)
handle | Handle to queried ThermoPhase object. |
h | Specific enthalpy (J/kg) |
p | Pressure (Pa) |
Definition at line 725 of file ctthermo.cpp.
int32_t thermo_setState_UV | ( | int32_t | handle, |
double | u, | ||
double | v | ||
) |
Set the specific internal energy (J/kg) and specific volume (m^3/kg).
Wraps C++ method: virtual void ThermoPhase::setState_UV(double, double)
handle | Handle to queried ThermoPhase object. |
u | specific internal energy (J/kg) |
v | specific volume (m^3/kg). |
Definition at line 736 of file ctthermo.cpp.
int32_t thermo_setState_SV | ( | int32_t | handle, |
double | s, | ||
double | v | ||
) |
Set the specific entropy (J/kg/K) and specific volume (m^3/kg).
Wraps C++ method: virtual void ThermoPhase::setState_SV(double, double)
handle | Handle to queried ThermoPhase object. |
s | specific entropy (J/kg/K) |
v | specific volume (m^3/kg). |
Definition at line 747 of file ctthermo.cpp.
int32_t thermo_setState_SP | ( | int32_t | handle, |
double | s, | ||
double | p | ||
) |
Set the specific entropy (J/kg/K) and pressure (Pa).
Wraps C++ method: virtual void ThermoPhase::setState_SP(double, double)
handle | Handle to queried ThermoPhase object. |
s | specific entropy (J/kg/K) |
p | specific pressure (Pa). |
Definition at line 758 of file ctthermo.cpp.
int32_t thermo_setState_ST | ( | int32_t | handle, |
double | s, | ||
double | t | ||
) |
Set the specific entropy (J/kg/K) and temperature (K).
Wraps C++ method: virtual void ThermoPhase::setState_ST(double, double)
handle | Handle to queried ThermoPhase object. |
s | specific entropy (J/kg/K) |
t | temperature (K) |
Definition at line 769 of file ctthermo.cpp.
int32_t thermo_setState_TV | ( | int32_t | handle, |
double | t, | ||
double | v | ||
) |
Set the temperature (K) and specific volume (m^3/kg).
Wraps C++ method: virtual void ThermoPhase::setState_TV(double, double)
handle | Handle to queried ThermoPhase object. |
t | temperature (K) |
v | specific volume (m^3/kg) |
Definition at line 780 of file ctthermo.cpp.
int32_t thermo_setState_PV | ( | int32_t | handle, |
double | p, | ||
double | v | ||
) |
Set the pressure (Pa) and specific volume (m^3/kg).
Wraps C++ method: virtual void ThermoPhase::setState_PV(double, double)
handle | Handle to queried ThermoPhase object. |
p | pressure (Pa) |
v | specific volume (m^3/kg) |
Definition at line 791 of file ctthermo.cpp.
int32_t thermo_setState_UP | ( | int32_t | handle, |
double | u, | ||
double | p | ||
) |
Set the specific internal energy (J/kg) and pressure (Pa).
Wraps C++ method: virtual void ThermoPhase::setState_UP(double, double)
handle | Handle to queried ThermoPhase object. |
u | specific internal energy (J/kg) |
p | pressure (Pa) |
Definition at line 802 of file ctthermo.cpp.
int32_t thermo_setState_VH | ( | int32_t | handle, |
double | v, | ||
double | h | ||
) |
Set the specific volume (m^3/kg) and the specific enthalpy (J/kg)
Wraps C++ method: virtual void ThermoPhase::setState_VH(double, double)
handle | Handle to queried ThermoPhase object. |
v | specific volume (m^3/kg) |
h | specific enthalpy (J/kg) |
Definition at line 813 of file ctthermo.cpp.
int32_t thermo_setState_TH | ( | int32_t | handle, |
double | t, | ||
double | h | ||
) |
Set the temperature (K) and the specific enthalpy (J/kg)
Wraps C++ method: virtual void ThermoPhase::setState_TH(double, double)
handle | Handle to queried ThermoPhase object. |
t | temperature (K) |
h | specific enthalpy (J/kg) |
Definition at line 824 of file ctthermo.cpp.
int32_t thermo_setState_SH | ( | int32_t | handle, |
double | s, | ||
double | h | ||
) |
Set the specific entropy (J/kg/K) and the specific enthalpy (J/kg)
Wraps C++ method: virtual void ThermoPhase::setState_SH(double, double)
handle | Handle to queried ThermoPhase object. |
s | specific entropy (J/kg/K) |
h | specific enthalpy (J/kg) |
Definition at line 835 of file ctthermo.cpp.
int32_t thermo_equilibrate | ( | int32_t | handle, |
const char * | XY, | ||
const char * | solver, | ||
double | rtol, | ||
int32_t | max_steps, | ||
int32_t | max_iter, | ||
int32_t | estimate_equil | ||
) |
Equilibrate a ThermoPhase object.
Wraps C++ method: void ThermoPhase::equilibrate(const string&, const string&, double, int, int, int)
handle | Handle to queried ThermoPhase object. |
XY | String representation of what two properties are being held constant |
solver | Name of the solver to be used to equilibrate the phase. If solver = 'element_potential', the ChemEquil element potential solver will be used. If solver = 'vcs', the VCS solver will be used. If solver = 'gibbs', the MultiPhaseEquil solver will be used. If solver = 'auto', the solvers will be tried in order if the initial solver(s) fail. |
rtol | Relative tolerance |
max_steps | Maximum number of steps to take to find the solution |
max_iter | For the 'gibbs' and 'vcs' solvers, this is the maximum number of outer temperature or pressure iterations to take when T and/or P is not held fixed. |
estimate_equil | For MultiPhaseEquil solver, an integer indicating whether the solver should estimate its own initial condition. If 0, the initial mole fraction vector in the ThermoPhase object is used as the initial condition. If 1, the initial mole fraction vector is used if the element abundances are satisfied. If -1, the initial mole fraction vector is thrown out, and an estimate is formulated. |
Definition at line 846 of file ctthermo.cpp.
double thermo_critTemperature | ( | int32_t | handle | ) |
Critical temperature (K).
Wraps C++ getter: virtual double ThermoPhase::critTemperature()
handle | Handle to queried ThermoPhase object. |
Definition at line 857 of file ctthermo.cpp.
double thermo_critPressure | ( | int32_t | handle | ) |
Critical pressure (Pa).
Wraps C++ getter: virtual double ThermoPhase::critPressure()
handle | Handle to queried ThermoPhase object. |
Definition at line 867 of file ctthermo.cpp.
double thermo_critDensity | ( | int32_t | handle | ) |
Critical density (kg/m3).
Wraps C++ getter: virtual double ThermoPhase::critDensity()
handle | Handle to queried ThermoPhase object. |
Definition at line 877 of file ctthermo.cpp.
double thermo_vaporFraction | ( | int32_t | handle | ) |
Return the fraction of vapor at the current conditions.
Wraps C++ getter: virtual double ThermoPhase::vaporFraction()
handle | Handle to queried ThermoPhase object. |
Definition at line 887 of file ctthermo.cpp.
double thermo_satTemperature | ( | int32_t | handle, |
double | p | ||
) |
Return the saturation temperature given the pressure.
Wraps C++ method: virtual double ThermoPhase::satTemperature(double)
handle | Handle to queried ThermoPhase object. |
p | Pressure (Pa) |
Definition at line 897 of file ctthermo.cpp.
double thermo_satPressure | ( | int32_t | handle, |
double | t | ||
) |
Return the saturation pressure given the temperature.
Wraps C++ method: virtual double ThermoPhase::satPressure(double)
handle | Handle to queried ThermoPhase object. |
t | Temperature (Kelvin) |
Definition at line 907 of file ctthermo.cpp.
int32_t thermo_setState_Psat | ( | int32_t | handle, |
double | p, | ||
double | x | ||
) |
Set the state to a saturated system at a particular pressure.
Wraps C++ method: virtual void ThermoPhase::setState_Psat(double, double)
handle | Handle to queried ThermoPhase object. |
p | Pressure (Pa) |
x | Fraction of vapor |
Definition at line 917 of file ctthermo.cpp.
int32_t thermo_setState_Tsat | ( | int32_t | handle, |
double | t, | ||
double | x | ||
) |
Set the state to a saturated system at a particular temperature.
Wraps C++ method: virtual void ThermoPhase::setState_Tsat(double, double)
handle | Handle to queried ThermoPhase object. |
t | Temperature (kelvin) |
x | Fraction of vapor |
Definition at line 928 of file ctthermo.cpp.
int32_t surf_getCoverages | ( | int32_t | handle, |
int32_t | thetaLen, | ||
double * | theta | ||
) |
Return a vector of surface coverages.
Wraps C++ getter: void SurfPhase::getCoverages(double*)
handle | Handle to queried SurfPhase object. | |
[in] | thetaLen | Length of array reserved for theta. |
theta | Array theta must be at least as long as the number of species. |
Definition at line 939 of file ctthermo.cpp.
int32_t surf_setCoverages | ( | int32_t | handle, |
int32_t | thetaLen, | ||
const double * | theta | ||
) |
Set the surface site fractions to a specified state.
Wraps C++ setter: void SurfPhase::setCoverages(const double*)
handle | Handle to queried SurfPhase object. | |
[in] | thetaLen | Length of array reserved for theta. |
theta | This is the surface site fraction for the kth species in the surface phase. This is a dimensionless quantity. |
Definition at line 952 of file ctthermo.cpp.
int32_t thermo_getConcentrations | ( | int32_t | handle, |
int32_t | cLen, | ||
double * | c | ||
) |
Get the species concentrations (kmol/m^3).
Wraps C++ getter: virtual void Phase::getConcentrations(double* const)
handle | Handle to queried Phase object. | |
[in] | cLen | Length of array reserved for c. |
[out] | c | The vector of species concentrations. Units are kmol/m^3. The length of the vector must be greater than or equal to the number of species within the phase. |
Definition at line 965 of file ctthermo.cpp.
int32_t thermo_setConcentrations | ( | int32_t | handle, |
int32_t | concLen, | ||
const double * | conc | ||
) |
Set the concentrations to the specified values within the phase.
Wraps C++ setter: virtual void Phase::setConcentrations(const double* const)
handle | Handle to queried Phase object. | |
[in] | concLen | Length of array reserved for conc. |
[in] | conc | Array of concentrations in dimensional units. For bulk phases c[k] is the concentration of the kth species in kmol/m3. For surface phases, c[k] is the concentration in kmol/m2. The length of the vector is the number of species in the phase. |
Definition at line 978 of file ctthermo.cpp.
double surf_siteDensity | ( | int32_t | handle | ) |
Returns the site density.
Wraps C++ getter: double SurfPhase::siteDensity()
handle | Handle to queried SurfPhase object. |
Definition at line 991 of file ctthermo.cpp.
int32_t surf_setSiteDensity | ( | int32_t | handle, |
double | n0 | ||
) |
Set the site density of the surface phase (kmol m-2)
Wraps C++ setter: void SurfPhase::setSiteDensity(double)
handle | Handle to queried SurfPhase object. |
n0 | Site density of the surface phase (kmol m-2) |
Definition at line 1001 of file ctthermo.cpp.
int32_t surf_setCoveragesByName | ( | int32_t | handle, |
const char * | cov | ||
) |
Set the coverages from a string of colon-separated name:value pairs.
Wraps C++ setter: void SurfPhase::setCoveragesByName(const string&)
handle | Handle to queried SurfPhase object. |
cov | String containing colon-separated name:value pairs |
Definition at line 1012 of file ctthermo.cpp.
int32_t thermo_report | ( | int32_t | handle, |
int32_t | show_thermo, | ||
double | threshold, | ||
int32_t | bufLen, | ||
char * | buf | ||
) |
returns a summary of the state of the phase as a string
Wraps C++ method: virtual string ThermoPhase::report(bool, double)
handle | Handle to queried ThermoPhase object. | |
show_thermo | If true, extra information is printed out about the thermodynamic state of the system. | |
threshold | Show information about species with mole fractions greater than | |
[in] | bufLen | Length of reserved array. |
[out] | buf | Returned string value. |
Definition at line 1023 of file ctthermo.cpp.
int32_t thermo_print | ( | int32_t | handle, |
int32_t | showThermo, | ||
double | threshold | ||
) |
Print a summary of the state of the phase to the logger.
Wraps C++ method: custom code
Uses:
virtual string ThermoPhase::report(bool, double)
handle | Handle to queried ThermoPhase object. |
showThermo | If true, extra information is printed out about the thermodynamic state of the system. |
threshold | Show information about species with mole fractions greater than |
Definition at line 1036 of file ctthermo.cpp.
int32_t thermo_del | ( | int32_t | handle | ) |
Destructor; required by some APIs although object is managed by Solution.
Wraps C++ noop: undefined
handle | Handle to ThermoPhase object. |
Definition at line 1050 of file ctthermo.cpp.
int32_t thermo_cabinetSize | ( | ) |
Return size of ThermoPhase storage.
Wraps C++ reserved CLib function: custom code
Definition at line 1056 of file ctthermo.cpp.
int32_t thermo_parentHandle | ( | int32_t | handle | ) |
Return handle to parent of ThermoPhase object.
Wraps C++ reserved CLib function: custom code
handle | Handle to queried ThermoPhase object. |
Definition at line 1068 of file ctthermo.cpp.