CoolProp.AbstractState module¶
This class is a one-to-one python wrapper of the AbstractState class
- CoolProp.AbstractState.Bvirial(self) double ¶
Get the B virial coefficient - wrapper of c++ function CoolProp::AbstractState::Bvirial(void)
- CoolProp.AbstractState.Cvirial(self) double ¶
Get the C virial coefficient - wrapper of c++ function CoolProp::AbstractState::Cvirial(void)
- CoolProp.AbstractState.PIP(self) double ¶
Get the phase identification parameter - wrapper of c++ function CoolProp::AbstractState::PIP()
- CoolProp.AbstractState.Prandtl(self) double ¶
Get the Prandtl number - wrapper of c++ function CoolProp::AbstractState::Prandtl(void)
- CoolProp.AbstractState.Q(self) double ¶
Get the vapor quality in mol/mol - wrapper of c++ function CoolProp::AbstractState::Q(void)
- CoolProp.AbstractState.T(self) double ¶
Get the temperature in K - wrapper of c++ function CoolProp::AbstractState::T(void)
- CoolProp.AbstractState.T_critical(self) double ¶
Gets the critical temperature in K - wrapper of c++ function CoolProp::AbstractState::T_critical()
- CoolProp.AbstractState.T_reducing(self) double ¶
Gets the reducing temperature in K - wrapper of c++ function CoolProp::AbstractState::T_reducing()
- CoolProp.AbstractState.Tmax(self) double ¶
Set the maximum temperature in K - wrapper of c++ function CoolProp::AbstractState::Tmax()
- CoolProp.AbstractState.Tmin(self) double ¶
Set the minimum temperature in K- wrapper of c++ function CoolProp::AbstractState::Tmin()
- CoolProp.AbstractState.Ttriple(self) double ¶
Set the triple point temperature in K - wrapper of c++ function CoolProp::AbstractState::Ttriple()
- CoolProp.AbstractState.acentric_factor(self) double ¶
Get the acentric factor - wrapper of c++ function CoolProp::AbstractState::acentric_factor(void)
- CoolProp.AbstractState.all_critical_points(self) list ¶
Calculate all the critical points - wrapper of c++ function CoolProp::AbstractState::all_critical_points()
- CoolProp.AbstractState.alpha0(self) CoolPropDbl ¶
Get the ideal-gas reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::alpha0()
- CoolProp.AbstractState.alphar(self) CoolPropDbl ¶
Get the residual reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::alphar()
- CoolProp.AbstractState.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()
- CoolProp.AbstractState.backend_name(self)¶
Get the backend name - wrapper of c++ function CoolProp::AbstractState::backend_name()
- CoolProp.AbstractState.build_phase_envelope(self, string type)¶
Build the phase envelope - wrapper of c++ function CoolProp::AbstractState::build_phase_envelope()
- CoolProp.AbstractState.build_spinodal(self)¶
Calculate the spinodal - wrapper of c++ function CoolProp::AbstractState::build_spinodal()
- CoolProp.AbstractState.change_EOS(self, size_t i, string EOS_name)¶
Change the EOS for one component - wrapper of c++ function CoolProp::AbstractState::change_EOS()
- CoolProp.AbstractState.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)
- CoolProp.AbstractState.compressibility_factor(self) double ¶
Get the compressibility factor Z=p/(rho*R*T) - wrapper of c++ function CoolProp::AbstractState::compressibility_factor(void)
- CoolProp.AbstractState.conductivity(self) double ¶
Get the thermal conductivity in W/m/K - wrapper of c++ function CoolProp::AbstractState::conductivity(void)
- CoolProp.AbstractState.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()
- CoolProp.AbstractState.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()
- CoolProp.AbstractState.cp0mass(self) double ¶
Get the ideal gas constant pressure specific heat in J/kg/K - wrapper of c++ function CoolProp::AbstractState::cp0mass(void)
- CoolProp.AbstractState.cp0molar(self) double ¶
Get the ideal gas constant pressure specific heat in J/mol/K - wrapper of c++ function CoolProp::AbstractState::cp0molar(void)
- CoolProp.AbstractState.cpmass(self) double ¶
Get the constant pressure specific heat in J/kg/K - wrapper of c++ function CoolProp::AbstractState::cpmass(void)
- CoolProp.AbstractState.cpmolar(self) double ¶
Get the constant pressure specific heat in J/mol/K - wrapper of c++ function CoolProp::AbstractState::cpmolar(void)
- CoolProp.AbstractState.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*)
- CoolProp.AbstractState.cvmass(self) double ¶
Get the constant volume specific heat in J/kg/K - wrapper of c++ function CoolProp::AbstractState::cvmass(void)
- CoolProp.AbstractState.cvmolar(self) double ¶
Get the constant volume specific heat in J/mol/K - wrapper of c++ function CoolProp::AbstractState::cvmolar(void)
- CoolProp.AbstractState.d2alpha0_dDelta2(self) CoolPropDbl ¶
Get the ideal-gas reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d2alpha0_dDelta2()
- CoolProp.AbstractState.d2alpha0_dDelta_dTau(self) CoolPropDbl ¶
Get the ideal-gas reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d2alpha0_dDelta_dTau()
- CoolProp.AbstractState.d2alpha0_dTau2(self) CoolPropDbl ¶
Get the ideal-gas reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d2alpha0_dTau2()
- CoolProp.AbstractState.d2alphar_dDelta2(self) CoolPropDbl ¶
Get the residual reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d2alphar_dDelta2()
- CoolProp.AbstractState.d2alphar_dDelta_dTau(self) CoolPropDbl ¶
Get the residual reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d2alphar_dDelta_dTau()
- CoolProp.AbstractState.d2alphar_dTau2(self) CoolPropDbl ¶
Get the residual reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d2alphar_dTau2()
- CoolProp.AbstractState.d3alpha0_dDelta2_dTau(self) CoolPropDbl ¶
Get the ideal-gas reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d3alpha0_dDelta2_dTau()
- CoolProp.AbstractState.d3alpha0_dDelta3(self) CoolPropDbl ¶
Get the ideal-gas reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d3alpha0_dDelta3()
- CoolProp.AbstractState.d3alpha0_dDelta_dTau2(self) CoolPropDbl ¶
Get the ideal-gas reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d3alpha0_dDelta_dTau2()
- CoolProp.AbstractState.d3alpha0_dTau3(self) CoolPropDbl ¶
Get the ideal-gas reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d3alpha0_dTau3()
- CoolProp.AbstractState.d3alphar_dDelta2_dTau(self) CoolPropDbl ¶
Get the residual reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d3alphar_dDelta2_dTau()
- CoolProp.AbstractState.d3alphar_dDelta3(self) CoolPropDbl ¶
Get the residual reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d3alphar_dDelta3()
- CoolProp.AbstractState.d3alphar_dDelta_dTau2(self) CoolPropDbl ¶
Get the residual reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d3alphar_dDelta_dTau2()
- CoolProp.AbstractState.d3alphar_dTau3(self) CoolPropDbl ¶
Get the residual reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d3alphar_dTau3()
- CoolProp.AbstractState.d4alphar_dDelta2_dTau2(self) CoolPropDbl ¶
Get the residual reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d4alphar_dDelta2_dTau2()
- CoolProp.AbstractState.d4alphar_dDelta3_dTau(self) CoolPropDbl ¶
Get the residual reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d4alphar_dDelta3_dTau()
- CoolProp.AbstractState.d4alphar_dDelta4(self) CoolPropDbl ¶
Get the residual reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d4alphar_dDelta4()
- CoolProp.AbstractState.d4alphar_dDelta_dTau3(self) CoolPropDbl ¶
Get the residual reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d4alphar_dDelta_dTau3()
- CoolProp.AbstractState.d4alphar_dTau4(self) CoolPropDbl ¶
Get the residual reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::d4alphar_dTau4()
- CoolProp.AbstractState.dalpha0_dDelta(self) CoolPropDbl ¶
Get the ideal-gas reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::dalpha0_dDelta()
- CoolProp.AbstractState.dalpha0_dTau(self) CoolPropDbl ¶
Get the ideal-gas reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::dalpha0_dTau()
- CoolProp.AbstractState.dalphar_dDelta(self) CoolPropDbl ¶
Get the residual reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::dalphar_dDelta()
- CoolProp.AbstractState.dalphar_dTau(self) CoolPropDbl ¶
Get the residual reduced Helmholtz energy - wrapper of c++ function CoolProp::AbstractState::dalphar_dTau()
- CoolProp.AbstractState.delta(self) double ¶
Get the reduced density - wrapper of c++ function CoolProp::AbstractState::delta(void)
- CoolProp.AbstractState.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()
- CoolProp.AbstractState.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()
- CoolProp.AbstractState.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()
- CoolProp.AbstractState.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()
- CoolProp.AbstractState.fluid_names(self)¶
Get the list of fluid names - wrapper of c++ function CoolProp::AbstractState::fluid_names()
- CoolProp.AbstractState.fluid_param_string(self, string key)¶
Get a fluid parameter string - wrapper of c++ function CoolProp::AbstractState::fluid_param_string()
- CoolProp.AbstractState.fugacity(self, size_t i) double ¶
Get the fugacity of the i-th component - wrapper of c++ function CoolProp::AbstractState::fugacity(std::size_t)
- CoolProp.AbstractState.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)
- CoolProp.AbstractState.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)
- CoolProp.AbstractState.gas_constant(self) double ¶
Get the gas constant in J/mol/K - wrapper of c++ function CoolProp::AbstractState::gas_constant(void)
- CoolProp.AbstractState.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()
- CoolProp.AbstractState.get_fluid_constant(self, size_t i, parameters param) double ¶
Get a constant for a fluid in the mixture CoolProp::AbstractState::get_fluid_constant()
- CoolProp.AbstractState.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()
- CoolProp.AbstractState.get_mass_fractions(self)¶
Get the mass fractions - wrapper of c++ function CoolProp::AbstractState::get_mass_fractions()
- CoolProp.AbstractState.get_mole_fractions(self)¶
Get the mole fractions - wrapper of c++ function CoolProp::AbstractState::get_mole_fractions()
- CoolProp.AbstractState.get_phase_envelope_data(self) PyPhaseEnvelopeData ¶
Get the phase envelope data - wrapper of c++ function CoolProp::AbstractState::get_phase_envelope_data()
- CoolProp.AbstractState.get_spinodal_data(self) PySpinodalData ¶
Get the data from the spinodal - wrapper of c++ function CoolProp::AbstractState::get_spinodal_data()
- CoolProp.AbstractState.gibbsmass(self) double ¶
Get the mass-specific Gibbs energy in J/kg - wrapper of c++ function CoolProp::AbstractState::gibbsmass(void)
- CoolProp.AbstractState.gibbsmass_excess(self) double ¶
Get the mass-specific excess Gibbs energy in J/kg - wrapper of c++ function CoolProp::AbstractState::gibbsmass_excess(void)
- CoolProp.AbstractState.gibbsmolar(self) double ¶
Get the mole-specific Gibbs energy in J/mol - wrapper of c++ function CoolProp::AbstractState::gibbsmolar(void)
- CoolProp.AbstractState.gibbsmolar_excess(self) double ¶
Get the mole-specific excess Gibbs energy in J/mol - wrapper of c++ function CoolProp::AbstractState::gibbsmolar_excess(void)
- CoolProp.AbstractState.gibbsmolar_residual(self) double ¶
Get the mole-specific residual Gibbs energy in J/mol - wrapper of c++ function CoolProp::AbstractState::gibbsmolar_residual(void)
- CoolProp.AbstractState.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()
- CoolProp.AbstractState.helmholtzmass(self) double ¶
Get the mass-specific Helmholtz energy in J/kg - wrapper of c++ function CoolProp::AbstractState::helmholtzmass(void)
- CoolProp.AbstractState.helmholtzmass_excess(self) double ¶
Get the mass-specific excess Helmholtz energy in J/kg - wrapper of c++ function CoolProp::AbstractState::helmholtzmass_excess(void)
- CoolProp.AbstractState.helmholtzmolar(self) double ¶
Get the mole-specific Helmholtz energy in J/mol - wrapper of c++ function CoolProp::AbstractState::helmholtzmolar(void)
- CoolProp.AbstractState.helmholtzmolar_excess(self) double ¶
Get the mole-specific excess Helmholtz energy in J/mol - wrapper of c++ function CoolProp::AbstractState::helmholtzmolar_excess(void)
- CoolProp.AbstractState.hmass(self) double ¶
Get the enthalpy in J/kg - wrapper of c++ function CoolProp::AbstractState::hmass(void)
- CoolProp.AbstractState.hmass_excess(self) double ¶
Get the mass-specific excess enthalpy in J/kg - wrapper of c++ function CoolProp::AbstractState::hmass_excess(void)
- CoolProp.AbstractState.hmolar(self) double ¶
Get the enthalpy in J/mol - wrapper of c++ function CoolProp::AbstractState::hmolar(void)
- CoolProp.AbstractState.hmolar_excess(self) double ¶
Get the mole-specific excess enthalpy in J/mol - wrapper of c++ function CoolProp::AbstractState::hmolar_excess(void)
- CoolProp.AbstractState.hmolar_residual(self) double ¶
Get the mole-specific residual enthalpy in J/mol - wrapper of c++ function CoolProp::AbstractState::hmolar_residual(void)
- CoolProp.AbstractState.ideal_curve(self, string type) tuple ¶
Get an ideal curve - wrapper of c++ function CoolProp::AbstractState::ideal_curve()
- CoolProp.AbstractState.isobaric_expansion_coefficient(self) double ¶
Get the isobaric expansion coefficient - wrapper of c++ function CoolProp::AbstractState::isobaric_expansion_coefficient(void)
- CoolProp.AbstractState.isothermal_compressibility(self) double ¶
Get the isothermal_compressibility - wrapper of c++ function CoolProp::AbstractState::isothermal_compressibility(void)
- CoolProp.AbstractState.keyed_output(self, parameters iOutput) double ¶
Get a keyed output CoolProp::AbstractState::keyed_output(parameters key)
- CoolProp.AbstractState.melting_line(self, int param, int given, double value) double ¶
Get values from the melting line - wrapper of c++ function CoolProp::AbstractState::melting_line()
- CoolProp.AbstractState.molar_mass(self) double ¶
Get the molar mass in kg/mol - wrapper of c++ function CoolProp::AbstractState::molar_mass(void)
- CoolProp.AbstractState.mole_fractions_liquid(self)¶
Get the mole fractions of the liquid phase - wrapper of c++ function CoolProp::AbstractState::mole_fractions_liquid(void)
- CoolProp.AbstractState.mole_fractions_vapor(self)¶
Get the mole fractions of the vapor phase - wrapper of c++ function CoolProp::AbstractState::mole_fractions_vapor(void)
- CoolProp.AbstractState.name(self)¶
Get the fluid name - wrapper of c++ function CoolProp::AbstractState::name()
- CoolProp.AbstractState.p(self) double ¶
Get the pressure in Pa - wrapper of c++ function CoolProp::AbstractState::p(void)
- CoolProp.AbstractState.p_critical(self) double ¶
Gets the critical pressure in Pa - wrapper of c++ function CoolProp::AbstractState::p_critical()
- CoolProp.AbstractState.phase(self) phases ¶
Get the phase as key value- wrapper of c++ function CoolProp::AbstractState::phase()
- CoolProp.AbstractState.pmax(self) double ¶
Set the maximum pressure in Pa - wrapper of c++ function CoolProp::AbstractState::pmax()
- CoolProp.AbstractState.rhomass(self) double ¶
Get the density in kg/m^3 - wrapper of c++ function CoolProp::AbstractState::rhomass(void)
- CoolProp.AbstractState.rhomass_critical(self) double ¶
Gets the critical density in kg/m^3 - wrapper of c++ function CoolProp::AbstractState::rhomass_critical()
- CoolProp.AbstractState.rhomass_reducing(self) double ¶
Gets the reducing density in kg/m^3 - wrapper of c++ function CoolProp::AbstractState::rhomass_reducing()
- CoolProp.AbstractState.rhomolar(self) double ¶
Get the density in mol/m^3 - wrapper of c++ function CoolProp::AbstractState::rhomolar(void)
- CoolProp.AbstractState.rhomolar_critical(self) double ¶
Gets the critical density in mol/m^3 - wrapper of c++ function CoolProp::AbstractState::rhomolar_critical()
- CoolProp.AbstractState.rhomolar_reducing(self) double ¶
Gets the reducing density in mol/m^3 - wrapper of c++ function CoolProp::AbstractState::rhomolar_reducing()
- CoolProp.AbstractState.saturated_liquid_keyed_output(self, parameters iOutput) double ¶
Get a trivial output for the saturated liquid CoolProp::AbstractState::saturated_liquid_keyed_output(parameters key)
- CoolProp.AbstractState.saturated_vapor_keyed_output(self, parameters iOutput) double ¶
Get a trivial output for the saturated vapor CoolProp::AbstractState::saturated_vapor_keyed_output(parameters key)
- CoolProp.AbstractState.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()
- CoolProp.AbstractState.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()
- CoolProp.AbstractState.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()
- CoolProp.AbstractState.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()
- CoolProp.AbstractState.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()
- CoolProp.AbstractState.set_mass_fractions(self, vector[double] z)¶
Set the mass fractions - wrapper of c++ function CoolProp::AbstractState::set_mass_fractions()
- CoolProp.AbstractState.set_mole_fractions(self, vector[double] z)¶
Set the mole fractions - wrapper of c++ function CoolProp::AbstractState::set_mole_fractions()
- CoolProp.AbstractState.set_volu_fractions(self, vector[double] z)¶
Set the volume fractions - wrapper of c++ function CoolProp::AbstractState::set_volu_fractions()
- CoolProp.AbstractState.smass(self) double ¶
Get the entropy in J/kg/K - wrapper of c++ function CoolProp::AbstractState::smass(void)
- CoolProp.AbstractState.smass_excess(self) double ¶
Get the mass-specific excess entropy in J/kg/K - wrapper of c++ function CoolProp::AbstractState::smass_excess(void)
- CoolProp.AbstractState.smolar(self) double ¶
Get the entropy in J/mol/K - wrapper of c++ function CoolProp::AbstractState::smolar(void)
- CoolProp.AbstractState.smolar_excess(self) double ¶
Get the mole-specific excess entropy in J/mol/K - wrapper of c++ function CoolProp::AbstractState::smolar_excess(void)
- CoolProp.AbstractState.smolar_residual(self) double ¶
Get the mole-specific residual entropy in J/mol/K - wrapper of c++ function CoolProp::AbstractState::smolar_residual(void)
- CoolProp.AbstractState.specify_phase(self, phases phase)¶
Specify the phase - wrapper of c++ function CoolProp::AbstractState::specify_phase()
- CoolProp.AbstractState.speed_sound(self) double ¶
Get the speed of sound in m/s - wrapper of c++ function CoolProp::AbstractState::speed_sound(void)
- CoolProp.AbstractState.surface_tension(self) double ¶
Get the surface tension N/m - wrapper of c++ function CoolProp::AbstractState::surface_tension(void)
- CoolProp.AbstractState.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()
- CoolProp.AbstractState.tau(self) double ¶
Get the reciprocal reduced temperature - wrapper of c++ function CoolProp::AbstractState::tau(void)
- CoolProp.AbstractState.trivial_keyed_output(self, parameters iOutput) double ¶
Get a trivial keyed output not requiring any iteration CoolProp::AbstractState::trivial_keyed_output(parameters key)
- CoolProp.AbstractState.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()
- CoolProp.AbstractState.umass(self) double ¶
Get the internal energy in J/kg - wrapper of c++ function CoolProp::AbstractState::umass(void)
- CoolProp.AbstractState.umass_excess(self) double ¶
Get the mass-specific excess internal energy in J/kg - wrapper of c++ function CoolProp::AbstractState::umass_excess(void)
- CoolProp.AbstractState.umolar(self) double ¶
Get the internal energy in J/mol - wrapper of c++ function CoolProp::AbstractState::umolar(void)
- CoolProp.AbstractState.umolar_excess(self) double ¶
Get the mole-specific excess internal energy in J/mol - wrapper of c++ function CoolProp::AbstractState::umolar_excess(void)
- CoolProp.AbstractState.unspecify_phase(self)¶
Unspecify the phase - wrapper of c++ function CoolProp::AbstractState::unspecify_phase()
- CoolProp.AbstractState.update(self, input_pairs ipair, double Value1, double Value2)¶
Update function - wrapper of c++ function CoolProp::AbstractState::update()
- CoolProp.AbstractState.update_with_guesses(self, input_pairs ipair, double Value1, double Value2, PyGuessesStructure guesses)¶
Update function - wrapper of c++ function CoolProp::AbstractState::update()
- CoolProp.AbstractState.viscosity(self) double ¶
Get the viscosity in Pa-s - wrapper of c++ function CoolProp::AbstractState::viscosity(void)
- CoolProp.AbstractState.viscosity_contributions(self) dict ¶
Retrieve each of the contributions to the viscosity, each in Pa-s - wrapper of c++ function CoolProp::AbstractState::viscosity_contributions()
- CoolProp.AbstractState.volumemass_excess(self) double ¶
Get the mass-specific excess volume in m^3/kg - wrapper of c++ function CoolProp::AbstractState::volumemass_excess(void)
- CoolProp.AbstractState.volumemolar_excess(self) double ¶
Get the mole-specific excess volume in m^3/mol - wrapper of c++ function CoolProp::AbstractState::volumemolar_excess(void)