Thermodynamic Properties¶
ThermoPhase¶
- class ThermoPhase(src, id)¶
ThermoPhase class constructor.
- Parameters
src – Input string of YAML, CTI, or XML file name.
id – ID of the phase to import as specified in the input file. (optional)
- Returns
Instance of class
ThermoPhase()
- atomicMasses(tp)¶
Get the atomic masses of the elements.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase).- Returns
Vector of element atomic masses. Units: kg/kmol
- charges(tp)¶
Get the array of species charges
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Vector of species charges. Units: elem. charge
- chemPotentials(tp)¶
Get the chemical potentials of the species.
The expressions used to compute the chemical potential depend on the model implemented by the underlying kernel thermo manager.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase).- Returns
Vector of species chemical potentials. Units: J/kmol
This method returns an array containing the species chemical potentials [J/kmol]. The expressions used to compute these depend on the model implemented by the underlying kernel thermo manager.”””
- cp_R(tp)¶
Get the non-dimensional specific heats at constant pressure.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Vector of specific heats of the species at constant pressure, non-dimensional basis
- cp_mass(tp)¶
Get the mass-basis specific heat at constant pressure.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Mass basis specific heat of the mixture at constant pressure. Units: J/kg-K
- cp_mole(tp)¶
Get the molar-basis specific heat at constant pressure.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Molar basis specific heat of the mixture at constant pressure. Units: J/kmol-K
- critDensity(tp)¶
Get the critical density.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Critical density. Units: kg/m**3
- critPressure(tp)¶
Get the critical pressure.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Critical pressure. Units: Pa
- critTemperature(tp)¶
Get the critical temperature.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Critical temperature. Units: K
- cv_mass(tp)¶
Get the mass-basis specific heat at constant volume.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Mass basis specific heat of the mixture at constant volume. Units: J/kg-K
- cv_mole(tp)¶
Get the molar-basis specific heat at constant volume.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Molar basis specific heat of the mixture at constant volume. Units: J/kmol-K
- density(tp)¶
Get the density.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Mass density. Units: kg/m**3
- electricPotential(tp)¶
Get the electric potential.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
The electric potential of the phase. Units: V
- elementIndex(tp, name)¶
Get the index of an element given its name.
The index is an integer assigned to each element in sequence as it is read in from the input file.
If
name
is a single string, the return value will be a integer containing the corresponding index. If it is an cell array of strings, the output will be an array of the same shape containing the indices.NOTE: In keeping with the conventions used by Matlab, this method returns 1 for the first element. In contrast, the corresponding method elementIndex in the Cantera C++ and Python interfaces returns 0 for the first element, 1 for the second one, etc.
>> ic = elementIndex(gas, 'C'); >> ih = elementIndex(gas, 'H');
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)name – String or cell array of strings of elements to look up
- Returns
Integer or vector of integers of element indices
- elementName(tp, m)¶
Get the name of an element given its index.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)m – Scalar or vector of integers of element indices
- Returns
If m is a scalar integer, the return value will be a string containing the name of the m^th species. If it is an array of integers, the output will be a cell array of the same shape containing the name strings.
- elementalMassFraction(tp, element)¶
Determine the elemental mass fraction in gas object.
- Parameters
tp – Object representing the gas, instance of class
Solution()
, and an ideal gas. The state of this object should be set to an estimate of the gas state before calling elementalMassFraction.element – String representing the element name.
- Returns
Elemental mass fraction within a gas object.
- enthalpies_RT(tp)¶
Get the non-dimensional enthalpies.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Vector of standard-state species enthalpies divided by RT, where R is the universal gas constant and T is the temperature. For gaseous species, these values are ideal gas enthalpies.
- enthalpy_mass(tp)¶
Get the mass specific enthalpy.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Mass specific enthalpy of the mixture. Units: J/kg
- enthalpy_mole(tp)¶
Get the mole specific enthalpy.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Molar specific enthalpy of the mixture. Units: J/kmol
- entropies_R(tp)¶
Get the non-dimensional entropy.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Vector of species non-dimensional entropies.
- entropy_mass(tp)¶
Get the mass specific entropy.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Mass specific entropy of the mixture. Units: J/kg-K
- entropy_mole(tp)¶
Get the mole specific entropy.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Molar specific entropy of the mixture. Units: J/kmol-K
- eosType(tp)¶
Get the type of the equation of state.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
An string identifying the equation of state.
- equilibrate(tp, xy, solver, rtol, maxsteps, maxiter, loglevel)¶
Set the phase to a state of chemical equilibrium.
- Parameters
XY – A two-letter string, which must be one of the set
['TP','TV','HP','SP','SV','UV','UP']
, indicating which pair of properties should be held constant. Not all of the properties to be held constant are available with all of the solvers.solver – Specifies the equilibrium solver to use. If solver = 0, a fast solver using the element potential method will be used. If solver = 1, a slower but more robust Gibbs minimization solver will be used. If solver >= 2, a version of the VCS algorithm will be used. If solver < 0 or is unspecified, the fast solver will be tried first, then if it fails the Gibbs minimization solver will be tried.
rtol – The relative error tolerance.
maxsteps – Maximum number of steps in composition to take to find a converged solution.
maxiter – For the Gibbs minimization solver only, this specifies the number of ‘outer’ iterations on T or P when some property pair other than TP is specified.
loglevel – Set to a value > 0 to write diagnostic output. Larger values generate more detailed information.
- gibbs_RT(tp)¶
Get the non-dimensional Gibbs function.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Vector of non-dimensional Gibbs functions of the species.
- gibbs_mass(tp)¶
Get the mass specific Gibbs function.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Mass specific Gibbs function of the mixture. Units: J/kg
- gibbs_mole(tp)¶
Get the mole specific Gibbs function.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Molar specific Gibbs function of the mixture. Units: J/kmol
- intEnergy_mass(tp)¶
Get the mass specific internal energy.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Mass specific internal energy of the mixture. Units: J/kg
- intEnergy_mole(tp)¶
Get the mole specific internal energy.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Molar specific internal energy of the mixture. Units: J/kmol
- isIdealGas(tp)¶
Get a flag indicating whether the phase is an ideal gas.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
True (1) if the phase is an ideal gas or ideal gas mixture, and false (0) otherwise.
- isothermalCompressibility(tp)¶
Get the isothermal compressibility.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Isothermal Compressibility. Units: 1/Pa
- massFraction(tp, species)¶
Get the mass fraction of a species.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)species – String or cell array of strings of species whose mass fraction is desired
- Returns
Scalar or vector double mass fractions
- massFractions(tp)¶
Get the mass fractions of all species.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Vector of species mass fractions for input phase. If no output argument is specified, a bar plot is produced.
- maxTemp(tp)¶
Get the maximum temperature of the parameter fits.
The parameterizations used to represent the temperature-dependent species thermodynamic properties are generally only valid in some finite temperature range, which may be different for each species in the phase.
See also:
minTemp()
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Vector of maximum temperatures of all species
- meanMolecularWeight(tp)¶
Get the mean molecular weight.
The mean molecular weight is the mole-fraction-weighted sum of the molar masses of the individual species in the phase.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Scalar double mean molecular weight. Units: kg/kmol
- minTemp(tp)¶
Get the minimum temperature of the parameter fits.
The parameterizations used to represent the temperature-dependent species thermodynamic properties are generally only valid in some finite temperature range, which may be different for each species in the phase.
See also:
maxTemp()
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Vector of minimum temperatures of all species
- molarDensity(tp)¶
Get the molar density.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Molar density. Units: kmol/m^3
- moleFraction(tp, species)¶
Get the mole fraction of a species.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)species – String or cell array of strings of species whose mole fraction is desired
- Returns
Scalar or vector double mole fractions
- moleFractions(tp)¶
Get the mole fractions of all species.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Vector of species mole fractions for input phase. If no output argument is specified, a bar plot is produced.
- molecularWeights(tp)¶
Get the molecular weights of the species.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Vector of species molecular weights. Units: kg/kmol
- nAtoms(tp, k, m)¶
Get the number of atoms of an element in a species.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)k – String species name or integer species number
m – String element name or integer element number
- Returns
Number of atoms of element
m
in speciesk
.
- nElements(tp)¶
Get the number of elements.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Number of elements in the phase.
- nSpecies(tp)¶
Get the number of species.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Number of species in the phase.
- name(tp)¶
Get the name of the phase.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
String name of the input phase
- pressure(tp)¶
Get the pressure.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Pressure. Units: Pa
- refPressure(tp)¶
Get the reference pressure.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)- Returns
Reference pressure in Pa. All standard-state thermodynamic properties are for this pressure.
- satPressure(tp, T)¶
Get the saturation pressure for a given temperature.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)T – Temperature Units: K
- Returns
Saturation pressure for temperature T. Units: Pa
- satTemperature(tp, p)¶
Get the saturation temperature for a given pressure.
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)p – Pressure. Units: Pa
- Returns
Saturation temperature for pressure p. Units: K
- set(tp, varargin)¶
Set properties of a phase.
The properties that may be set are
Temperature (T)
Density (Rho)
Volume (V)
Pressure (P)
Enthalpy (H)
Entropy (S)
Mole Fractions (X)
Mass Fractions (Y)
Vapor Fraction (Vapor)
Liquid Fraction (Liquid)
Either the full property name or the symbol may be specified. For the extensive properties (V,H,U,S), the values must be given per unit mass. H, U, and S must be set in conjunction with pressure (for H,S) or volume (for U,S). Either (specific) volume or density may be specified. Mole and mass fractions must be input as vectors (either row or column) with length equal to the number of species. Two properties may be specified in a single call to
set()
, plus one of mass fractions or mole fractions.Examples:
>> set(gas,'Temperature',600.0); >> set(gas,'T',600.0); >> set(gas,'T',600.0,'P',2*oneatm,'Y',massfracs); >> set(gas,'H',0.5*enthalpy_mass(gas),'P',pressure(gas)); >> set(gas,'S',entropy_mass(gas),'P',0.5*pressure(gas)); >> set(gas,'X',ones(nSpecies(gas),1)); >> set(gas,'T',500.0,'Vapor',0.8)
Alternatively, individual methods to set properties may be called (setTemperature, setMoleFractions, etc.)
See also:
setDensity()
,setMassFractions()
,setMoleFractions()
,setPressure()
,setState_HP()
,setState_Psat()
,setState_SP()
,setState_SV()
,setState_Tsat()
,setState_UV()
,setTemperature()
- Parameters
tp – Instance of class
ThermoPhase()
(or another object that derives from ThermoPhase)varargin – Comma separated list of
property, value
pairs to be set
- setDensity(tp, rho)¶
Set the density.
- Parameters
tp – Instance of class
ThermoPhase()
(or another class derived from ThermoPhase)rho – Density. Units: kg/m**3
- setElectricPotential(tp, phi)¶
Set the electric potential.
- Parameters
tp – Instance of class
ThermoPhase()
(or another class derived from ThermoPhase)phi – Electric potential. Units: V
- setMassFractions(tp, y, norm)¶
Set the species mass fractions.
Note that calling
setMassFractions()
leaves the temperature and density unchanged, and therefore the pressure changes if the new composition has a molar mass that is different than the old composition. If it is desired to change the composition and hold the pressure fixed, use methodset()
and specify the mass fractions and the pressure, or callsetPressure()
after callingsetMassFractions()
.- Parameters
tp – Instance of class
ThermoPhase()
(or another class derived from ThermoPhase)y – Vector of mass fractions whose length must be the same as the number of species. Alternatively, a string in the format
'SPEC:Y,SPEC2:Y2'
that specifies the mass fraction of specific species.norm – If
'nonorm'
is specified,y
will be normalized. This only works ify
is a vector, not a string. Since unnormalized mass fractions can lead to unphysical results,'nonorm'
should be used only in rare cases, such as computing partial derivatives with respect to a species mass fraction.
- setMoleFractions(tp, x, norm)¶
Set the species mole fractions.
Note that calling
setMoleFractions()
leaves the temperature and density unchanged, and therefore the pressure changes if the new composition has a molar mass that is different than the old composition. If it is desired to change the composition and hold the pressure fixed, use methodset()
and specify the mole fractions and the pressure, or callsetPressure()
after callingsetmoleFractions()
.- Parameters
tp – Instance of class
ThermoPhase()
(or another class derived from ThermoPhase)y – Vector of mole fractions whose length must be the same as the number of species. Alternatively, a string in the format
'SPEC:Y,SPEC2:Y2'
that specifies the mole fraction of specific species.norm – If
'nonorm'
is specified,y
will be normalized. This only works ify
is a vector, not a string. Since unnormalized mole fractions can lead to unphysical results,'nonorm'
should be used only in rare cases, such as computing partial derivatives with respect to a species mole fraction.
- setName(tp, name)¶
Set the name of the phase.
- Parameters
tp – Instance of class
ThermoPhase()
(or another class derived from ThermoPhase)name – String, name of the phase
- setPressure(tp, p)¶
Set the pressure.
The pressure is set by changing the density holding the temperature and chemical composition fixed.
- Parameters
tp – Instance of class
ThermoPhase()
(or another class derived from ThermoPhase)p – Pressure. Units: Pa
- setState_HP(tp, hp)¶
Set the specific enthalpy and pressure.
- Parameters
tp – Instance of class
ThermoPhase()
(or another class derived from ThermoPhase)hp – Vector of length 2 containing the desired values for the specific enthalpy (J/kg) and pressure (Pa).
- setState_PV(tp, pv)¶
Set the pressure and specific volume.
- Parameters
tp – Instance of class
ThermoPhase()
(or another class derived from ThermoPhase)pv – Vector of length 2 containing the desired values for the pressure (Pa) and specific volume (m^3/kg).
- setState_Psat(tp, px)¶
Set the fluid state using the given pressure and quality.
The fluid state will be set to a saturated liquid-vapor state using the input pressure and vapor fraction (quality) as the independent, intensive variables.
- Parameters
tp – Instance of class
ThermoPhase()
(or another class derived from ThermoPhase)px – Vector of length 2 containing the desired values for the pressure (Pa) and the vapor fraction
- setState_RP(tp, rp)¶
Set the density and pressure.
The density is set first, then the pressure is set by changing the temperature holding the density and chemical composition fixed.
- Parameters
tp – Instance of class
ThermoPhase()
(or another class derived from ThermoPhase)rp – Vector of length 2 containing the desired values for the density (kg/m^3) and pressure (Pa)
- setState_SH(tp, sh)¶
Set the specific entropy and specific enthalpy.
- Parameters
tp – Instance of class
ThermoPhase()
(or another class derived from ThermoPhase)sh – Vector of length 2 containing the desired values for the specific entropy (J/kg/K) and specific enthalpy (J/kg).
- setState_SP(tp, sp)¶
Set the specific entropy and pressure.
- Parameters
tp – Instance of class
ThermoPhase()
(or another class derived from ThermoPhase)sp – Vector of length 2 containing the desired values for the specific entropy (J/kg-K) and pressure (Pa).
- setState_ST(tp, st)¶
Set the specific entropy and temperature.
- Parameters
tp – Instance of class
ThermoPhase()
(or another class derived from ThermoPhase)st – Vector of length 2 containing the desired values for the specific entropy (J/kg-K) and temperature (K).
- setState_SV(tp, sv)¶
Set the specific entropy and specific volume.
- Parameters
tp – Instance of class
ThermoPhase()
(or another class derived from ThermoPhase)sv – Vector of length 2 containing the desired values for the specific entropy (J/kg-K) and specific volume (m**3/kg).
- setState_TH(tp, th)¶
Set the temperature and specific enthalpy.
- Parameters
tp – Instance of class
ThermoPhase()
(or another class derived from ThermoPhase)th – Vector of length 2 containing the desired values for the temperature (K) and specific enthalpy (J/kg).
- setState_TV(tp, tv)¶
Set the temperature and specific volume.
- Parameters
tp – Instance of class
ThermoPhase()
(or another class derived from ThermoPhase)tv – Vector of length 2 containing the desired values for the temperature (K) and specific volume (m^3/kg).
- setState_Tsat(tp, tx)¶
Set the fluid state using the given temperature and quality.
The fluid state will be set to a saturated liquid-vapor state using the input temperature and vapor fraction (quality) as the independent, intensive variables.
- Parameters
tp – Instance of class
ThermoPhase()
(or another class derived from ThermoPhase)tx – Vector of length 2 containing the desired values for the temperature (K) and the vapor fraction (quality)
- setState_UP(tp, up)¶
Set the specific internal energy and pressure.
- Parameters
tp – Instance of class
ThermoPhase()
(or another class derived from ThermoPhase)up – Vector of length 2 containing the desired values for the specific internal energy (J/kg) and pressure (Pa).
- setState_UV(tp, uv)¶
Set the specific internal energy and specific volume.
- Parameters
tp – Instance of class
ThermoPhase()
(or another class derived from ThermoPhase)uv – Vector of length 2 containing the desired values for the specific internal energy (J/kg) and specific volume (m**3/kg).
- setState_VH(tp, vh)¶
Set the specific volume and specific enthalpy.
- Parameters
tp – Instance of class
ThermoPhase()
(or another class derived from ThermoPhase)vh – Vector of length 2 containing the desired values for the specific volume (m^3/kg) and specific enthalpy (J/kg).
- setState_satLiquid(tp)¶
Set the fluid to the saturated liquid state at the current temperature.
- Parameters
tp – Instance of class
ThermoPhase()
(or another class derived from ThermoPhase)
- setState_satVapor(tp)¶
Set the fluid to the saturated vapor state at the current temperature.
- Parameters
tp – Instance of class
ThermoPhase()
(or another class derived from ThermoPhase)
- setTemperature(tp, t)¶
Set the temperature.
- Parameters
tp – Instance of class
ThermoPhase()
(or another class derived from ThermoPhase)t – Temperature. Units: K
- soundspeed(tp)¶
Get the speed of sound.
If the phase is an ideal gas, the speed of sound is calculated by:
\[c = \sqrt{\gamma * R * T} \]where \(\gamma\) is the ratio of specific heats, \(R\) is the specific gas constant, and \(T\) is the temperature. If the phase is not an ideal gas, the speed of sound is calculated by
\[c = \sqrt{\left(\frac{\partial p}{\partial \rho}\right)_s} \]where \(p\) is the pressure and \(\rho\) is the density, and the subscript \(s\) indicates constant entropy. This is approximated by slightly increasing the density at constant entropy and computing the change in pressure.
\[c = \sqrt{\frac{p_1 - p_0}{\rho_1-\rho_0}} \]- Parameters
tp – Instance of class
ThermoPhase()
(or another class derived from ThermoPhase)- Returns
The speed of sound. Units: m/s
- speciesIndex(tp, name)¶
Get the index of a species given the name.
The index is an integer assigned to each species in sequence as it is read in from the input file.
NOTE: In keeping with the conventions used by Matlab, this method returns 1 for the first species, 2 for the second, etc. In contrast, the corresponding method speciesIndex in the Cantera C++ and Python interfaces returns 0 for the first species, 1 for the second one, etc.
>> ich4 = speciesIndex(gas, 'CH4'); >> iho2 = speciesIndex(gas, 'HO2');
- Parameters
tp – Instance of class
ThermoPhase()
(or another class derived from ThermoPhase)name – If name is a single string, the return value will be a integer containing the corresponding index. If it is an cell array of strings, the output will be an array of the same shape containing the indices.
- Returns
Scalar or array of integers
- speciesName(tp, k)¶
Get the name of a species given the index.
- Parameters
tp – Instance of class
ThermoPhase()
(or another class derived from ThermoPhase)k – Scalar or array of integer species numbers
- Returns
Cell array of strings
- speciesNames(tp)¶
Get the species names.
- Parameters
tp – Instance of class
ThermoPhase()
(or another class derived from ThermoPhase)- Returns
Cell array of strings of all of the species names
- temperature(tp)¶
Get the temperature.
- Parameters
tp – Instance of class
ThermoPhase()
(or another class derived from ThermoPhase)- Returns
Temperature. Units: K
- thermalExpansionCoeff(tp)¶
Get the thermal expansion coefficient.
- Parameters
tp – Instance of class
ThermoPhase()
(or another class derived from ThermoPhase)- Returns
Thermal Expansion Coefficient. Units: 1/K
- vaporFraction(tp)¶
Get the vapor fraction.
- Parameters
tp – Instance of class
ThermoPhase()
(or another class derived from ThermoPhase)- Returns
Vapor fraction.