CoolProp.CoolProp module#

class CoolProp.CoolProp.AbstractState#

Bases: object

This class is a one-to-one python wrapper of the AbstractState class

Bvirial(self) double#

Get the B virial coefficient - wrapper of c++ function CoolProp::AbstractState::Bvirial(void)

Cvirial(self) double#

Get the C virial coefficient - wrapper of c++ function CoolProp::AbstractState::Cvirial(void)

PIP(self) double#

Get the phase identification parameter - wrapper of c++ function CoolProp::AbstractState::PIP()

Prandtl(self) double#

Get the Prandtl number - wrapper of c++ function CoolProp::AbstractState::Prandtl(void)

Q(self) double#

Get the vapor quality in mol/mol - wrapper of c++ function CoolProp::AbstractState::Q(void)

T(self) double#

Get the temperature in K - wrapper of c++ function CoolProp::AbstractState::T(void)

T_critical(self) double#

Gets the critical temperature in K - wrapper of c++ function CoolProp::AbstractState::T_critical()

T_reducing(self) double#

Gets the reducing temperature in K - wrapper of c++ function CoolProp::AbstractState::T_reducing()

Tmax(self) double#

Set the maximum temperature in K - wrapper of c++ function CoolProp::AbstractState::Tmax()

Tmin(self) double#

Set the minimum temperature in K- wrapper of c++ function CoolProp::AbstractState::Tmin()

Ttriple(self) double#

Set the triple point temperature in K - wrapper of c++ function CoolProp::AbstractState::Ttriple()

acentric_factor(self) double#

Get the acentric factor - wrapper of c++ function CoolProp::AbstractState::acentric_factor(void)

all_critical_points(self) list#

Calculate all the critical points - wrapper of c++ function CoolProp::AbstractState::all_critical_points()

alpha0(self) CoolPropDbl#

Get the ideal-gas reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::alpha0()

alphar(self) CoolPropDbl#

Get the residual reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::alphar()

apply_simple_mixing_rule(self, size_t i, size_t j, string model)#

Apply a simple mixing rule - wrapper of c++ function CoolProp::AbstractState::apply_simple_mixing_rule()

backend_name(self)#

Get the backend name - wrapper of c++ function CoolProp::AbstractState::backend_name()

build_phase_envelope(self, string type)#

Build the phase envelope - wrapper of c++ function CoolProp::AbstractState::build_phase_envelope()

build_spinodal(self)#

Calculate the spinodal - wrapper of c++ function CoolProp::AbstractState::build_spinodal()

change_EOS(self, size_t i, string EOS_name)#

Change the EOS for one component - wrapper of c++ function CoolProp::AbstractState::change_EOS()

chemical_potential(self, size_t i) double#

Get the chemical potential of the i-th component - wrapper of c++ function CoolProp::AbstractState::chemical_potential(std::size_t)

compressibility_factor(self) double#

Get the compressibility factor Z=p/(rho*R*T) - wrapper of c++ function CoolProp::AbstractState::compressibility_factor(void)

conductivity(self) double#

Get the thermal conductivity in W/m/K - wrapper of c++ function CoolProp::AbstractState::conductivity(void)

conductivity_contributions(self) dict#

Retrieve each of the contributions to the conductivity, each in W/m/K - wrapper of c++ function CoolProp::AbstractState::conductivity_contributions()

conformal_state(self, string reference_fluid, CoolPropDbl T, CoolPropDbl rho) dict#

Solve for conformal state used in extended corresponding states - wrapper of c++ function CoolProp::AbstractState::conformal_state()

cp0mass(self) double#

Get the ideal gas constant pressure specific heat in J/kg/K - wrapper of c++ function CoolProp::AbstractState::cp0mass(void)

cp0molar(self) double#

Get the ideal gas constant pressure specific heat in J/mol/K - wrapper of c++ function CoolProp::AbstractState::cp0molar(void)

cpmass(self) double#

Get the constant pressure specific heat in J/kg/K - wrapper of c++ function CoolProp::AbstractState::cpmass(void)

cpmolar(self) double#

Get the constant pressure specific heat in J/mol/K - wrapper of c++ function CoolProp::AbstractState::cpmolar(void)

criticality_contour_values(self) tuple#

Gets the criticality matrix values L1* and M1* - wrapper of c++ function CoolProp::AbstractState::criticality_contour_values() Returns a tuple of (L1*, M1*)

cvmass(self) double#

Get the constant volume specific heat in J/kg/K - wrapper of c++ function CoolProp::AbstractState::cvmass(void)

cvmolar(self) double#

Get the constant volume specific heat in J/mol/K - wrapper of c++ function CoolProp::AbstractState::cvmolar(void)

d2alpha0_dDelta2(self) CoolPropDbl#

Get the ideal-gas reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d2alpha0_dDelta2()

d2alpha0_dDelta_dTau(self) CoolPropDbl#

Get the ideal-gas reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d2alpha0_dDelta_dTau()

d2alpha0_dTau2(self) CoolPropDbl#

Get the ideal-gas reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d2alpha0_dTau2()

d2alphar_dDelta2(self) CoolPropDbl#

Get the residual reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d2alphar_dDelta2()

d2alphar_dDelta_dTau(self) CoolPropDbl#

Get the residual reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d2alphar_dDelta_dTau()

d2alphar_dTau2(self) CoolPropDbl#

Get the residual reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d2alphar_dTau2()

d3alpha0_dDelta2_dTau(self) CoolPropDbl#

Get the ideal-gas reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d3alpha0_dDelta2_dTau()

d3alpha0_dDelta3(self) CoolPropDbl#

Get the ideal-gas reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d3alpha0_dDelta3()

d3alpha0_dDelta_dTau2(self) CoolPropDbl#

Get the ideal-gas reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d3alpha0_dDelta_dTau2()

d3alpha0_dTau3(self) CoolPropDbl#

Get the ideal-gas reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d3alpha0_dTau3()

d3alphar_dDelta2_dTau(self) CoolPropDbl#

Get the residual reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d3alphar_dDelta2_dTau()

d3alphar_dDelta3(self) CoolPropDbl#

Get the residual reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d3alphar_dDelta3()

d3alphar_dDelta_dTau2(self) CoolPropDbl#

Get the residual reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d3alphar_dDelta_dTau2()

d3alphar_dTau3(self) CoolPropDbl#

Get the residual reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d3alphar_dTau3()

d4alphar_dDelta2_dTau2(self) CoolPropDbl#

Get the residual reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d4alphar_dDelta2_dTau2()

d4alphar_dDelta3_dTau(self) CoolPropDbl#

Get the residual reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d4alphar_dDelta3_dTau()

d4alphar_dDelta4(self) CoolPropDbl#

Get the residual reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d4alphar_dDelta4()

d4alphar_dDelta_dTau3(self) CoolPropDbl#

Get the residual reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d4alphar_dDelta_dTau3()

d4alphar_dTau4(self) CoolPropDbl#

Get the residual reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d4alphar_dTau4()

dalpha0_dDelta(self) CoolPropDbl#

Get the ideal-gas reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::dalpha0_dDelta()

dalpha0_dTau(self) CoolPropDbl#

Get the ideal-gas reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::dalpha0_dTau()

dalphar_dDelta(self) CoolPropDbl#

Get the residual reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::dalphar_dDelta()

dalphar_dTau(self) CoolPropDbl#

Get the residual reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::dalphar_dTau()

delta(self) double#

Get the reduced density - wrapper of c++ function CoolProp::AbstractState::delta(void)

first_partial_deriv(self, parameters OF, parameters WRT, parameters CONSTANT) CoolPropDbl#

Get the first partial derivative - wrapper of c++ function CoolProp::AbstractState::first_partial_deriv()

first_saturation_deriv(self, parameters OF, parameters WRT) CoolPropDbl#

Get the first derivative along the saturation curve - wrapper of c++ function CoolProp::AbstractState::first_saturation_deriv()

first_two_phase_deriv(self, parameters Of, parameters Wrt, parameters Constant) double#

Get the first two-phase derivative - wrapper of C++ function CoolProp::AbstractState::first_two_phase_deriv()

first_two_phase_deriv_splined(self, parameters Of, parameters Wrt, parameters Constant, double x_end) double#

Get the first two-phase derivative using splines - wrapper of C++ function CoolProp::AbstractState::first_two_phase_deriv_splined()

fluid_names(self)#

Get the list of fluid names - wrapper of c++ function CoolProp::AbstractState::fluid_names()

fluid_param_string(self, string key)#

Get a fluid parameter string - wrapper of c++ function CoolProp::AbstractState::fluid_param_string()

fugacity(self, size_t i) double#

Get the fugacity of the i-th component - wrapper of c++ function CoolProp::AbstractState::fugacity(std::size_t)

fugacity_coefficient(self, size_t i) double#

Get the fugacity coefficient of the i-th component - wrapper of c++ function CoolProp::AbstractState::fugacity_coefficient(std::size_t)

fundamental_derivative_of_gas_dynamics(self) double#

Get the fundamental derivative of gas dynamics - wrapper of c++ function CoolProp::AbstractState::fundamental_derivative_of_gas_dynamics(void)

gas_constant(self) double#

Get the gas constant in J/mol/K - wrapper of c++ function CoolProp::AbstractState::gas_constant(void)

get_binary_interaction_double#

Get a double precision interaction parameter - wrapper of c++ function CoolProp::AbstractState::get_binary_interaction_double()

get_binary_interaction_string(self, string CAS1, string CAS2, string parameter) string#

Get a string interaction parameter - wrapper of c++ function CoolProp::AbstractState::get_binary_interaction_string()

get_fluid_constant(self, size_t i, parameters param) double#

Get a constant for a fluid in the mixture CoolProp::AbstractState::get_fluid_constant()

get_fluid_parameter_double(self, size_t i, string parameter) double#

Get a fluid parameter that is a double-precision number - wrapper of c++ function CoolProp::AbstractState::get_fluid_parameter_double()

get_mass_fractions(self)#

Get the mass fractions - wrapper of c++ function CoolProp::AbstractState::get_mass_fractions()

get_mole_fractions(self)#

Get the mole fractions - wrapper of c++ function CoolProp::AbstractState::get_mole_fractions()

get_phase_envelope_data(self) PyPhaseEnvelopeData#

Get the phase envelope data - wrapper of c++ function CoolProp::AbstractState::get_phase_envelope_data()

get_spinodal_data(self) PySpinodalData#

Get the data from the spinodal - wrapper of c++ function CoolProp::AbstractState::get_spinodal_data()

gibbsmass(self) double#

Get the mass-specific Gibbs energy in J/kg - wrapper of c++ function CoolProp::AbstractState::gibbsmass(void)

gibbsmass_excess(self) double#

Get the mass-specific excess Gibbs energy in J/kg - wrapper of c++ function CoolProp::AbstractState::gibbsmass_excess(void)

gibbsmolar(self) double#

Get the mole-specific Gibbs energy in J/mol - wrapper of c++ function CoolProp::AbstractState::gibbsmolar(void)

gibbsmolar_excess(self) double#

Get the mole-specific excess Gibbs energy in J/mol - wrapper of c++ function CoolProp::AbstractState::gibbsmolar_excess(void)

gibbsmolar_residual(self) double#

Get the mole-specific residual Gibbs energy in J/mol - wrapper of c++ function CoolProp::AbstractState::gibbsmolar_residual(void)

has_melting_line(self) bool#

Check if the fluid has a melting line - True if is does, False otherwise - wrapper of c++ function CoolProp::AbstractState::has_melting_line()

helmholtzmass(self) double#

Get the mass-specific Helmholtz energy in J/kg - wrapper of c++ function CoolProp::AbstractState::helmholtzmass(void)

helmholtzmass_excess(self) double#

Get the mass-specific excess Helmholtz energy in J/kg - wrapper of c++ function CoolProp::AbstractState::helmholtzmass_excess(void)

helmholtzmolar(self) double#

Get the mole-specific Helmholtz energy in J/mol - wrapper of c++ function CoolProp::AbstractState::helmholtzmolar(void)

helmholtzmolar_excess(self) double#

Get the mole-specific excess Helmholtz energy in J/mol - wrapper of c++ function CoolProp::AbstractState::helmholtzmolar_excess(void)

hmass(self) double#

Get the enthalpy in J/kg - wrapper of c++ function CoolProp::AbstractState::hmass(void)

hmass_excess(self) double#

Get the mass-specific excess enthalpy in J/kg - wrapper of c++ function CoolProp::AbstractState::hmass_excess(void)

hmolar(self) double#

Get the enthalpy in J/mol - wrapper of c++ function CoolProp::AbstractState::hmolar(void)

hmolar_excess(self) double#

Get the mole-specific excess enthalpy in J/mol - wrapper of c++ function CoolProp::AbstractState::hmolar_excess(void)

hmolar_residual(self) double#

Get the mole-specific residual enthalpy in J/mol - wrapper of c++ function CoolProp::AbstractState::hmolar_residual(void)

ideal_curve(self, string type) tuple#

Get an ideal curve - wrapper of c++ function CoolProp::AbstractState::ideal_curve()

isobaric_expansion_coefficient(self) double#

Get the isobaric expansion coefficient - wrapper of c++ function CoolProp::AbstractState::isobaric_expansion_coefficient(void)

isothermal_compressibility(self) double#

Get the isothermal_compressibility - wrapper of c++ function CoolProp::AbstractState::isothermal_compressibility(void)

keyed_output(self, parameters iOutput) double#

Get a keyed output CoolProp::AbstractState::keyed_output(parameters key)

melting_line(self, int param, int given, double value) double#

Get values from the melting line - wrapper of c++ function CoolProp::AbstractState::melting_line()

molar_mass(self) double#

Get the molar mass in kg/mol - wrapper of c++ function CoolProp::AbstractState::molar_mass(void)

mole_fractions_liquid(self)#

Get the mole fractions of the liquid phase - wrapper of c++ function CoolProp::AbstractState::mole_fractions_liquid(void)

mole_fractions_vapor(self)#

Get the mole fractions of the vapor phase - wrapper of c++ function CoolProp::AbstractState::mole_fractions_vapor(void)

name(self)#

Get the fluid name - wrapper of c++ function CoolProp::AbstractState::name()

neff(self) double#

Get the effective hardness of interaction - wrapper of c++ function CoolProp::AbstractState::neff(void)

p(self) double#

Get the pressure in Pa - wrapper of c++ function CoolProp::AbstractState::p(void)

p_critical(self) double#

Gets the critical pressure in Pa - wrapper of c++ function CoolProp::AbstractState::p_critical()

phase(self) phases#

Get the phase as key value- wrapper of c++ function CoolProp::AbstractState::phase()

pmax(self) double#

Set the maximum pressure in Pa - wrapper of c++ function CoolProp::AbstractState::pmax()

rhomass(self) double#

Get the density in kg/m^3 - wrapper of c++ function CoolProp::AbstractState::rhomass(void)

rhomass_critical(self) double#

Gets the critical density in kg/m^3 - wrapper of c++ function CoolProp::AbstractState::rhomass_critical()

rhomass_reducing(self) double#

Gets the reducing density in kg/m^3 - wrapper of c++ function CoolProp::AbstractState::rhomass_reducing()

rhomolar(self) double#

Get the density in mol/m^3 - wrapper of c++ function CoolProp::AbstractState::rhomolar(void)

rhomolar_critical(self) double#

Gets the critical density in mol/m^3 - wrapper of c++ function CoolProp::AbstractState::rhomolar_critical()

rhomolar_reducing(self) double#

Gets the reducing density in mol/m^3 - wrapper of c++ function CoolProp::AbstractState::rhomolar_reducing()

saturated_liquid_keyed_output(self, parameters iOutput) double#

Get a trivial output for the saturated liquid CoolProp::AbstractState::saturated_liquid_keyed_output(parameters key)

saturated_vapor_keyed_output(self, parameters iOutput) double#

Get a trivial output for the saturated vapor CoolProp::AbstractState::saturated_vapor_keyed_output(parameters key)

saturation_ancillary(self, parameters param, int Q, parameters given, double value) double#

Get values from the saturation_ancillary - wrapper of c++ function CoolProp::AbstractState::saturation_ancillary()

second_partial_deriv(self, parameters OF, parameters WRT1, parameters CONSTANT1, parameters WRT2, parameters CONSTANT2) CoolPropDbl#

Get the second partial derivative - wrapper of c++ function CoolProp::AbstractState::second_partial_deriv()

second_saturation_deriv(self, parameters OF1, parameters WRT1, parameters WRT2) CoolPropDbl#

Get the second derivative along the saturation curve - wrapper of c++ function CoolProp::AbstractState::second_saturation_deriv()

second_two_phase_deriv(self, parameters Of1, parameters Wrt1, parameters Constant1, parameters Wrt2, parameters Constant2) double#

Get the second two-phase derivative - wrapper of C++ function CoolProp::AbstractState::second_two_phase_deriv()

set_binary_interaction_double#

Set a double precision interaction parameter - wrapper of c++ function CoolProp::AbstractState::set_binary_interaction_double()

set_binary_interaction_string#

Set a string interaction parameter - wrapper of c++ function CoolProp::AbstractState::set_binary_interaction_string()

set_cubic_alpha_C(self, size_t i, string parameter, double c1, double c2, double c3)#

Set alpha function (MC or Twu) and fluid specific parameters - wrapper of c++ function CoolProp::AbstractCubicBackend::set_cubic_alpha_C()

set_fluid_parameter_double(self, size_t i, string parameter, double val)#

Set a fluid parameter that is a double-precision number - wrapper of c++ function CoolProp::AbstractState::set_fluid_parameter_double()

set_mass_fractions(self, vector[double] z)#

Set the mass fractions - wrapper of c++ function CoolProp::AbstractState::set_mass_fractions()

set_mole_fractions(self, vector[double] z)#

Set the mole fractions - wrapper of c++ function CoolProp::AbstractState::set_mole_fractions()

set_volu_fractions(self, vector[double] z)#

Set the volume fractions - wrapper of c++ function CoolProp::AbstractState::set_volu_fractions()

smass(self) double#

Get the entropy in J/kg/K - wrapper of c++ function CoolProp::AbstractState::smass(void)

smass_excess(self) double#

Get the mass-specific excess entropy in J/kg/K - wrapper of c++ function CoolProp::AbstractState::smass_excess(void)

smolar(self) double#

Get the entropy in J/mol/K - wrapper of c++ function CoolProp::AbstractState::smolar(void)

smolar_excess(self) double#

Get the mole-specific excess entropy in J/mol/K - wrapper of c++ function CoolProp::AbstractState::smolar_excess(void)

smolar_residual(self) double#

Get the mole-specific residual entropy in J/mol/K - wrapper of c++ function CoolProp::AbstractState::smolar_residual(void)

specify_phase(self, phases phase)#

Specify the phase - wrapper of c++ function CoolProp::AbstractState::specify_phase()

speed_sound(self) double#

Get the speed of sound in m/s - wrapper of c++ function CoolProp::AbstractState::speed_sound(void)

surface_tension(self) double#

Get the surface tension N/m - wrapper of c++ function CoolProp::AbstractState::surface_tension(void)

tangent_plane_distance(self, double T, double p, vector[double] w, double rhomolar_guess=-1) double#

Gets the tangent_plane_distance - wrapper of c++ function CoolProp::AbstractState::tangent_plane_distance()

tau(self) double#

Get the reciprocal reduced temperature - wrapper of c++ function CoolProp::AbstractState::tau(void)

trivial_keyed_output(self, parameters iOutput) double#

Get a trivial keyed output not requiring any iteration CoolProp::AbstractState::trivial_keyed_output(parameters key)

true_critical_point(self) tuple#

Get the “true” critical point where dp/drho|T = 0 & d2p/drho^2|T = 0 - wrapper of c++ function CoolProp::AbstractState::true_critical_point()

umass(self) double#

Get the internal energy in J/kg - wrapper of c++ function CoolProp::AbstractState::umass(void)

umass_excess(self) double#

Get the mass-specific excess internal energy in J/kg - wrapper of c++ function CoolProp::AbstractState::umass_excess(void)

umolar(self) double#

Get the internal energy in J/mol - wrapper of c++ function CoolProp::AbstractState::umolar(void)

umolar_excess(self) double#

Get the mole-specific excess internal energy in J/mol - wrapper of c++ function CoolProp::AbstractState::umolar_excess(void)

unspecify_phase(self)#

Unspecify the phase - wrapper of c++ function CoolProp::AbstractState::unspecify_phase()

update(self, input_pairs ipair, double Value1, double Value2)#

Update function - wrapper of c++ function CoolProp::AbstractState::update()

update_QT_pure_superanc(self, double Q, double T)#

Update function - wrapper of c++ function CoolProp::AbstractState::update()

update_with_guesses(self, input_pairs ipair, double Value1, double Value2, PyGuessesStructure guesses)#

Update function - wrapper of c++ function CoolProp::AbstractState::update()

viscosity(self) double#

Get the viscosity in Pa-s - wrapper of c++ function CoolProp::AbstractState::viscosity(void)

viscosity_contributions(self) dict#

Retrieve each of the contributions to the viscosity, each in Pa-s - wrapper of c++ function CoolProp::AbstractState::viscosity_contributions()

volumemass_excess(self) double#

Get the mass-specific excess volume in m^3/kg - wrapper of c++ function CoolProp::AbstractState::volumemass_excess(void)

volumemolar_excess(self) double#

Get the mole-specific excess volume in m^3/mol - wrapper of c++ function CoolProp::AbstractState::volumemolar_excess(void)

class CoolProp.CoolProp.ChebyshevApproximation1D#

Bases: object

count_x_for_y_many(self, double[::1] y, unsigned int bits, size_t max_iter, double boundstytol, double[::1] counts)#
eval_many(self, double[::1] x, double[::1] y)#
get_x_for_y(self, double y, unsigned int bits, size_t max_iter, double boundstytol)#
is_monotonic(self)#
monotonic_intervals(self)#
xmax(self)#
xmin(self)#
class CoolProp.CoolProp.ChebyshevExpansion#

Bases: object

coeff(self)#
eval_many(self, double[::1] x, double[::1] y)#
solve_for_x(self, double y, double a, double b, unsigned int bits, size_t max_iter, double boundstytol)#
solve_for_x_many(self, double[::1] y, double a, double b, unsigned int bits, size_t max_iter, double boundstytol, double[::1] x, double[::1] counts)#
xmax(self)#
xmin(self)#
CoolProp.CoolProp.FluidsList() list#

Return a list of strings of all fluid names

Returns:

FluidsList – All the fluids that are included in CoolProp

Return type:

list of strings of fluid names

Notes

Here is an example:

In [0]: from CoolProp.CoolProp import FluidsList

In [1]: FluidsList()
CoolProp.CoolProp.HAProps(string OutputName, string Input1Name, Input1, string Input2Name, Input2, string Input3Name, Input3)#

DEPRECATED!!! Used the HAPropsSI function

CoolProp.CoolProp.HAPropsSI(string OutputName, string Input1Name, Input1, string Input2Name, Input2, string Input3Name, Input3)#

Copyright Ian Bell, 2011 email: ian.h.bell@gmail.com

The function is called like

HAPropsSI(‘H’,’T’,298.15,’P’,101325,’R’,0.5)

which will return the enthalpy of the air for a set of inputs of dry bulb temperature of 25C, atmospheric pressure, and a relative humidity of 50%.

This function implements humid air properties based on the analysis in ASHRAE RP-1845 which is available online: http://rp.ashrae.biz/page/ASHRAE-D-RP-1485-20091216.pdf

It employs real gas properties for both air and water, as well as the most accurate interaction parameters and enhancement factors. The IAPWS-95 formulation for the properties of water is used throughout in preference to the industrial formulation. It is unclear why the industrial formulation is used in the first place.

Since humid air is nominally a binary mixture, three variables are needed to fix the state. At least one of the input parameters must be dry-bulb temperature, relative humidity, dew-point temperature, or humidity ratio. The others will be calculated. If the output variable is a transport property (conductivity or viscosity), the state must be able to be found directly - i.e. make sure you give temperature and relative humidity or humidity ratio. The list of possible input variables are

String

Aliases

Description

T

Tdb

Dry-Bulb Temperature [K]

B

Twb

Wet-Bulb Temperature [K]

D

Tdp

Dew-Point Temperature [K]

P

Pressure [Pa]

V

Vda

Mixture volume [m3/kg dry air]

R

RH

Relative humidity in (0,1) [-]

W

Omega

Humidity Ratio [kg water/kg dry air]

H

Hda

Mixture enthalpy [J/kg dry air]

S

Sda

Mixture entropy [J/kg dry air/K]

C

cp

Mixture specific heat [J/kg dry air/K]

M

Visc

Mixture viscosity [Pa-s]

K

Mixture thermal conductivity [W/m/K]

There are also strings for the mixture volume and mixture enthalpy that will return the properties on a total humid air flow rate basis, they are given by ‘Vha’ [units of m^3/kg humid air] and ‘Cha’ [units of kJ/kg humid air/K] and ‘Hha’ [units of kJ/kg humid air] respectively.

For more information, go to http://www.coolprop.org

CoolProp.CoolProp.HAProps_Aux(str OutputName, double T, double p, double w) tuple#

Allows low-level access to some of the routines employed in HumidAirProps

Returns tuples of the form (Value, Units) where Value is the actual value and Units is a string that describes the units

The list of possible inputs is

  • Baa [First virial air-air coefficient]

  • Caaa [Second virial air coefficient]

  • Bww [First virial water-water coefficient]

  • Cwww [Second virial water coefficient]

  • Baw [First cross virial coefficient]

  • Caww [Second air-water-water virial coefficient]

  • Caaw [Second air-air-water virial coefficient]

  • beta_H

  • kT

  • vbar_ws [Molar saturated volume of water vapor]

  • p_ws [Saturated vapor pressure of pure water (>=0.01C) or ice (<0.01 C)]

  • f [Enhancement factor]

CoolProp.CoolProp.PhaseSI(in1, in2, in3, in4, in5)#

A Python wrapper of C++ function CoolProp::PhaseSI()

Does not support vectorization of the inputs like PropsSI

CoolProp.CoolProp.Props(in1, in2, in3=None, in4=None, in5=None, in6=None)#

This function is deprecated, use PropsSI instead

CoolProp.CoolProp.PropsSI(in1, in2, in3=None, in4=None, in5=None, in6=None, in7=None)#

A Python wrapper of C++ function CoolProp::PropsSI() .

class CoolProp.CoolProp.PyCriticalState#

Bases: object

T#

‘double’

Type:

T

hmolar#

‘double’

Type:

hmolar

p#

‘double’

Type:

p

rhomolar#

‘double’

Type:

rhomolar

smolar#

‘double’

Type:

smolar

stable#

‘bool’

Type:

stable

class CoolProp.CoolProp.PyGuessesStructure#

Bases: object

T#

‘double’

Type:

T

hmolar#

‘double’

Type:

hmolar

p#

‘double’

Type:

p

rhomolar#

‘double’

Type:

rhomolar

rhomolar_liq#

‘double’

Type:

rhomolar_liq

rhomolar_vap#

‘double’

Type:

rhomolar_vap

smolar#

‘double’

Type:

smolar

x#

list

Type:

x

y#

list

Type:

y

class CoolProp.CoolProp.PyPhaseEnvelopeData#

Bases: object

K#

list

Type:

K

Q#

list

Type:

Q

T#

list

Type:

T

TypeI#

‘bool’

Type:

TypeI

hmolar_liq#

list

Type:

hmolar_liq

hmolar_vap#

list

Type:

hmolar_vap

iTsat_max#

‘size_t’

Type:

iTsat_max

icrit#

‘size_t’

Type:

icrit

ipsat_max#

‘size_t’

Type:

ipsat_max

lnT#

list

Type:

lnT

lnp#

list

Type:

lnp

lnrhomolar_liq#

list

Type:

lnrhomolar_liq

lnrhomolar_vap#

list

Type:

lnrhomolar_vap

p#

list

Type:

p

rhomolar_liq#

list

Type:

rhomolar_liq

rhomolar_vap#

list

Type:

rhomolar_vap

smolar_liq#

list

Type:

smolar_liq

smolar_vap#

list

Type:

smolar_vap

x#

list

Type:

x

y#

list

Type:

y

class CoolProp.CoolProp.PySpinodalData#

Bases: object

M1#

‘vector[double]’

Type:

M1

delta#

‘vector[double]’

Type:

delta

tau#

‘vector[double]’

Type:

tau

class CoolProp.CoolProp.State(_Fluid, dict StateDict, phase=None, backend=None)#

Bases: object

A class that contains all the code that represents a thermodynamic state

Warning

This class is deprecated. You should use CoolProp.AbstractState instead

The motivation for this class is that it is useful to be able to define the state once using whatever state inputs you like and then be able to calculate other thermodynamic properties with the minimum of computational work.

Let’s suppose that you have inputs of pressure and temperature and you want to calculate the enthalpy and pressure. Since the Equations of State are all explicit in temperature and density, each time you call something like:

h = PropsSI('H','T',T','P',P,Fluid)
s = PropsSI('S','T',T','P',P,Fluid)

the solver is used to carry out the T-P flash calculation. And if you wanted entropy as well you could either intermediately calculate T, rho and then use T, rho in the EOS in a manner like:

rho = PropsSI('D','T',T','P',P,Fluid)
h = PropsSI('H','T',T','D',rho,Fluid)
s = PropsSI('S','T',T','D',rho,Fluid)

Instead in this class all that is handled internally. So the call to update sets the internal variables in the most computationally efficient way possible

Parameters:
  • Fluid (string)

  • StateDict (dictionary) – The state of the fluid - passed to the update function; if None, does not do a state update

  • phase (string) – DEPRECATED : this input is ignored

  • backend (string) – The CoolProp backend that should be used, one of “HEOS” (default), “REFPROP”, “INCOMP”, “BRINE”, etc.

Fluid#
MM#

The molar mass [kg/kmol] or [g/mol]

Phase(self) long#

Returns an integer flag for the phase of the fluid, where the flag value is one of iLiquid, iSupercritical, iGas, iTwoPhase

These constants are defined in the phase_constants module, and are imported into this module

Prandtl#

The Prandtl number (cp*mu/k) [-]

Props(self, parameters iOutput) double#
Q#

The quality [-]

T#

The temperature [K]

Tsat#

The saturation temperature (dew) for the given pressure, in [K]

copy(self) State#

Make a copy of this State class

cp#

The specific heat at constant pressure [kJ/kg/K]

cp0#

The ideal-gas specific heat at constant pressure [kJ/kg/K]

cv#

The specific heat at constant volume [kJ/kg/K]

dpdT#
get_MM(self) double#

Get the mole mass [kg/kmol] or [g/mol]

get_Q(self) double#

Get the quality [-]

get_T(self) double#

Get the temperature [K]

get_Tsat(self, double Q=1)#

Get the saturation temperature, in [K]

Returns None if pressure is not within the two-phase pressure range

get_cond(self) double#

Get the thermal conductivity, in [kW/m/K]

get_cp(self) double#

Get the specific heat at constant pressure [kJ/kg/K]

get_cp0(self) double#

Get the specific heat at constant pressure for the ideal gas [kJ/kg/K]

get_cv(self) double#

Get the specific heat at constant volume [kJ/kg/K]

get_dpdT(self) double#
get_h(self) double#

Get the specific enthalpy [kJ/kg]

get_p(self) double#

Get the pressure [kPa]

get_rho(self) double#

Get the density [kg/m^3]

get_s(self) double#

Get the specific enthalpy [kJ/kg/K]

get_speed_sound(self) double#

Get the speed of sound [m/s]

get_subcooling(self)#

Get the amount of subcooling below the saturation temperature corresponding to the pressure, in [K]

Returns None if pressure is not within the two-phase pressure range

get_superheat(self)#

Get the amount of superheat above the saturation temperature corresponding to the pressure, in [K]

Returns None if pressure is not within the two-phase pressure range

get_u(self) double#

Get the specific internal energy [kJ/kg]

get_visc(self) double#

Get the viscosity, in [Pa-s]

h#

The specific enthalpy [kJ/kg]

k#

The thermal conductivity, in [kW/m/K]

p#

The pressure [kPa]

pAS#

CoolProp.CoolProp.AbstractState

Type:

pAS

phase#
rho#

The density [kg/m^3]

s#

The specific enthalpy [kJ/kg/K]

set_Fluid(self, Fluid, backend)#
speed_test(self, int N)#
subcooling#

The amount of subcooling below the saturation temperature corresponding to the pressure, in [K]

Returns None if pressure is not within the two-phase pressure range

superheat#

The amount of superheat above the saturation temperature corresponding to the pressure, in [K]

Returns None if pressure is not within the two-phase pressure range

u#

The internal energy [kJ/kg]

update(self, dict params)#

Parameters params, dictionary

A dictionary of terms to be updated, with keys equal to single-char inputs to the Props function, for instance dict(T=298, P = 101.325) would be one standard atmosphere

update_Trho(self, double T, double rho)#

Just use the temperature and density directly for speed

Parameters:
  • T (float) – Temperature [K]

  • rho (float) – Density [kg/m^3]

update_ph(self, double p, double h)#

Use the pressure and enthalpy directly

Parameters:
  • p (float) – Pressure (absolute) [kPa]

  • h (float) – Enthalpy [kJ/kg]

visc#

The viscosity, in [Pa-s]

class CoolProp.CoolProp.SuperAncillary#

Bases: object

eval_sat(self, double T, str prop, short Q)#
eval_sat_many(self, double[::1] T, str prop, short Q, double[::1] y)#
CoolProp.CoolProp.add_fluids_as_JSON(backend, JSONstring)#

Add fluids in a JSON-formatted string format. Python wrapper of C++ function CoolProp::add_fluids_as_JSON()

CoolProp.CoolProp.apply_simple_mixing_rule(CAS1, CAS2, rule)#

Apply simple mixing rule. Currently linear or Lorentz-Berthelot. Python wrapper of C++ function CoolProp::apply_simple_mixing_rule()

CoolProp.CoolProp.cair_sat(double T) double#

The derivative of the saturation enthalpy cair_sat = d(hsat)/dT

CoolProp.CoolProp.config_key_description(string key) string#

Obtain the string description for a configuration key. Python wrapper of C++ function CoolProp::config_key_description()

CoolProp.CoolProp.extract_backend(in_str)#

A Python wrapper of C++ function CoolProp::extract_backend() .

CoolProp.CoolProp.extract_fractions(flds)#

A Python wrapper of C++ function CoolProp::extract_fractions() .

CoolProp.CoolProp.generate_update_pair(parameters key1, double value1, parameters key2, double value2) tuple#

This function will generate an input pair to the update() function given the key, value pairs for both inputs

CoolProp.CoolProp.get_BibTeXKey(Fluid, key) string#

Return the BibTeX key for the given fluid.

The possible keys are

  • EOS

  • CP0

  • VISCOSITY

  • CONDUCTIVITY

  • ECS_LENNARD_JONES

  • ECS_FITS

  • SURFACE_TENSION

  • MELTING_LINE

BibTeX keys refer to the BibTeX file in the trunk/CoolProp folder

Returns:

empty string if Fluid not in CoolProp, “Bad key” if key is invalid

Return type:

key, string

CoolProp.CoolProp.get_REFPROPname(Fluid) string#

Return the REFPROP compatible name for the fluid

Some fluids do not use the REFPROP name. For instance, ammonia is R717, and propane is R290. You can still can still call CoolProp using the name ammonia or R717, but REFPROP requires that you use a limited subset of names. Therefore, this function that returns the REFPROP compatible name. To then use this to call REFPROP, you would do something like:

In [0]: from CoolProp.CoolProp import get_REFPROPname, PropsSI

In [1]: get_REFPROPname('R290')

In [2]: PropsSI('D', 'T', 300, 'P', 300, Fluid)
CoolProp.CoolProp.get_aliases(Fluid)#

Return a comma separated string of aliases for the given fluid

CoolProp.CoolProp.get_config_as_json_string() string#

Obtain a json formulation of the internal configuration in CoolProp

Values can be set by passing a modified json library (converted to string) to set_config_as_json_string

CoolProp.CoolProp.get_config_bool(configuration_keys key) bool#

Get a configuration key that is a boolean; wrapper of wrapper of C++ function CoolProp::get_config_bool()

CoolProp.CoolProp.get_config_double(configuration_keys key) double#

Get a configuration key that is a double-precision float; wrapper of wrapper of C++ function CoolProp::get_config_double()

CoolProp.CoolProp.get_config_int(configuration_keys key) int#

Get a configuration key that is an integer; wrapper of wrapper of C++ function CoolProp::get_config_int()

CoolProp.CoolProp.get_config_string(configuration_keys key) string#

Get a configuration key that is a string; wrapper of wrapper of C++ function CoolProp::get_config_string()

CoolProp.CoolProp.get_debug_level() int#

Return the current debug level as integer

Returns:

level – If level is 0, no output will be written to screen, if >0, some output will be written to screen. The larger level is, the more verbose the output will be

Return type:

int

CoolProp.CoolProp.get_errstr() string#

Return the current error string

CoolProp.CoolProp.get_fluid_param_string(fluid, param)#
CoolProp.CoolProp.get_global_param_string(param)#
CoolProp.CoolProp.get_mixture_binary_pair_data(CAS1, CAS2, key) string#

Obtain mixture interaction parameter. Python wrapper of C++ function CoolProp::get_mixture_binary_pair_data()

CoolProp.CoolProp.get_mixture_binary_pair_pcsaft(CAS1, CAS2, key) string#

Obtain mixture PC-SAFT interaction parameter. Python wrapper of C++ function CoolProp::get_mixture_binary_pair_pcsaft()

CoolProp.CoolProp.get_parameter_index(string key) int#
CoolProp.CoolProp.get_parameter_information(int key, string info) string#
CoolProp.CoolProp.get_phase_index(string key) int#
CoolProp.CoolProp.is_trivial_parameter(int key)#
CoolProp.CoolProp.rebuildState(d)#
CoolProp.CoolProp.saturation_ancillary(string name, string output, int Q, string input, double value) double#

Return a value from the saturation ancillary equations; python wrapper of CoolProp::saturation_ancillary()

CoolProp.CoolProp.set_config_as_json_string(string s)#

Set the internal configuration in CoolProp from a json data string

Current state can be obtained by calling get_config_as_json_string

CoolProp.CoolProp.set_config_bool(configuration_keys key, bool value)#

Set a configuration key that is a boolean; wrapper of wrapper of C++ function CoolProp::set_config_bool()

CoolProp.CoolProp.set_config_double(configuration_keys key, double value)#

Set configuration key that is a double-precision float; wrapper of wrapper of C++ function CoolProp::set_config_double()

CoolProp.CoolProp.set_config_int(configuration_keys key, int value)#

Set a configuration key that is an integer; wrapper of wrapper of C++ function CoolProp::set_config_int()

CoolProp.CoolProp.set_config_string(configuration_keys key, string value)#

Set a configuration key that is a string; wrapper of wrapper of C++ function CoolProp::set_config_string()

CoolProp.CoolProp.set_debug_level(int level)#

Set the current debug level as integer in the range [0,10]

Parameters:

level (int) – If level is 0, no output will be written to screen, if >0, some output will be written to screen. The larger level is, the more verbose the output will be

CoolProp.CoolProp.set_departure_functions(functions)#

Specify the departure terms as JSON. Python wrapper of C++ function CoolProp::set_departure_functions()

CoolProp.CoolProp.set_interaction_parameters(data)#

Specify the binary interaction terms as JSON. Python wrapper of C++ function CoolProp::set_interaction_parameters()

CoolProp.CoolProp.set_mixture_binary_pair_data(CAS1, CAS2, key, val)#

Set mixture interaction parameter. Python wrapper of C++ function CoolProp::set_mixture_binary_pair_data()

CoolProp.CoolProp.set_mixture_binary_pair_pcsaft(CAS1, CAS2, key, val)#

Set mixture PC-SAFT interaction parameter. Python wrapper of C++ function CoolProp::set_mixture_binary_pair_pcsaft()

CoolProp.CoolProp.set_predefined_mixtures(data)#

Specify predefined mixtures as JSON. Python wrapper of C++ function CoolProp::set_predefined_mixtures()

CoolProp.CoolProp.set_reference_state(FluidName, *args)#

Accepts one of two signatures:

Type #1 (A Python wrapper of CoolProp::set_reference_stateS()):

set_reference_state(FluidName,reference_state)

FluidName The name of the fluid param reference_state The reference state to use, one of

IIR

(h=200 kJ/kg, s=1 kJ/kg/K at 0C sat. liq.)

ASHRAE

(h=0,s=0 @ -40C sat liq)

NBP

(h=0,s=0 @ 1.0 bar sat liq.)

Type #2 (A Python wrapper of CoolProp::set_reference_stateD()):

set_reference_state(FluidName,T0,rhomolar,hmolar0,smolar0)

Note

Only supported for internal backend currently

FluidName The name of the fluid

T0 The temperature at the reference point [K]

rhomolar The density at the reference point [mol/m^3]

hmolar0 The enthalpy at the reference point [J/mol]

smolar0 The entropy at the reference point [J/mol/K]