|
CoolProp 7.2.0
An open-source fluid property and humid air property database
|
The mother of all state classes.
This class provides the basic properties based on interrelations of the properties, their derivatives and the Helmholtz energy terms. It does not provide the mechanism to update the values. This has to be implemented in a subclass. Most functions are defined as virtual functions allowing us redefine them later, for example to implement the TTSE technique. The functions defined here are always used as a fall-back.
This base class does not perform any checks on the two-phase conditions and alike. Most of the functions defined here only apply to compressible single state substances. Make sure you are aware of all the assumptions we made when using this class.
Add build table function to Abstract State Interpolator inherit AS implemented by TTSE BICUBIC
Definition at line 78 of file AbstractState.h.
#include <AbstractState.h>
Public Member Functions | |
| AbstractState () | |
| virtual | ~AbstractState () |
| void | set_T (CoolPropDbl T) |
| Set the internal variable T without a flash call (expert use only!) More... | |
| virtual std::string | backend_name (void)=0 |
| virtual bool | using_mole_fractions (void)=0 |
| virtual bool | using_mass_fractions (void)=0 |
| virtual bool | using_volu_fractions (void)=0 |
| virtual void | set_mole_fractions (const std::vector< CoolPropDbl > &mole_fractions)=0 |
| virtual void | set_mass_fractions (const std::vector< CoolPropDbl > &mass_fractions)=0 |
| virtual void | set_volu_fractions (const std::vector< CoolPropDbl > &mass_fractions) |
| virtual void | set_reference_stateS (const std::string &reference_state) |
| Set the reference state based on a string representation. More... | |
| virtual void | set_reference_stateD (double T, double rhomolar, double hmolar0, double smolar0) |
| std::vector< CoolPropDbl > | mole_fractions_liquid (void) |
| Get the mole fractions of the equilibrium liquid phase. More... | |
| std::vector< double > | mole_fractions_liquid_double (void) |
| Get the mole fractions of the equilibrium liquid phase (but as a double for use in SWIG wrapper) More... | |
| std::vector< CoolPropDbl > | mole_fractions_vapor (void) |
| Get the mole fractions of the equilibrium vapor phase. More... | |
| std::vector< double > | mole_fractions_vapor_double (void) |
| Get the mole fractions of the equilibrium vapor phase (but as a double for use in SWIG wrapper) More... | |
| virtual const std::vector< CoolPropDbl > & | get_mole_fractions (void)=0 |
| Get the mole fractions of the fluid. More... | |
| virtual const std::vector< CoolPropDbl > | get_mass_fractions (void) |
| Get the mass fractions of the fluid. More... | |
| virtual void | update (CoolProp::input_pairs input_pair, double Value1, double Value2)=0 |
| Update the state using two state variables. More... | |
| virtual void | update_QT_pure_superanc (double Q, double T) |
| Update the state for QT inputs for pure fluids when using the superancillary functions. More... | |
| virtual void | update_with_guesses (CoolProp::input_pairs input_pair, double Value1, double Value2, const GuessesStructure &guesses) |
| virtual bool | available_in_high_level (void) |
| virtual std::string | fluid_param_string (const std::string &) |
| Return a string from the backend for the mixture/fluid - backend dependent - could be CAS #, name, etc. More... | |
| std::vector< std::string > | fluid_names (void) |
| Return a vector of strings of the fluid names that are in use. More... | |
| virtual const double | get_fluid_constant (std::size_t i, parameters param) const |
| virtual void | set_binary_interaction_double (const std::string &CAS1, const std::string &CAS2, const std::string ¶meter, const double value) |
| Set binary mixture floating point parameter (EXPERT USE ONLY!!!) More... | |
| virtual void | set_binary_interaction_double (const std::size_t i, const std::size_t j, const std::string ¶meter, const double value) |
| Set binary mixture floating point parameter (EXPERT USE ONLY!!!) More... | |
| virtual void | set_binary_interaction_string (const std::string &CAS1, const std::string &CAS2, const std::string ¶meter, const std::string &value) |
| Set binary mixture string parameter (EXPERT USE ONLY!!!) More... | |
| virtual void | set_binary_interaction_string (const std::size_t i, const std::size_t j, const std::string ¶meter, const std::string &value) |
| Set binary mixture string parameter (EXPERT USE ONLY!!!) More... | |
| virtual double | get_binary_interaction_double (const std::string &CAS1, const std::string &CAS2, const std::string ¶meter) |
| Get binary mixture double value (EXPERT USE ONLY!!!) More... | |
| virtual double | get_binary_interaction_double (const std::size_t i, const std::size_t j, const std::string ¶meter) |
| Get binary mixture double value (EXPERT USE ONLY!!!) More... | |
| virtual std::string | get_binary_interaction_string (const std::string &CAS1, const std::string &CAS2, const std::string ¶meter) |
| Get binary mixture string value (EXPERT USE ONLY!!!) More... | |
| virtual void | apply_simple_mixing_rule (std::size_t i, std::size_t j, const std::string &model) |
| Apply a simple mixing rule (EXPERT USE ONLY!!!) More... | |
| virtual void | set_cubic_alpha_C (const size_t i, const std::string ¶meter, const double c1, const double c2, const double c3) |
| Set the cubic alpha function's constants: More... | |
| virtual void | set_fluid_parameter_double (const size_t i, const std::string ¶meter, const double value) |
| Set fluid parameter (currently the volume translation parameter for cubic) More... | |
| virtual double | get_fluid_parameter_double (const size_t i, const std::string ¶meter) |
| Double fluid parameter (currently the volume translation parameter for cubic) More... | |
| virtual bool | clear () |
| Clear all the cached values. More... | |
| virtual bool | clear_comp_change () |
| When the composition changes, clear all cached values that are only dependent on composition, but not the thermodynamic state. More... | |
| virtual const CoolProp::SimpleState & | get_reducing_state () |
| const CoolProp::SimpleState & | get_state (const std::string &state) |
| Get a desired state point - backend dependent. More... | |
| double | Tmin (void) |
| Get the minimum temperature in K. More... | |
| double | Tmax (void) |
| Get the maximum temperature in K. More... | |
| double | pmax (void) |
| Get the maximum pressure in Pa. More... | |
| double | Ttriple (void) |
| Get the triple point temperature in K. More... | |
| phases | phase (void) |
| Get the phase of the state. More... | |
| void | specify_phase (phases phase) |
| Specify the phase for all further calculations with this state class. More... | |
| void | unspecify_phase (void) |
| Unspecify the phase and go back to calculating it based on the inputs. More... | |
| double | T_critical (void) |
| Return the critical temperature in K. More... | |
| double | p_critical (void) |
| Return the critical pressure in Pa. More... | |
| double | rhomolar_critical (void) |
| Return the critical molar density in mol/m^3. More... | |
| double | rhomass_critical (void) |
| Return the critical mass density in kg/m^3. More... | |
| std::vector< CriticalState > | all_critical_points (void) |
| Return the vector of critical points, including points that are unstable or correspond to negative pressure. More... | |
| void | build_spinodal () |
| Construct the spinodal curve for the mixture (or pure fluid) More... | |
| SpinodalData | get_spinodal_data () |
| Get the data from the spinodal curve constructed in the call to build_spinodal() More... | |
| void | criticality_contour_values (double &L1star, double &M1star) |
| Calculate the criticality contour values \(\mathcal{L}_1^*\) and \(\mathcal{M}_1^*\). More... | |
| double | tangent_plane_distance (const double T, const double p, const std::vector< double > &w, const double rhomolar_guess=-1) |
| double | T_reducing (void) |
| Return the reducing point temperature in K. More... | |
| double | rhomolar_reducing (void) |
| Return the molar density at the reducing point in mol/m^3. More... | |
| double | rhomass_reducing (void) |
| Return the mass density at the reducing point in kg/m^3. More... | |
| double | p_triple (void) |
| Return the triple point pressure in Pa. More... | |
| std::string | name () |
| Return the name - backend dependent. More... | |
| std::string | description () |
| Return the description - backend dependent. More... | |
| double | dipole_moment () |
| Return the dipole moment in C-m (1 D = 3.33564e-30 C-m) More... | |
| double | keyed_output (parameters key) |
| Retrieve a value by key. More... | |
| double | trivial_keyed_output (parameters key) |
| A trivial keyed output like molar mass that does not depend on the state. More... | |
| double | saturated_liquid_keyed_output (parameters key) |
| Get an output from the saturated liquid state by key. More... | |
| double | saturated_vapor_keyed_output (parameters key) |
| Get an output from the saturated vapor state by key. More... | |
| double | T (void) |
| Return the temperature in K. More... | |
| double | rhomolar (void) |
| Return the molar density in mol/m^3. More... | |
| double | rhomass (void) |
| Return the mass density in kg/m^3. More... | |
| double | p (void) |
| Return the pressure in Pa. More... | |
| double | Q (void) |
| Return the vapor quality (mol/mol); Q = 0 for saturated liquid. More... | |
| double | tau (void) |
| Return the reciprocal of the reduced temperature ( \(\tau = T_c/T\)) More... | |
| double | delta (void) |
| Return the reduced density ( \(\delta = \rho/\rho_c\)) More... | |
| double | molar_mass (void) |
| Return the molar mass in kg/mol. More... | |
| double | acentric_factor (void) |
| Return the acentric factor. More... | |
| double | gas_constant (void) |
| Return the mole-fraction weighted gas constant in J/mol/K. More... | |
| double | Bvirial (void) |
| Return the B virial coefficient. More... | |
| double | dBvirial_dT (void) |
| Return the derivative of the B virial coefficient with respect to temperature. More... | |
| double | Cvirial (void) |
| Return the C virial coefficient. More... | |
| double | dCvirial_dT (void) |
| Return the derivative of the C virial coefficient with respect to temperature. More... | |
| double | compressibility_factor (void) |
| Return the compressibility factor \( Z = p/(rho R T) \). More... | |
| double | hmolar (void) |
| Return the molar enthalpy in J/mol. More... | |
| double | hmolar_residual (void) |
| Return the residual molar enthalpy in J/mol. More... | |
| double | hmolar_idealgas (void) |
| Return the ideal gas molar enthalpy in J/mol. More... | |
| double | hmass_idealgas (void) |
| Return the ideal gas specific enthalpy in J/kg. More... | |
| double | hmass (void) |
| Return the mass enthalpy in J/kg. More... | |
| double | hmolar_excess (void) |
| Return the excess molar enthalpy in J/mol. More... | |
| double | hmass_excess (void) |
| Return the excess mass enthalpy in J/kg. More... | |
| double | smolar (void) |
| Return the molar entropy in J/mol/K. More... | |
| double | smolar_residual (void) |
| Return the residual molar entropy (as a function of temperature and density) in J/mol/K. More... | |
| double | smolar_idealgas (void) |
| Return the ideal gas molar entropy in J/mol/K. More... | |
| double | smass_idealgas (void) |
| Return the ideal gas specific entropy in J/kg/K. More... | |
| double | neff (void) |
| Return the effective hardness of interaction. More... | |
| double | smass (void) |
| Return the molar entropy in J/kg/K. More... | |
| double | smolar_excess (void) |
| Return the molar entropy in J/mol/K. More... | |
| double | smass_excess (void) |
| Return the molar entropy in J/kg/K. More... | |
| double | umolar (void) |
| Return the molar internal energy in J/mol. More... | |
| double | umass (void) |
| Return the mass internal energy in J/kg. More... | |
| double | umolar_excess (void) |
| Return the excess internal energy in J/mol. More... | |
| double | umass_excess (void) |
| Return the excess internal energy in J/kg. More... | |
| double | umolar_idealgas (void) |
| Return the ideal gas molar internal energy in J/mol. More... | |
| double | umass_idealgas (void) |
| Return the ideal gas specific internal energy in J/kg. More... | |
| double | cpmolar (void) |
| Return the molar constant pressure specific heat in J/mol/K. More... | |
| double | cpmass (void) |
| Return the mass constant pressure specific heat in J/kg/K. More... | |
| double | cp0molar (void) |
| Return the molar constant pressure specific heat for ideal gas part only in J/mol/K. More... | |
| double | cp0mass (void) |
| Return the mass constant pressure specific heat for ideal gas part only in J/kg/K. More... | |
| double | cvmolar (void) |
| Return the molar constant volume specific heat in J/mol/K. More... | |
| double | cvmass (void) |
| Return the mass constant volume specific heat in J/kg/K. More... | |
| double | gibbsmolar (void) |
| Return the Gibbs energy in J/mol. More... | |
| double | gibbsmolar_residual (void) |
| Return the residual Gibbs energy in J/mol. More... | |
| double | gibbsmass (void) |
| Return the Gibbs energy in J/kg. More... | |
| double | gibbsmolar_excess (void) |
| Return the excess Gibbs energy in J/mol. More... | |
| double | gibbsmass_excess (void) |
| Return the excess Gibbs energy in J/kg. More... | |
| double | helmholtzmolar (void) |
| Return the Helmholtz energy in J/mol. More... | |
| double | helmholtzmass (void) |
| Return the Helmholtz energy in J/kg. More... | |
| double | helmholtzmolar_excess (void) |
| Return the excess Helmholtz energy in J/mol. More... | |
| double | helmholtzmass_excess (void) |
| Return the excess Helmholtz energy in J/kg. More... | |
| double | volumemolar_excess (void) |
| Return the excess volume in m^3/mol. More... | |
| double | volumemass_excess (void) |
| Return the excess volume in m^3/kg. More... | |
| double | speed_sound (void) |
| Return the speed of sound in m/s. More... | |
| double | isothermal_compressibility (void) |
| Return the isothermal compressibility \( \kappa = -\frac{1}{v}\left.\frac{\partial v}{\partial p}\right|_T=\frac{1}{\rho}\left.\frac{\partial \rho}{\partial p}\right|_T\) in 1/Pa. More... | |
| double | isobaric_expansion_coefficient (void) |
| Return the isobaric expansion coefficient \( \beta = \frac{1}{v}\left.\frac{\partial v}{\partial T}\right|_p = -\frac{1}{\rho}\left.\frac{\partial \rho}{\partial T}\right|_p\) in 1/K. More... | |
| double | isentropic_expansion_coefficient (void) |
| Return the isentropic expansion coefficient \( \kappa_s = -\frac{c_p}{c_v}\frac{v}{p}\left.\frac{\partial p}{\partial v}\right|_T = \frac{\rho}{p}\left.\frac{\partial p}{\partial \rho}\right|_s\). More... | |
| double | fugacity_coefficient (std::size_t i) |
| Return the fugacity coefficient of the i-th component of the mixture. More... | |
| std::vector< double > | fugacity_coefficients () |
| Return a vector of the fugacity coefficients for all components in the mixture. More... | |
| double | fugacity (std::size_t i) |
| Return the fugacity of the i-th component of the mixture. More... | |
| double | chemical_potential (std::size_t i) |
| Return the chemical potential of the i-th component of the mixture. More... | |
| double | fundamental_derivative_of_gas_dynamics (void) |
| Return the fundamental derivative of gas dynamics \( \Gamma \). More... | |
| double | PIP () |
| Return the phase identification parameter (PIP) of G. Venkatarathnam and L.R. Oellrich, "Identification of the phase of a fluid using partial derivatives of pressure, volume, and temperature without reference to saturation properties: Applications in phase equilibria calculations". More... | |
| void | true_critical_point (double &T, double &rho) |
| Calculate the "true" critical point for pure fluids where dpdrho|T and d2p/drho2|T are equal to zero. More... | |
| void | ideal_curve (const std::string &type, std::vector< double > &T, std::vector< double > &p) |
| Calculate an ideal curve for a pure fluid. More... | |
| CoolPropDbl | first_partial_deriv (parameters Of, parameters Wrt, parameters Constant) |
| The first partial derivative in homogeneous phases. More... | |
| CoolPropDbl | second_partial_deriv (parameters Of1, parameters Wrt1, parameters Constant1, parameters Wrt2, parameters Constant2) |
| The second partial derivative in homogeneous phases. More... | |
| CoolPropDbl | first_saturation_deriv (parameters Of1, parameters Wrt1) |
| The first partial derivative along the saturation curve. More... | |
| CoolPropDbl | second_saturation_deriv (parameters Of1, parameters Wrt1, parameters Wrt2) |
| The second partial derivative along the saturation curve. More... | |
| double | first_two_phase_deriv (parameters Of, parameters Wrt, parameters Constant) |
| Calculate the first "two-phase" derivative as described by Thorade and Sadaat, EAS, 2013. More... | |
| double | second_two_phase_deriv (parameters Of, parameters Wrt1, parameters Constant1, parameters Wrt2, parameters Constant2) |
| Calculate the second "two-phase" derivative as described by Thorade and Sadaat, EAS, 2013. More... | |
| double | first_two_phase_deriv_splined (parameters Of, parameters Wrt, parameters Constant, double x_end) |
| Calculate the first "two-phase" derivative as described by Thorade and Sadaat, EAS, 2013. More... | |
| void | build_phase_envelope (const std::string &type="") |
| Construct the phase envelope for a mixture. More... | |
| const CoolProp::PhaseEnvelopeData & | get_phase_envelope_data () |
| After having calculated the phase envelope, return the phase envelope data. More... | |
| virtual bool | has_melting_line (void) |
| Return true if the fluid has a melting line - default is false, but can be re-implemented by derived class. More... | |
| double | melting_line (int param, int given, double value) |
| double | saturation_ancillary (parameters param, int Q, parameters given, double value) |
| double | viscosity (void) |
| Return the viscosity in Pa-s. More... | |
| void | viscosity_contributions (CoolPropDbl &dilute, CoolPropDbl &initial_density, CoolPropDbl &residual, CoolPropDbl &critical) |
| Return the viscosity contributions, each in Pa-s. More... | |
| double | conductivity (void) |
| Return the thermal conductivity in W/m/K. More... | |
| void | conductivity_contributions (CoolPropDbl &dilute, CoolPropDbl &initial_density, CoolPropDbl &residual, CoolPropDbl &critical) |
| Return the thermal conductivity contributions, each in W/m/K. More... | |
| double | surface_tension (void) |
| Return the surface tension in N/m. More... | |
| double | Prandtl (void) |
| Return the Prandtl number (dimensionless) More... | |
| void | conformal_state (const std::string &reference_fluid, CoolPropDbl &T, CoolPropDbl &rhomolar) |
| Find the conformal state needed for ECS. More... | |
| void | change_EOS (const std::size_t i, const std::string &EOS_name) |
| Change the equation of state for a given component to a specified EOS. More... | |
| CoolPropDbl | alpha0 (void) |
| Return the term \( \alpha^0 \). More... | |
| CoolPropDbl | dalpha0_dDelta (void) |
| Return the term \( \alpha^0_{\delta} \). More... | |
| CoolPropDbl | dalpha0_dTau (void) |
| Return the term \( \alpha^0_{\tau} \). More... | |
| CoolPropDbl | d2alpha0_dDelta2 (void) |
| Return the term \( \alpha^0_{\delta\delta} \). More... | |
| CoolPropDbl | d2alpha0_dDelta_dTau (void) |
| Return the term \( \alpha^0_{\delta\tau} \). More... | |
| CoolPropDbl | d2alpha0_dTau2 (void) |
| Return the term \( \alpha^0_{\tau\tau} \). More... | |
| CoolPropDbl | d3alpha0_dTau3 (void) |
| Return the term \( \alpha^0_{\tau\tau\tau} \). More... | |
| CoolPropDbl | d3alpha0_dDelta_dTau2 (void) |
| Return the term \( \alpha^0_{\delta\tau\tau} \). More... | |
| CoolPropDbl | d3alpha0_dDelta2_dTau (void) |
| Return the term \( \alpha^0_{\delta\delta\tau} \). More... | |
| CoolPropDbl | d3alpha0_dDelta3 (void) |
| Return the term \( \alpha^0_{\delta\delta\delta} \). More... | |
| CoolPropDbl | alphar (void) |
| Return the term \( \alpha^r \). More... | |
| CoolPropDbl | dalphar_dDelta (void) |
| Return the term \( \alpha^r_{\delta} \). More... | |
| CoolPropDbl | dalphar_dTau (void) |
| Return the term \( \alpha^r_{\tau} \). More... | |
| CoolPropDbl | d2alphar_dDelta2 (void) |
| Return the term \( \alpha^r_{\delta\delta} \). More... | |
| CoolPropDbl | d2alphar_dDelta_dTau (void) |
| Return the term \( \alpha^r_{\delta\tau} \). More... | |
| CoolPropDbl | d2alphar_dTau2 (void) |
| Return the term \( \alpha^r_{\tau\tau} \). More... | |
| CoolPropDbl | d3alphar_dDelta3 (void) |
| Return the term \( \alpha^r_{\delta\delta\delta} \). More... | |
| CoolPropDbl | d3alphar_dDelta2_dTau (void) |
| Return the term \( \alpha^r_{\delta\delta\tau} \). More... | |
| CoolPropDbl | d3alphar_dDelta_dTau2 (void) |
| Return the term \( \alpha^r_{\delta\tau\tau} \). More... | |
| CoolPropDbl | d3alphar_dTau3 (void) |
| Return the term \( \alpha^r_{\tau\tau\tau} \). More... | |
| CoolPropDbl | d4alphar_dDelta4 (void) |
| Return the term \( \alpha^r_{\delta\delta\delta\delta} \). More... | |
| CoolPropDbl | d4alphar_dDelta3_dTau (void) |
| Return the term \( \alpha^r_{\delta\delta\delta\tau} \). More... | |
| CoolPropDbl | d4alphar_dDelta2_dTau2 (void) |
| Return the term \( \alpha^r_{\delta\delta\tau\tau} \). More... | |
| CoolPropDbl | d4alphar_dDelta_dTau3 (void) |
| Return the term \( \alpha^r_{\delta\tau\tau\tau} \). More... | |
| CoolPropDbl | d4alphar_dTau4 (void) |
| Return the term \( \alpha^r_{\tau\tau\tau\tau} \). More... | |
Static Public Member Functions | |
| static AbstractState * | factory (const std::string &backend, const std::string &fluid_names) |
| A factory function to return a pointer to a new-allocated instance of one of the backends. More... | |
| static AbstractState * | factory (const std::string &backend, const std::vector< std::string > &fluid_names) |
| A factory function to return a pointer to a new-allocated instance of one of the backends. More... | |
Protected Types | |
| using | CAE = CacheArrayElement< double > |
Protected Member Functions | |
| bool | isSupercriticalPhase (void) |
| bool | isHomogeneousPhase (void) |
| bool | isTwoPhase (void) |
| virtual CoolPropDbl | calc_hmolar (void) |
| Using this backend, calculate the molar enthalpy in J/mol. More... | |
| virtual CoolPropDbl | calc_hmolar_residual (void) |
| Using this backend, calculate the residual molar enthalpy in J/mol. More... | |
| virtual CoolPropDbl | calc_smolar (void) |
| Using this backend, calculate the molar entropy in J/mol/K. More... | |
| virtual CoolPropDbl | calc_smolar_residual (void) |
| Using this backend, calculate the residual molar entropy in J/mol/K. More... | |
| virtual CoolPropDbl | calc_neff (void) |
| Using this backend, calculate effective hardness of interaction. More... | |
| virtual CoolPropDbl | calc_umolar (void) |
| Using this backend, calculate the molar internal energy in J/mol. More... | |
| virtual CoolPropDbl | calc_cpmolar (void) |
| Using this backend, calculate the molar constant-pressure specific heat in J/mol/K. More... | |
| virtual CoolPropDbl | calc_cpmolar_idealgas (void) |
| Using this backend, calculate the ideal gas molar constant-pressure specific heat in J/mol/K. More... | |
| virtual CoolPropDbl | calc_cvmolar (void) |
| Using this backend, calculate the molar constant-volume specific heat in J/mol/K. More... | |
| virtual CoolPropDbl | calc_gibbsmolar (void) |
| Using this backend, calculate the molar Gibbs function in J/mol. More... | |
| virtual CoolPropDbl | calc_gibbsmolar_residual (void) |
| Using this backend, calculate the residual molar Gibbs function in J/mol. More... | |
| virtual CoolPropDbl | calc_helmholtzmolar (void) |
| Using this backend, calculate the molar Helmholtz energy in J/mol. More... | |
| virtual CoolPropDbl | calc_speed_sound (void) |
| Using this backend, calculate the speed of sound in m/s. More... | |
| virtual CoolPropDbl | calc_isothermal_compressibility (void) |
| Using this backend, calculate the isothermal compressibility \( \kappa = -\frac{1}{v}\left.\frac{\partial v}{\partial p}\right|_T=\frac{1}{\rho}\left.\frac{\partial \rho}{\partial p}\right|_T\) in 1/Pa. More... | |
| virtual CoolPropDbl | calc_isobaric_expansion_coefficient (void) |
| Using this backend, calculate the isobaric expansion coefficient \( \beta = \frac{1}{v}\left.\frac{\partial v}{\partial T}\right|_p = -\frac{1}{\rho}\left.\frac{\partial \rho}{\partial T}\right|_p\) in 1/K. More... | |
| virtual CoolPropDbl | calc_isentropic_expansion_coefficient (void) |
| Using this backend, calculate the isentropic expansion coefficient \( \kappa_s = -\frac{c_p}{c_v}\frac{v}{p}\left.\frac{\partial p}{\partial v}\right|_T = \frac{\rho}{p}\left.\frac{\partial p}{\partial \rho}\right|_s\). More... | |
| virtual CoolPropDbl | calc_viscosity (void) |
| Using this backend, calculate the viscosity in Pa-s. More... | |
| virtual CoolPropDbl | calc_conductivity (void) |
| Using this backend, calculate the thermal conductivity in W/m/K. More... | |
| virtual CoolPropDbl | calc_surface_tension (void) |
| Using this backend, calculate the surface tension in N/m. More... | |
| virtual CoolPropDbl | calc_molar_mass (void) |
| Using this backend, calculate the molar mass in kg/mol. More... | |
| virtual CoolPropDbl | calc_acentric_factor (void) |
| Using this backend, calculate the acentric factor. More... | |
| virtual CoolPropDbl | calc_pressure (void) |
| Using this backend, calculate the pressure in Pa. More... | |
| virtual CoolPropDbl | calc_gas_constant (void) |
| Using this backend, calculate the universal gas constant \(R_u\) in J/mol/K. More... | |
| virtual CoolPropDbl | calc_fugacity_coefficient (std::size_t i) |
| Using this backend, calculate the fugacity coefficient (dimensionless) More... | |
| virtual std::vector< CoolPropDbl > | calc_fugacity_coefficients () |
| Using this backend, calculate the fugacity in Pa. More... | |
| virtual CoolPropDbl | calc_fugacity (std::size_t i) |
| Using this backend, calculate the fugacity in Pa. More... | |
| virtual CoolPropDbl | calc_chemical_potential (std::size_t i) |
| Using this backend, calculate the chemical potential in J/mol. More... | |
| virtual CoolPropDbl | calc_PIP (void) |
| Using this backend, calculate the phase identification parameter (PIP) More... | |
| virtual void | calc_excess_properties (void) |
| Using this backend, calculate and cache the excess properties. More... | |
| virtual CoolPropDbl | calc_alphar (void) |
| Using this backend, calculate the residual Helmholtz energy term \(\alpha^r\) (dimensionless) More... | |
| virtual CoolPropDbl | calc_dalphar_dDelta (void) |
| Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta}\) (dimensionless) More... | |
| virtual CoolPropDbl | calc_dalphar_dTau (void) |
| Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\tau}\) (dimensionless) More... | |
| virtual CoolPropDbl | calc_d2alphar_dDelta2 (void) |
| Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\delta}\) (dimensionless) More... | |
| virtual CoolPropDbl | calc_d2alphar_dDelta_dTau (void) |
| Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\tau}\) (dimensionless) More... | |
| virtual CoolPropDbl | calc_d2alphar_dTau2 (void) |
| Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\tau\tau}\) (dimensionless) More... | |
| virtual CoolPropDbl | calc_d3alphar_dDelta3 (void) |
| Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\delta\delta}\) (dimensionless) More... | |
| virtual CoolPropDbl | calc_d3alphar_dDelta2_dTau (void) |
| Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\delta\tau}\) (dimensionless) More... | |
| virtual CoolPropDbl | calc_d3alphar_dDelta_dTau2 (void) |
| Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\tau\tau}\) (dimensionless) More... | |
| virtual CoolPropDbl | calc_d3alphar_dTau3 (void) |
| Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\tau\tau\tau}\) (dimensionless) More... | |
| virtual CoolPropDbl | calc_d4alphar_dDelta4 (void) |
| Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\delta\delta\delta}\) (dimensionless) More... | |
| virtual CoolPropDbl | calc_d4alphar_dDelta3_dTau (void) |
| Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\delta\delta\tau}\) (dimensionless) More... | |
| virtual CoolPropDbl | calc_d4alphar_dDelta2_dTau2 (void) |
| Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\delta\tau\tau}\) (dimensionless) More... | |
| virtual CoolPropDbl | calc_d4alphar_dDelta_dTau3 (void) |
| Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\tau\tau\tau}\) (dimensionless) More... | |
| virtual CoolPropDbl | calc_d4alphar_dTau4 (void) |
| Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\tau\tau\tau\tau}\) (dimensionless) More... | |
| virtual CoolPropDbl | calc_alpha0 (void) |
| Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0\) (dimensionless) More... | |
| virtual CoolPropDbl | calc_dalpha0_dDelta (void) |
| Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\delta}\) (dimensionless) More... | |
| virtual CoolPropDbl | calc_dalpha0_dTau (void) |
| Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\tau}\) (dimensionless) More... | |
| virtual CoolPropDbl | calc_d2alpha0_dDelta_dTau (void) |
| Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\delta\delta}\) (dimensionless) More... | |
| virtual CoolPropDbl | calc_d2alpha0_dDelta2 (void) |
| Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\delta\tau}\) (dimensionless) More... | |
| virtual CoolPropDbl | calc_d2alpha0_dTau2 (void) |
| Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\tau\tau}\) (dimensionless) More... | |
| virtual CoolPropDbl | calc_d3alpha0_dDelta3 (void) |
| Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\delta\delta\delta}\) (dimensionless) More... | |
| virtual CoolPropDbl | calc_d3alpha0_dDelta2_dTau (void) |
| Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\delta\delta\tau}\) (dimensionless) More... | |
| virtual CoolPropDbl | calc_d3alpha0_dDelta_dTau2 (void) |
| Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\delta\tau\tau}\) (dimensionless) More... | |
| virtual CoolPropDbl | calc_d3alpha0_dTau3 (void) |
| Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\tau\tau\tau}\) (dimensionless) More... | |
| virtual void | calc_reducing_state (void) |
| virtual CoolPropDbl | calc_Tmax (void) |
| Using this backend, calculate the maximum temperature in K. More... | |
| virtual CoolPropDbl | calc_Tmin (void) |
| Using this backend, calculate the minimum temperature in K. More... | |
| virtual CoolPropDbl | calc_pmax (void) |
| Using this backend, calculate the maximum pressure in Pa. More... | |
| virtual CoolPropDbl | calc_GWP20 (void) |
| Using this backend, calculate the 20-year global warming potential (GWP) More... | |
| virtual CoolPropDbl | calc_GWP100 (void) |
| Using this backend, calculate the 100-year global warming potential (GWP) More... | |
| virtual CoolPropDbl | calc_GWP500 (void) |
| Using this backend, calculate the 500-year global warming potential (GWP) More... | |
| virtual CoolPropDbl | calc_ODP (void) |
| Using this backend, calculate the ozone depletion potential (ODP) More... | |
| virtual CoolPropDbl | calc_flame_hazard (void) |
| Using this backend, calculate the flame hazard. More... | |
| virtual CoolPropDbl | calc_health_hazard (void) |
| Using this backend, calculate the health hazard. More... | |
| virtual CoolPropDbl | calc_physical_hazard (void) |
| Using this backend, calculate the physical hazard. More... | |
| virtual CoolPropDbl | calc_dipole_moment (void) |
| Using this backend, calculate the dipole moment in C-m (1 D = 3.33564e-30 C-m) More... | |
| virtual CoolPropDbl | calc_first_partial_deriv (parameters Of, parameters Wrt, parameters Constant) |
| Calculate the first partial derivative for the desired derivative. More... | |
| virtual CoolPropDbl | calc_second_partial_deriv (parameters Of1, parameters Wrt1, parameters Constant1, parameters Wrt2, parameters Constant2) |
| Calculate the second partial derivative using the given backend. More... | |
| virtual CoolPropDbl | calc_reduced_density (void) |
| Using this backend, calculate the reduced density (rho/rhoc) More... | |
| virtual CoolPropDbl | calc_reciprocal_reduced_temperature (void) |
| Using this backend, calculate the reciprocal reduced temperature (Tc/T) More... | |
| virtual CoolPropDbl | calc_Bvirial (void) |
| Using this backend, calculate the second virial coefficient. More... | |
| virtual CoolPropDbl | calc_Cvirial (void) |
| Using this backend, calculate the third virial coefficient. More... | |
| virtual CoolPropDbl | calc_dBvirial_dT (void) |
| Using this backend, calculate the derivative dB/dT. More... | |
| virtual CoolPropDbl | calc_dCvirial_dT (void) |
| Using this backend, calculate the derivative dC/dT. More... | |
| virtual CoolPropDbl | calc_compressibility_factor (void) |
| Using this backend, calculate the compressibility factor Z \( Z = p/(\rho R T) \). More... | |
| virtual std::string | calc_name (void) |
| Using this backend, get the name of the fluid. More... | |
| virtual std::string | calc_description (void) |
| Using this backend, get the description of the fluid. More... | |
| virtual CoolPropDbl | calc_Ttriple (void) |
| Using this backend, get the triple point temperature in K. More... | |
| virtual CoolPropDbl | calc_p_triple (void) |
| Using this backend, get the triple point pressure in Pa. More... | |
| virtual CoolPropDbl | calc_T_critical (void) |
| Using this backend, get the critical point temperature in K. More... | |
| virtual CoolPropDbl | calc_T_reducing (void) |
| Using this backend, get the reducing point temperature in K. More... | |
| virtual CoolPropDbl | calc_p_critical (void) |
| Using this backend, get the critical point pressure in Pa. More... | |
| virtual CoolPropDbl | calc_p_reducing (void) |
| Using this backend, get the reducing point pressure in Pa. More... | |
| virtual CoolPropDbl | calc_rhomolar_critical (void) |
| Using this backend, get the critical point molar density in mol/m^3. More... | |
| virtual CoolPropDbl | calc_rhomass_critical (void) |
| Using this backend, get the critical point mass density in kg/m^3 - Added for IF97Backend which is mass based. More... | |
| virtual CoolPropDbl | calc_rhomolar_reducing (void) |
| Using this backend, get the reducing point molar density in mol/m^3. More... | |
| virtual void | calc_phase_envelope (const std::string &type) |
| Using this backend, construct the phase envelope, the variable type describes the type of phase envelope to be built. More... | |
| virtual CoolPropDbl | calc_rhomass (void) |
| virtual CoolPropDbl | calc_hmass (void) |
| virtual CoolPropDbl | calc_hmass_excess (void) |
| virtual CoolPropDbl | calc_smass (void) |
| virtual CoolPropDbl | calc_smass_excess (void) |
| virtual CoolPropDbl | calc_cpmass (void) |
| virtual CoolPropDbl | calc_cp0mass (void) |
| virtual CoolPropDbl | calc_cvmass (void) |
| virtual CoolPropDbl | calc_umass (void) |
| virtual CoolPropDbl | calc_umass_excess (void) |
| virtual CoolPropDbl | calc_gibbsmass (void) |
| virtual CoolPropDbl | calc_gibbsmass_excess (void) |
| virtual CoolPropDbl | calc_helmholtzmass (void) |
| virtual CoolPropDbl | calc_helmholtzmass_excess (void) |
| virtual CoolPropDbl | calc_volumemass_excess (void) |
| virtual void | update_states (void) |
| Update the states after having changed the reference state for enthalpy and entropy. More... | |
| virtual CoolPropDbl | calc_melting_line (int param, int given, CoolPropDbl value) |
| virtual CoolPropDbl | calc_saturation_ancillary (parameters param, int Q, parameters given, double value) |
| virtual phases | calc_phase (void) |
| Using this backend, calculate the phase. More... | |
| virtual void | calc_specify_phase (phases phase) |
| Using this backend, specify the phase to be used for all further calculations. More... | |
| virtual void | calc_unspecify_phase (void) |
| Using this backend, unspecify the phase. More... | |
| virtual std::vector< std::string > | calc_fluid_names (void) |
| Using this backend, get a vector of fluid names. More... | |
| virtual const CoolProp::SimpleState & | calc_state (const std::string &state) |
| virtual const CoolProp::PhaseEnvelopeData & | calc_phase_envelope_data (void) |
| virtual std::vector< CoolPropDbl > | calc_mole_fractions_liquid (void) |
| virtual std::vector< CoolPropDbl > | calc_mole_fractions_vapor (void) |
| virtual const std::vector< CoolPropDbl > | calc_mass_fractions (void) |
| virtual CoolPropDbl | calc_fraction_min (void) |
| Get the minimum fraction (mole, mass, volume) for incompressible fluid. More... | |
| virtual CoolPropDbl | calc_fraction_max (void) |
| Get the maximum fraction (mole, mass, volume) for incompressible fluid. More... | |
| virtual CoolPropDbl | calc_T_freeze (void) |
| virtual CoolPropDbl | calc_first_saturation_deriv (parameters Of1, parameters Wrt1) |
| virtual CoolPropDbl | calc_second_saturation_deriv (parameters Of1, parameters Wrt1, parameters Wrt2) |
| virtual CoolPropDbl | calc_first_two_phase_deriv (parameters Of, parameters Wrt, parameters Constant) |
| virtual CoolPropDbl | calc_second_two_phase_deriv (parameters Of, parameters Wrt, parameters Constant, parameters Wrt2, parameters Constant2) |
| virtual CoolPropDbl | calc_first_two_phase_deriv_splined (parameters Of, parameters Wrt, parameters Constant, CoolPropDbl x_end) |
| virtual CoolPropDbl | calc_saturated_liquid_keyed_output (parameters key) |
| virtual CoolPropDbl | calc_saturated_vapor_keyed_output (parameters key) |
| virtual void | calc_ideal_curve (const std::string &type, std::vector< double > &T, std::vector< double > &p) |
| virtual CoolPropDbl | calc_T (void) |
| Using this backend, get the temperature. More... | |
| virtual CoolPropDbl | calc_rhomolar (void) |
| Using this backend, get the molar density in mol/m^3. More... | |
| virtual double | calc_tangent_plane_distance (const double T, const double p, const std::vector< double > &w, const double rhomolar_guess) |
| Using this backend, calculate the tangent plane distance for a given trial composition. More... | |
| virtual void | calc_true_critical_point (double &T, double &rho) |
| Using this backend, return true critical point where dp/drho|T = 0 and d2p/drho^2|T = 0. More... | |
| virtual void | calc_conformal_state (const std::string &reference_fluid, CoolPropDbl &T, CoolPropDbl &rhomolar) |
| virtual void | calc_viscosity_contributions (CoolPropDbl &dilute, CoolPropDbl &initial_density, CoolPropDbl &residual, CoolPropDbl &critical) |
| virtual void | calc_conductivity_contributions (CoolPropDbl &dilute, CoolPropDbl &initial_density, CoolPropDbl &residual, CoolPropDbl &critical) |
| virtual std::vector< CriticalState > | calc_all_critical_points (void) |
| virtual void | calc_build_spinodal () |
| virtual SpinodalData | calc_get_spinodal_data () |
| virtual void | calc_criticality_contour_values (double &L1star, double &M1star) |
| virtual void | mass_to_molar_inputs (CoolProp::input_pairs &input_pair, CoolPropDbl &value1, CoolPropDbl &value2) |
| Convert mass-based input pair to molar-based input pair; If molar-based, do nothing. More... | |
| virtual void | calc_change_EOS (const std::size_t i, const std::string &EOS_name) |
| Change the equation of state for a given component to a specified EOS. More... | |
Protected Attributes | |
| long | _fluid_type |
| Some administrative variables. More... | |
| phases | _phase |
| The key for the phase from CoolProp::phases enum. More... | |
| phases | imposed_phase_index |
| If the phase is imposed, the imposed phase index. More... | |
| CacheArray< 70 > | cache |
| SimpleState | _critical |
| Two important points. More... | |
| SimpleState | _reducing |
| CAE | _molar_mass = cache.next() |
| Molar mass [mol/kg]. More... | |
| CAE | _gas_constant = cache.next() |
| Universal gas constant [J/mol/K]. More... | |
| double | _rhomolar |
| Bulk values. More... | |
| double | _T |
| double | _p |
| double | _Q |
| CAE | _tau = cache.next() |
| CAE | _delta = cache.next() |
| CAE | _viscosity = cache.next() |
| Transport properties. More... | |
| CAE | _conductivity = cache.next() |
| CAE | _surface_tension = cache.next() |
| CAE | _hmolar = cache.next() |
| CAE | _smolar = cache.next() |
| CAE | _umolar = cache.next() |
| CAE | _logp = cache.next() |
| CAE | _logrhomolar = cache.next() |
| CAE | _cpmolar = cache.next() |
| CAE | _cp0molar = cache.next() |
| CAE | _cvmolar = cache.next() |
| CAE | _speed_sound = cache.next() |
| CAE | _gibbsmolar = cache.next() |
| CAE | _helmholtzmolar = cache.next() |
| CAE | _hmolar_residual = cache.next() |
| Residual properties. More... | |
| CAE | _smolar_residual = cache.next() |
| CAE | _gibbsmolar_residual = cache.next() |
| CAE | _hmolar_excess = cache.next() |
| Excess properties. More... | |
| CAE | _smolar_excess = cache.next() |
| CAE | _gibbsmolar_excess = cache.next() |
| CAE | _umolar_excess = cache.next() |
| CAE | _volumemolar_excess = cache.next() |
| CAE | _helmholtzmolar_excess = cache.next() |
| CAE | _rhoLanc = cache.next() |
| Ancillary values. More... | |
| CAE | _rhoVanc = cache.next() |
| CAE | _pLanc = cache.next() |
| CAE | _pVanc = cache.next() |
| CAE | _TLanc = cache.next() |
| CAE | _TVanc = cache.next() |
| CachedElement | _fugacity_coefficient |
| CAE | _rho_spline = cache.next() |
| Smoothing values. More... | |
| CAE | _drho_spline_dh__constp = cache.next() |
| CAE | _drho_spline_dp__consth = cache.next() |
| CAE | _alpha0 = cache.next() |
| Cached low-level elements for in-place calculation of other properties. More... | |
| CAE | _dalpha0_dTau = cache.next() |
| CAE | _dalpha0_dDelta = cache.next() |
| CAE | _d2alpha0_dTau2 = cache.next() |
| CAE | _d2alpha0_dDelta_dTau = cache.next() |
| CAE | _d2alpha0_dDelta2 = cache.next() |
| CAE | _d3alpha0_dTau3 = cache.next() |
| CAE | _d3alpha0_dDelta_dTau2 = cache.next() |
| CAE | _d3alpha0_dDelta2_dTau = cache.next() |
| CAE | _d3alpha0_dDelta3 = cache.next() |
| CAE | _alphar = cache.next() |
| CAE | _dalphar_dTau = cache.next() |
| CAE | _dalphar_dDelta = cache.next() |
| CAE | _d2alphar_dTau2 = cache.next() |
| CAE | _d2alphar_dDelta_dTau = cache.next() |
| CAE | _d2alphar_dDelta2 = cache.next() |
| CAE | _d3alphar_dTau3 = cache.next() |
| CAE | _d3alphar_dDelta_dTau2 = cache.next() |
| CAE | _d3alphar_dDelta2_dTau = cache.next() |
| CAE | _d3alphar_dDelta3 = cache.next() |
| CAE | _d4alphar_dTau4 = cache.next() |
| CAE | _d4alphar_dDelta_dTau3 = cache.next() |
| CAE | _d4alphar_dDelta2_dTau2 = cache.next() |
| CAE | _d4alphar_dDelta3_dTau = cache.next() |
| CAE | _d4alphar_dDelta4 = cache.next() |
| CAE | _dalphar_dDelta_lim = cache.next() |
| CAE | _d2alphar_dDelta2_lim = cache.next() |
| CAE | _d2alphar_dDelta_dTau_lim = cache.next() |
| CAE | _d3alphar_dDelta2_dTau_lim = cache.next() |
| CAE | _rhoLmolar = cache.next() |
| Two-Phase variables. More... | |
| CAE | _rhoVmolar = cache.next() |
|
protected |
Definition at line 99 of file AbstractState.h.
|
inline |
Definition at line 704 of file AbstractState.h.
|
inlinevirtual |
Definition at line 707 of file AbstractState.h.
| double CoolProp::AbstractState::acentric_factor | ( | void | ) |
Return the acentric factor.
Definition at line 641 of file AbstractState.cpp.
|
inline |
Return the vector of critical points, including points that are unstable or correspond to negative pressure.
Definition at line 978 of file AbstractState.h.
|
inline |
Return the term \( \alpha^0 \).
Definition at line 1487 of file AbstractState.h.
|
inline |
Return the term \( \alpha^r \).
Definition at line 1538 of file AbstractState.h.
|
inlinevirtual |
Apply a simple mixing rule (EXPERT USE ONLY!!!)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 913 of file AbstractState.h.
|
inlinevirtual |
A function that says whether the backend instance can be instantiated in the high-level interface In general this should be true, except for some other backends (especially the tabular backends) To disable use in high-level interface, implement this function and return false
Reimplemented in CoolProp::TabularBackend.
Definition at line 863 of file AbstractState.h.
|
pure virtual |
Get a string representation of the backend - for instance "HelmholtzEOSMixtureBackend" for the core mixture model in CoolProp
Must be overloaded by the backend to provide the backend's name
Implemented in CoolProp::SRKBackend, CoolProp::PengRobinsonBackend, CoolProp::VTPRBackend, CoolProp::HelmholtzEOSBackend, CoolProp::HelmholtzEOSMixtureBackend, CoolProp::IF97Backend, CoolProp::IncompressibleBackend, CoolProp::PCSAFTBackend, CoolProp::REFPROPBackend, CoolProp::REFPROPMixtureBackend, CoolProp::BicubicBackend, and CoolProp::TTSEBackend.
| void CoolProp::AbstractState::build_phase_envelope | ( | const std::string & | type = "" | ) |
Construct the phase envelope for a mixture.
| type | currently a dummy variable that is not used |
Definition at line 675 of file AbstractState.cpp.
|
inline |
Construct the spinodal curve for the mixture (or pure fluid)
Definition at line 983 of file AbstractState.h.
| double CoolProp::AbstractState::Bvirial | ( | void | ) |
Return the B virial coefficient.
Definition at line 687 of file AbstractState.cpp.
|
inlineprotectedvirtual |
Using this backend, calculate the acentric factor.
Reimplemented in CoolProp::AbstractCubicBackend, CoolProp::HelmholtzEOSMixtureBackend, CoolProp::IF97Backend, and CoolProp::REFPROPMixtureBackend.
Definition at line 237 of file AbstractState.h.
|
inlineprotectedvirtual |
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 682 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0\) (dimensionless)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 340 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r\) (dimensionless)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, CoolProp::PCSAFTBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 277 of file AbstractState.h.
|
inlineprotectedvirtual |
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 685 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the second virial coefficient.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 445 of file AbstractState.h.
|
inlineprotectedvirtual |
Change the equation of state for a given component to a specified EOS.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 699 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the chemical potential in J/mol.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 261 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the compressibility factor Z \( Z = p/(\rho R T) \).
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, CoolProp::PCSAFTBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 461 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the thermal conductivity in W/m/K.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, CoolProp::IF97Backend, CoolProp::IncompressibleBackend, CoolProp::REFPROPMixtureBackend, and CoolProp::TabularBackend.
Definition at line 225 of file AbstractState.h.
|
inlineprotectedvirtual |
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 679 of file AbstractState.h.
|
inlineprotectedvirtual |
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 672 of file AbstractState.h.
|
inlineprotectedvirtual |
Definition at line 534 of file AbstractState.h.
|
inlineprotectedvirtual |
Reimplemented in CoolProp::IF97Backend, and CoolProp::IncompressibleBackend.
Definition at line 531 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the molar constant-pressure specific heat in J/mol/K.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, CoolProp::IF97Backend, and CoolProp::TabularBackend.
Definition at line 181 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the ideal gas molar constant-pressure specific heat in J/mol/K.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, CoolProp::REFPROPMixtureBackend, and CoolProp::TabularBackend.
Definition at line 185 of file AbstractState.h.
|
inlineprotectedvirtual |
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 691 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the third virial coefficient.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 449 of file AbstractState.h.
|
inlineprotectedvirtual |
Reimplemented in CoolProp::IF97Backend, and CoolProp::IncompressibleBackend.
Definition at line 537 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the molar constant-volume specific heat in J/mol/K.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, CoolProp::IF97Backend, and CoolProp::TabularBackend.
Definition at line 189 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\delta\tau}\) (dimensionless)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 356 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\delta\delta}\) (dimensionless)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 352 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\tau\tau}\) (dimensionless)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 360 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\delta}\) (dimensionless)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 289 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\tau}\) (dimensionless)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 293 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\tau\tau}\) (dimensionless)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 297 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\delta\delta\tau}\) (dimensionless)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 368 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\delta\delta\delta}\) (dimensionless)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 364 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\delta\tau\tau}\) (dimensionless)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 372 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\tau\tau\tau}\) (dimensionless)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 376 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\delta\tau}\) (dimensionless)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 305 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\delta\delta}\) (dimensionless)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 301 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\tau\tau}\) (dimensionless)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 309 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\tau\tau\tau}\) (dimensionless)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 313 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\delta\tau\tau}\) (dimensionless)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 326 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\delta\delta\tau}\) (dimensionless)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 322 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\delta\delta\delta}\) (dimensionless)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 318 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\tau\tau\tau}\) (dimensionless)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 330 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\tau\tau\tau\tau}\) (dimensionless)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 334 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\delta}\) (dimensionless)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 344 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\tau}\) (dimensionless)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 348 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta}\) (dimensionless)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 281 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\tau}\) (dimensionless)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 285 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the derivative dB/dT.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 453 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the derivative dC/dT.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 457 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, get the description of the fluid.
Reimplemented in CoolProp::IncompressibleBackend.
Definition at line 470 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the dipole moment in C-m (1 D = 3.33564e-30 C-m)
Reimplemented in CoolProp::REFPROPMixtureBackend.
Definition at line 426 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate and cache the excess properties.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 271 of file AbstractState.h.
|
protectedvirtual |
Calculate the first partial derivative for the desired derivative.
Reimplemented in CoolProp::IncompressibleBackend, and CoolProp::TabularBackend.
Definition at line 967 of file AbstractState.cpp.
|
inlineprotectedvirtual |
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::TabularBackend.
Definition at line 627 of file AbstractState.h.
|
inlineprotectedvirtual |
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::TabularBackend.
Definition at line 633 of file AbstractState.h.
|
inlineprotectedvirtual |
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::TabularBackend.
Definition at line 639 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the flame hazard.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 414 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, get a vector of fluid names.
Reimplemented in CoolProp::REFPROPMixtureBackend, CoolProp::AbstractCubicBackend, CoolProp::VTPRBackend, CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::TabularBackend.
Definition at line 592 of file AbstractState.h.
|
inlineprotectedvirtual |
Get the maximum fraction (mole, mass, volume) for incompressible fluid.
Reimplemented in CoolProp::IncompressibleBackend.
Definition at line 620 of file AbstractState.h.
|
inlineprotectedvirtual |
Get the minimum fraction (mole, mass, volume) for incompressible fluid.
Reimplemented in CoolProp::IncompressibleBackend.
Definition at line 616 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the fugacity in Pa.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 257 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the fugacity coefficient (dimensionless)
Reimplemented in CoolProp::VTPRBackend, CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 249 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the fugacity in Pa.
Reimplemented in CoolProp::PCSAFTBackend.
Definition at line 253 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the universal gas constant \(R_u\) in J/mol/K.
Reimplemented in CoolProp::AbstractCubicBackend, CoolProp::HelmholtzEOSMixtureBackend, CoolProp::IF97Backend, and CoolProp::REFPROPMixtureBackend.
Definition at line 245 of file AbstractState.h.
|
inlineprotectedvirtual |
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 688 of file AbstractState.h.
|
inlineprotectedvirtual |
Definition at line 546 of file AbstractState.h.
|
inlineprotectedvirtual |
Definition at line 549 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the molar Gibbs function in J/mol.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 193 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the residual molar Gibbs function in J/mol.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::PCSAFTBackend.
Definition at line 197 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the 100-year global warming potential (GWP)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 402 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the 20-year global warming potential (GWP)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 398 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the 500-year global warming potential (GWP)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 406 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the health hazard.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 418 of file AbstractState.h.
|
inlineprotectedvirtual |
Definition at line 552 of file AbstractState.h.
|
inlineprotectedvirtual |
Definition at line 555 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the molar Helmholtz energy in J/mol.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 201 of file AbstractState.h.
|
inlineprotectedvirtual |
Reimplemented in CoolProp::IF97Backend, and CoolProp::IncompressibleBackend.
Definition at line 519 of file AbstractState.h.
|
inlineprotectedvirtual |
Definition at line 522 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the molar enthalpy in J/mol.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, CoolProp::IF97Backend, and CoolProp::TabularBackend.
Definition at line 157 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the residual molar enthalpy in J/mol.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::PCSAFTBackend.
Definition at line 161 of file AbstractState.h.
|
inlineprotectedvirtual |
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 649 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the isentropic expansion coefficient \( \kappa_s = -\frac{c_p}{c_v}\frac{v}{p}\left.\frac{\partial p}{\partial v}\right|_T = \frac{\rho}{p}\left.\frac{\partial p}{\partial \rho}\right|_s\).
Definition at line 217 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the isobaric expansion coefficient \( \beta = \frac{1}{v}\left.\frac{\partial v}{\partial T}\right|_p = -\frac{1}{\rho}\left.\frac{\partial \rho}{\partial T}\right|_p\) in 1/K.
Definition at line 213 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the isothermal compressibility \( \kappa = -\frac{1}{v}\left.\frac{\partial v}{\partial p}\right|_T=\frac{1}{\rho}\left.\frac{\partial \rho}{\partial p}\right|_T\) in 1/Pa.
Definition at line 209 of file AbstractState.h.
|
inlineprotectedvirtual |
Reimplemented in CoolProp::REFPROPMixtureBackend, CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::TabularBackend.
Definition at line 611 of file AbstractState.h.
|
inlineprotectedvirtual |
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, CoolProp::IncompressibleBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 567 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the molar mass in kg/mol.
Reimplemented in CoolProp::AbstractCubicBackend, CoolProp::VTPRBackend, CoolProp::HelmholtzEOSMixtureBackend, CoolProp::IF97Backend, CoolProp::PCSAFTBackend, CoolProp::REFPROPMixtureBackend, and CoolProp::TabularBackend.
Definition at line 233 of file AbstractState.h.
|
inlineprotectedvirtual |
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 605 of file AbstractState.h.
|
inlineprotectedvirtual |
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 608 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, get the name of the fluid.
Reimplemented in CoolProp::REFPROPMixtureBackend, CoolProp::HelmholtzEOSMixtureBackend, CoolProp::IncompressibleBackend, and CoolProp::TabularBackend.
Definition at line 466 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate effective hardness of interaction.
Definition at line 173 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the ozone depletion potential (ODP)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 410 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, get the critical point pressure in Pa.
Reimplemented in CoolProp::AbstractCubicBackend, CoolProp::HelmholtzEOSMixtureBackend, CoolProp::IF97Backend, CoolProp::REFPROPMixtureBackend, and CoolProp::TabularBackend.
Definition at line 492 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, get the reducing point pressure in Pa.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 496 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, get the triple point pressure in Pa.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, CoolProp::IF97Backend, CoolProp::REFPROPMixtureBackend, and CoolProp::TabularBackend.
Definition at line 479 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the phase.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, CoolProp::IF97Backend, CoolProp::PCSAFTBackend, CoolProp::REFPROPMixtureBackend, and CoolProp::TabularBackend.
Definition at line 580 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, construct the phase envelope, the variable type describes the type of phase envelope to be built.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 512 of file AbstractState.h.
|
inlineprotectedvirtual |
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 601 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the physical hazard.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 422 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the phase identification parameter (PIP)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 265 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the maximum pressure in Pa.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, CoolProp::IF97Backend, CoolProp::REFPROPMixtureBackend, and CoolProp::TabularBackend.
Definition at line 393 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the pressure in Pa.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, CoolProp::IF97Backend, and CoolProp::PCSAFTBackend.
Definition at line 241 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the reciprocal reduced temperature (Tc/T)
Reimplemented in CoolProp::AbstractCubicBackend.
Definition at line 440 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the reduced density (rho/rhoc)
Reimplemented in CoolProp::AbstractCubicBackend.
Definition at line 436 of file AbstractState.h.
|
inlineprotectedvirtual |
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 380 of file AbstractState.h.
|
inlineprotectedvirtual |
Reimplemented in CoolProp::IF97Backend, and CoolProp::IncompressibleBackend.
Definition at line 516 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, get the critical point mass density in kg/m^3 - Added for IF97Backend which is mass based.
Reimplemented in CoolProp::IF97Backend.
Definition at line 504 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, get the molar density in mol/m^3.
Reimplemented in CoolProp::IF97Backend, and CoolProp::TabularBackend.
Definition at line 658 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, get the critical point molar density in mol/m^3.
Reimplemented in CoolProp::AbstractCubicBackend, CoolProp::HelmholtzEOSMixtureBackend, CoolProp::IF97Backend, CoolProp::REFPROPMixtureBackend, and CoolProp::TabularBackend.
Definition at line 500 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, get the reducing point molar density in mol/m^3.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 508 of file AbstractState.h.
|
inlineprotectedvirtual |
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, CoolProp::REFPROPMixtureBackend, and CoolProp::TabularBackend.
Definition at line 643 of file AbstractState.h.
|
inlineprotectedvirtual |
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, CoolProp::REFPROPMixtureBackend, and CoolProp::TabularBackend.
Definition at line 646 of file AbstractState.h.
|
inlineprotectedvirtual |
| param | The key for the parameter to be returned |
| Q | The quality for the parameter that is given (0 = saturated liquid, 1 = saturated vapor) |
| given | The key for the parameter that is given |
| value | The value for the parameter that is given |
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 575 of file AbstractState.h.
|
protectedvirtual |
Calculate the second partial derivative using the given backend.
Definition at line 976 of file AbstractState.cpp.
|
inlineprotectedvirtual |
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 630 of file AbstractState.h.
|
inlineprotectedvirtual |
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 636 of file AbstractState.h.
|
inlineprotectedvirtual |
Reimplemented in CoolProp::IF97Backend, and CoolProp::IncompressibleBackend.
Definition at line 525 of file AbstractState.h.
|
inlineprotectedvirtual |
Definition at line 528 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the molar entropy in J/mol/K.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, CoolProp::IF97Backend, and CoolProp::TabularBackend.
Definition at line 165 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the residual molar entropy in J/mol/K.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, CoolProp::PCSAFTBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 169 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, specify the phase to be used for all further calculations.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, CoolProp::PCSAFTBackend, CoolProp::REFPROPMixtureBackend, and CoolProp::TabularBackend.
Definition at line 584 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the speed of sound in m/s.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, CoolProp::IF97Backend, and CoolProp::TabularBackend.
Definition at line 205 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate a phase given by the state string
| state | A string that describes the state desired, one of "hs_anchor", "critical"/"crit", "reducing" |
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 597 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the surface tension in N/m.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, CoolProp::IF97Backend, CoolProp::REFPROPMixtureBackend, and CoolProp::TabularBackend.
Definition at line 229 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, get the temperature.
Reimplemented in CoolProp::TabularBackend.
Definition at line 654 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, get the critical point temperature in K.
Reimplemented in CoolProp::AbstractCubicBackend, CoolProp::HelmholtzEOSMixtureBackend, CoolProp::IF97Backend, CoolProp::REFPROPMixtureBackend, and CoolProp::TabularBackend.
Definition at line 484 of file AbstractState.h.
|
inlineprotectedvirtual |
Reimplemented in CoolProp::IncompressibleBackend.
Definition at line 623 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, get the reducing point temperature in K.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 488 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the tangent plane distance for a given trial composition.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 663 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the maximum temperature in K.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, CoolProp::IF97Backend, CoolProp::IncompressibleBackend, CoolProp::REFPROPMixtureBackend, and CoolProp::TabularBackend.
Definition at line 385 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the minimum temperature in K.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, CoolProp::IF97Backend, CoolProp::IncompressibleBackend, CoolProp::REFPROPMixtureBackend, and CoolProp::TabularBackend.
Definition at line 389 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, return true critical point where dp/drho|T = 0 and d2p/drho^2|T = 0.
Reimplemented in CoolProp::REFPROPMixtureBackend.
Definition at line 668 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, get the triple point temperature in K.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, CoolProp::IF97Backend, CoolProp::REFPROPMixtureBackend, and CoolProp::TabularBackend.
Definition at line 475 of file AbstractState.h.
|
inlineprotectedvirtual |
Reimplemented in CoolProp::IF97Backend, and CoolProp::IncompressibleBackend.
Definition at line 540 of file AbstractState.h.
|
inlineprotectedvirtual |
Definition at line 543 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the molar internal energy in J/mol.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, CoolProp::IF97Backend, and CoolProp::TabularBackend.
Definition at line 177 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, unspecify the phase.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, CoolProp::PCSAFTBackend, CoolProp::REFPROPMixtureBackend, and CoolProp::TabularBackend.
Definition at line 588 of file AbstractState.h.
|
inlineprotectedvirtual |
Using this backend, calculate the viscosity in Pa-s.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, CoolProp::IF97Backend, CoolProp::IncompressibleBackend, CoolProp::REFPROPMixtureBackend, and CoolProp::TabularBackend.
Definition at line 221 of file AbstractState.h.
|
inlineprotectedvirtual |
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 676 of file AbstractState.h.
|
inlineprotectedvirtual |
Definition at line 558 of file AbstractState.h.
|
inline |
Change the equation of state for a given component to a specified EOS.
| i | Index of the component to change (if a pure fluid, i=0) |
| EOS_name | Name of the EOS to use (something like "SRK", "PR", "XiangDeiters", but backend-specific) |
Definition at line 1479 of file AbstractState.h.
| double CoolProp::AbstractState::chemical_potential | ( | std::size_t | i | ) |
Return the chemical potential of the i-th component of the mixture.
Definition at line 671 of file AbstractState.cpp.
|
virtual |
Clear all the cached values.
Bulk values
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, CoolProp::IF97Backend, and CoolProp::IncompressibleBackend.
Definition at line 186 of file AbstractState.cpp.
|
virtual |
When the composition changes, clear all cached values that are only dependent on composition, but not the thermodynamic state.
Definition at line 176 of file AbstractState.cpp.
| double CoolProp::AbstractState::compressibility_factor | ( | void | ) |
Return the compressibility factor \( Z = p/(rho R T) \).
Definition at line 699 of file AbstractState.cpp.
| double CoolProp::AbstractState::conductivity | ( | void | ) |
Return the thermal conductivity in W/m/K.
Definition at line 634 of file AbstractState.cpp.
|
inline |
Return the thermal conductivity contributions, each in W/m/K.
Definition at line 1456 of file AbstractState.h.
|
inline |
Find the conformal state needed for ECS.
| reference_fluid | The reference fluid for which the conformal state will be calculated |
| T | Temperature (initial guess must be provided, or < 0 to start with unity shape factors) |
| rhomolar | Molar density (initial guess must be provided, or < 0 to start with unity shape factors) |
Definition at line 1471 of file AbstractState.h.
|
inline |
Return the mass constant pressure specific heat for ideal gas part only in J/kg/K.
Definition at line 1168 of file AbstractState.h.
| double CoolProp::AbstractState::cp0molar | ( | void | ) |
Return the molar constant pressure specific heat for ideal gas part only in J/mol/K.
Definition at line 619 of file AbstractState.cpp.
|
inline |
Return the mass constant pressure specific heat in J/kg/K.
Definition at line 1162 of file AbstractState.h.
| double CoolProp::AbstractState::cpmolar | ( | void | ) |
Return the molar constant pressure specific heat in J/mol/K.
Definition at line 615 of file AbstractState.cpp.
|
inline |
Calculate the criticality contour values \(\mathcal{L}_1^*\) and \(\mathcal{M}_1^*\).
Definition at line 993 of file AbstractState.h.
| double CoolProp::AbstractState::Cvirial | ( | void | ) |
Return the C virial coefficient.
Definition at line 690 of file AbstractState.cpp.
|
inline |
Return the mass constant volume specific heat in J/kg/K.
Definition at line 1174 of file AbstractState.h.
| double CoolProp::AbstractState::cvmolar | ( | void | ) |
Return the molar constant volume specific heat in J/mol/K.
Definition at line 622 of file AbstractState.cpp.
|
inline |
Return the term \( \alpha^0_{\delta\delta} \).
Definition at line 1502 of file AbstractState.h.
|
inline |
Return the term \( \alpha^0_{\delta\tau} \).
Definition at line 1507 of file AbstractState.h.
|
inline |
Return the term \( \alpha^0_{\tau\tau} \).
Definition at line 1512 of file AbstractState.h.
|
inline |
Return the term \( \alpha^r_{\delta\delta} \).
Definition at line 1553 of file AbstractState.h.
|
inline |
Return the term \( \alpha^r_{\delta\tau} \).
Definition at line 1558 of file AbstractState.h.
|
inline |
Return the term \( \alpha^r_{\tau\tau} \).
Definition at line 1563 of file AbstractState.h.
|
inline |
Return the term \( \alpha^0_{\delta\delta\tau} \).
Definition at line 1527 of file AbstractState.h.
|
inline |
Return the term \( \alpha^0_{\delta\delta\delta} \).
Definition at line 1532 of file AbstractState.h.
|
inline |
Return the term \( \alpha^0_{\delta\tau\tau} \).
Definition at line 1522 of file AbstractState.h.
|
inline |
Return the term \( \alpha^0_{\tau\tau\tau} \).
Definition at line 1517 of file AbstractState.h.
|
inline |
Return the term \( \alpha^r_{\delta\delta\tau} \).
Definition at line 1573 of file AbstractState.h.
|
inline |
Return the term \( \alpha^r_{\delta\delta\delta} \).
Definition at line 1568 of file AbstractState.h.
|
inline |
Return the term \( \alpha^r_{\delta\tau\tau} \).
Definition at line 1578 of file AbstractState.h.
|
inline |
Return the term \( \alpha^r_{\tau\tau\tau} \).
Definition at line 1583 of file AbstractState.h.
|
inline |
Return the term \( \alpha^r_{\delta\delta\tau\tau} \).
Definition at line 1598 of file AbstractState.h.
|
inline |
Return the term \( \alpha^r_{\delta\delta\delta\tau} \).
Definition at line 1593 of file AbstractState.h.
|
inline |
Return the term \( \alpha^r_{\delta\delta\delta\delta} \).
Definition at line 1588 of file AbstractState.h.
|
inline |
Return the term \( \alpha^r_{\delta\tau\tau\tau} \).
Definition at line 1603 of file AbstractState.h.
|
inline |
Return the term \( \alpha^r_{\tau\tau\tau\tau} \).
Definition at line 1608 of file AbstractState.h.
|
inline |
Return the term \( \alpha^0_{\delta} \).
Definition at line 1492 of file AbstractState.h.
|
inline |
Return the term \( \alpha^0_{\tau} \).
Definition at line 1497 of file AbstractState.h.
|
inline |
Return the term \( \alpha^r_{\delta} \).
Definition at line 1543 of file AbstractState.h.
|
inline |
Return the term \( \alpha^r_{\tau} \).
Definition at line 1548 of file AbstractState.h.
| double CoolProp::AbstractState::dBvirial_dT | ( | void | ) |
Return the derivative of the B virial coefficient with respect to temperature.
Definition at line 693 of file AbstractState.cpp.
| double CoolProp::AbstractState::dCvirial_dT | ( | void | ) |
Return the derivative of the C virial coefficient with respect to temperature.
Definition at line 696 of file AbstractState.cpp.
| double CoolProp::AbstractState::delta | ( | void | ) |
Return the reduced density ( \(\delta = \rho/\rho_c\))
Definition at line 487 of file AbstractState.cpp.
|
inline |
Return the description - backend dependent.
Definition at line 1039 of file AbstractState.h.
|
inline |
Return the dipole moment in C-m (1 D = 3.33564e-30 C-m)
Definition at line 1044 of file AbstractState.h.
|
inlinestatic |
A factory function to return a pointer to a new-allocated instance of one of the backends.
This is a convenience function to allow for the use of '&' delimited fluid names. Slightly less computationally efficient than the
| backend | The backend in use, one of "HEOS", "REFPROP", etc. |
| fluid_names | Fluid names as a '&' delimited string |
Definition at line 716 of file AbstractState.h.
|
static |
A factory function to return a pointer to a new-allocated instance of one of the backends.
| backend | The backend in use, "HEOS", "REFPROP", etc. |
| fluid_names | A vector of strings of the fluid names |
Several backends are possible:
Very Important!! : Use a smart pointer to manage the pointer returned. In older versions of C++, you can use std::tr1::smart_ptr. In C++2011 you can use std::shared_ptr
Definition at line 126 of file AbstractState.cpp.
|
inline |
The first partial derivative in homogeneous phases.
\[ \left(\frac{\partial A}{\partial B}\right)_C = \frac{\left(\frac{\partial A}{\partial \tau}\right)_\delta\left(\frac{\partial C}{\partial \delta}\right)_\tau-\left(\frac{\partial A}{\partial \delta}\right)_\tau\left(\frac{\partial C}{\partial \tau}\right)_\delta}{\left(\frac{\partial B}{\partial \tau}\right)_\delta\left(\frac{\partial C}{\partial \delta}\right)_\tau-\left(\frac{\partial B}{\partial \delta}\right)_\tau\left(\frac{\partial C}{\partial \tau}\right)_\delta} = \frac{N}{D}\]
Definition at line 1263 of file AbstractState.h.
|
inline |
The first partial derivative along the saturation curve.
Implementing the algorithms and ideas of: Matthis Thorade, Ali Saadat, "Partial derivatives of thermodynamic state properties for dynamic simulation", Environmental Earth Sciences, December 2013, Volume 70, Issue 8, pp 3497-3503
Basically the idea is that the p-T derivative is given by Clapeyron relations:
\[ \left(\frac{\partial T}{\partial p}\right)_{\sigma} = T\left(\frac{v'' - v'}{h'' - h'}\right)_{\sigma} \]
and then other derivatives can be obtained along the saturation curve from
\[ \left(\frac{\partial y}{\partial p}\right)_{\sigma} = \left(\frac{\partial y}{\partial p}\right)+\left(\frac{\partial y}{\partial T}\right)\left(\frac{\partial T}{\partial p}\right)_{\sigma} \]
\[ \left(\frac{\partial y}{\partial T}\right)_{\sigma} = \left(\frac{\partial y}{\partial T}\right)+\left(\frac{\partial y}{\partial p}\right)\left(\frac{\partial p}{\partial T}\right)_{\sigma} \]
where derivatives without the \( \sigma \) are homogeneous (conventional) derivatives.
| Of1 | The parameter that the derivative is taken of |
| Wrt1 | The parameter that the derivative is taken with respect to |
Definition at line 1315 of file AbstractState.h.
|
inline |
Calculate the first "two-phase" derivative as described by Thorade and Sadaat, EAS, 2013.
Implementing the algorithms and ideas of: Matthis Thorade, Ali Saadat, "Partial derivatives of thermodynamic state properties for dynamic simulation", Environmental Earth Sciences, December 2013, Volume 70, Issue 8, pp 3497-3503
Spline evaluation is as described in: S Quoilin, I Bell, A Desideri, P Dewallef, V Lemort, "Methods to increase the robustness of finite-volume flow models in thermodynamic systems", Energies 7 (3), 1621-1640
| Of | The parameter to be derived |
| Wrt | The parameter that the derivative is taken with respect to |
| Constant | The parameter that is held constant |
Definition at line 1359 of file AbstractState.h.
|
inline |
Calculate the first "two-phase" derivative as described by Thorade and Sadaat, EAS, 2013.
Implementing the algorithms and ideas of: Matthis Thorade, Ali Saadat, "Partial derivatives of thermodynamic state properties for dynamic simulation", Environmental Earth Sciences, December 2013, Volume 70, Issue 8, pp 3497-3503
Spline evaluation is as described in: S Quoilin, I Bell, A Desideri, P Dewallef, V Lemort, "Methods to increase the robustness of finite-volume flow models in thermodynamic systems", Energies 7 (3), 1621-1640
| Of | The parameter to be derived |
| Wrt | The parameter that the derivative is taken with respect to |
| Constant | The parameter that is held constant |
| x_end | The end vapor quality at which the spline is defined (spline is active in [0, x_end]) |
Definition at line 1403 of file AbstractState.h.
| std::vector< std::string > CoolProp::AbstractState::fluid_names | ( | void | ) |
Return a vector of strings of the fluid names that are in use.
Definition at line 173 of file AbstractState.cpp.
|
inlinevirtual |
Return a string from the backend for the mixture/fluid - backend dependent - could be CAS #, name, etc.
Reimplemented in CoolProp::AbstractCubicBackend, CoolProp::HelmholtzEOSMixtureBackend, CoolProp::IncompressibleBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 868 of file AbstractState.h.
| double CoolProp::AbstractState::fugacity | ( | std::size_t | i | ) |
Return the fugacity of the i-th component of the mixture.
Definition at line 667 of file AbstractState.cpp.
| double CoolProp::AbstractState::fugacity_coefficient | ( | std::size_t | i | ) |
Return the fugacity coefficient of the i-th component of the mixture.
Definition at line 659 of file AbstractState.cpp.
| std::vector< double > CoolProp::AbstractState::fugacity_coefficients | ( | ) |
Return a vector of the fugacity coefficients for all components in the mixture.
Definition at line 663 of file AbstractState.cpp.
| double CoolProp::AbstractState::fundamental_derivative_of_gas_dynamics | ( | void | ) |
Return the fundamental derivative of gas dynamics \( \Gamma \).
see also Colonna et al, FPE, 2010
\[ \Gamma = 1+\frac{\rho}{c}\left(\frac{\partial c}{\partial \rho}\right)_{s} = 1+\frac{\rho}{2c^2}\left(\frac{\partial^2 p}{\partial \rho^2}\right)_{s} = \frac{v^3}{2c^2}\left(\frac{\partial^2 p}{\partial v^2}\right)_{s}\]
Note: densities are mass-based densities, not mole-based densities
Definition at line 703 of file AbstractState.cpp.
| double CoolProp::AbstractState::gas_constant | ( | void | ) |
Return the mole-fraction weighted gas constant in J/mol/K.
Definition at line 655 of file AbstractState.cpp.
|
inlinevirtual |
Get binary mixture double value (EXPERT USE ONLY!!!)
Reimplemented in CoolProp::VTPRBackend, CoolProp::HelmholtzEOSMixtureBackend, CoolProp::REFPROPMixtureBackend, and CoolProp::AbstractCubicBackend.
Definition at line 905 of file AbstractState.h.
|
inlinevirtual |
Get binary mixture double value (EXPERT USE ONLY!!!)
Reimplemented in CoolProp::AbstractCubicBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 901 of file AbstractState.h.
|
inlinevirtual |
Get binary mixture string value (EXPERT USE ONLY!!!)
Reimplemented in CoolProp::REFPROPMixtureBackend.
Definition at line 909 of file AbstractState.h.
|
inlinevirtual |
Get a constant for one of the fluids forming this mixture
| i | Index (0-based) of the fluid |
| param | parameter you want to obtain (probably one that is a trivial parameter) |
Reimplemented in CoolProp::AbstractCubicBackend, CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::PCSAFTBackend.
Definition at line 879 of file AbstractState.h.
|
inlinevirtual |
Double fluid parameter (currently the volume translation parameter for cubic)
Reimplemented in CoolProp::AbstractCubicBackend, and CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 925 of file AbstractState.h.
|
inlinevirtual |
Get the mass fractions of the fluid.
Definition at line 842 of file AbstractState.h.
|
pure virtual |
Get the mole fractions of the fluid.
Implemented in CoolProp::HelmholtzEOSMixtureBackend, CoolProp::REFPROPMixtureBackend, CoolProp::TabularBackend, CoolProp::AbstractCubicBackend, CoolProp::IF97Backend, CoolProp::IncompressibleBackend, and CoolProp::PCSAFTBackend.
|
inline |
After having calculated the phase envelope, return the phase envelope data.
Definition at line 1420 of file AbstractState.h.
|
inlinevirtual |
Get the state that is used in the equation of state or mixture model to reduce the state. For pure fluids this is usually, but not always, the critical point. For mixture models, it is usually composition dependent
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 937 of file AbstractState.h.
|
inline |
Get the data from the spinodal curve constructed in the call to build_spinodal()
Definition at line 988 of file AbstractState.h.
|
inline |
Get a desired state point - backend dependent.
Definition at line 942 of file AbstractState.h.
|
inline |
Return the Gibbs energy in J/kg.
Definition at line 1182 of file AbstractState.h.
|
inline |
Return the excess Gibbs energy in J/kg.
Definition at line 1188 of file AbstractState.h.
| double CoolProp::AbstractState::gibbsmolar | ( | void | ) |
Return the Gibbs energy in J/mol.
Definition at line 591 of file AbstractState.cpp.
| double CoolProp::AbstractState::gibbsmolar_excess | ( | void | ) |
Return the excess Gibbs energy in J/mol.
Definition at line 599 of file AbstractState.cpp.
| double CoolProp::AbstractState::gibbsmolar_residual | ( | void | ) |
Return the residual Gibbs energy in J/mol.
Definition at line 595 of file AbstractState.cpp.
|
inlinevirtual |
Return true if the fluid has a melting line - default is false, but can be re-implemented by derived class.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 1429 of file AbstractState.h.
|
inline |
Return the Helmholtz energy in J/kg.
Definition at line 1194 of file AbstractState.h.
|
inline |
Return the excess Helmholtz energy in J/kg.
Definition at line 1200 of file AbstractState.h.
| double CoolProp::AbstractState::helmholtzmolar | ( | void | ) |
Return the Helmholtz energy in J/mol.
Definition at line 603 of file AbstractState.cpp.
| double CoolProp::AbstractState::helmholtzmolar_excess | ( | void | ) |
Return the excess Helmholtz energy in J/mol.
Definition at line 607 of file AbstractState.cpp.
|
inline |
Return the mass enthalpy in J/kg.
Definition at line 1114 of file AbstractState.h.
|
inline |
Return the excess mass enthalpy in J/kg.
Definition at line 1120 of file AbstractState.h.
| double CoolProp::AbstractState::hmass_idealgas | ( | void | ) |
Return the ideal gas specific enthalpy in J/kg.
Definition at line 544 of file AbstractState.cpp.
| double CoolProp::AbstractState::hmolar | ( | void | ) |
Return the molar enthalpy in J/mol.
Definition at line 533 of file AbstractState.cpp.
| double CoolProp::AbstractState::hmolar_excess | ( | void | ) |
Return the excess molar enthalpy in J/mol.
Definition at line 547 of file AbstractState.cpp.
| double CoolProp::AbstractState::hmolar_idealgas | ( | void | ) |
Return the ideal gas molar enthalpy in J/mol.
Definition at line 541 of file AbstractState.cpp.
| double CoolProp::AbstractState::hmolar_residual | ( | void | ) |
Return the residual molar enthalpy in J/mol.
Definition at line 537 of file AbstractState.cpp.
|
inline |
Calculate an ideal curve for a pure fluid.
| type | The type of ideal curve you would like to calculate - "Ideal", "Boyle", "Joule-Thomson", "Joule Inversion", etc. |
| T | The temperatures along the curve in K |
| p | The pressures along the curve in Pa |
Definition at line 1251 of file AbstractState.h.
| double CoolProp::AbstractState::isentropic_expansion_coefficient | ( | void | ) |
Return the isentropic expansion coefficient \( \kappa_s = -\frac{c_p}{c_v}\frac{v}{p}\left.\frac{\partial p}{\partial v}\right|_T = \frac{\rho}{p}\left.\frac{\partial p}{\partial \rho}\right|_s\).
Definition at line 684 of file AbstractState.cpp.
|
inlineprotected |
Definition at line 90 of file AbstractState.h.
| double CoolProp::AbstractState::isobaric_expansion_coefficient | ( | void | ) |
Return the isobaric expansion coefficient \( \beta = \frac{1}{v}\left.\frac{\partial v}{\partial T}\right|_p = -\frac{1}{\rho}\left.\frac{\partial \rho}{\partial T}\right|_p\) in 1/K.
Definition at line 681 of file AbstractState.cpp.
| double CoolProp::AbstractState::isothermal_compressibility | ( | void | ) |
Return the isothermal compressibility \( \kappa = -\frac{1}{v}\left.\frac{\partial v}{\partial p}\right|_T=\frac{1}{\rho}\left.\frac{\partial \rho}{\partial p}\right|_T\) in 1/Pa.
Definition at line 678 of file AbstractState.cpp.
|
inlineprotected |
Definition at line 86 of file AbstractState.h.
|
inlineprotected |
Definition at line 94 of file AbstractState.h.
| double CoolProp::AbstractState::keyed_output | ( | parameters | key | ) |
Retrieve a value by key.
Definition at line 352 of file AbstractState.cpp.
|
protectedvirtual |
Convert mass-based input pair to molar-based input pair; If molar-based, do nothing.
< Mass density in kg/m^3, Temperature in K
< Entropy in J/kg/K, Temperature in K
< Mass density in kg/m^3, Pressure in Pa
< Mass density in kg/m^3, molar quality
< Enthalpy in J/kg, Pressure in Pa
< Pressure in Pa, Entropy in J/kg/K
< Pressure in Pa, Internal energy in J/kg
< Enthalpy in J/kg, Entropy in J/kg/K
< Entropy in J/kg/K, Internal energy in J/kg
< Mass density in kg/m^3, Enthalpy in J/kg
< Mass density in kg/m^3, Entropy in J/kg/K
< Mass density in kg/m^3, Internal energy in J/kg
Definition at line 200 of file AbstractState.cpp.
| double CoolProp::AbstractState::melting_line | ( | int | param, |
| int | given, | ||
| double | value | ||
| ) |
Return a value from the melting line
| param | The key for the parameter to be returned |
| given | The key for the parameter that is given |
| value | The value for the parameter that is given |
Definition at line 638 of file AbstractState.cpp.
| double CoolProp::AbstractState::molar_mass | ( | void | ) |
Return the molar mass in kg/mol.
Definition at line 651 of file AbstractState.cpp.
|
inline |
Get the mole fractions of the equilibrium liquid phase.
Definition at line 820 of file AbstractState.h.
|
inline |
Get the mole fractions of the equilibrium liquid phase (but as a double for use in SWIG wrapper)
Definition at line 824 of file AbstractState.h.
|
inline |
Get the mole fractions of the equilibrium vapor phase.
Definition at line 830 of file AbstractState.h.
|
inline |
Get the mole fractions of the equilibrium vapor phase (but as a double for use in SWIG wrapper)
Definition at line 834 of file AbstractState.h.
|
inline |
Return the name - backend dependent.
Definition at line 1035 of file AbstractState.h.
| double CoolProp::AbstractState::neff | ( | void | ) |
Return the effective hardness of interaction.
Definition at line 565 of file AbstractState.cpp.
|
inline |
Return the pressure in Pa.
Definition at line 1078 of file AbstractState.h.
| double CoolProp::AbstractState::p_critical | ( | void | ) |
Return the critical pressure in Pa.
Definition at line 512 of file AbstractState.cpp.
| double CoolProp::AbstractState::p_triple | ( | void | ) |
Return the triple point pressure in Pa.
Definition at line 515 of file AbstractState.cpp.
|
inline |
Get the phase of the state.
Definition at line 956 of file AbstractState.h.
|
inline |
Return the phase identification parameter (PIP) of G. Venkatarathnam and L.R. Oellrich, "Identification of the phase of a fluid using partial derivatives of pressure, volume, and temperature without reference to saturation properties: Applications in phase equilibria calculations".
Definition at line 1235 of file AbstractState.h.
| double CoolProp::AbstractState::pmax | ( | void | ) |
Get the maximum pressure in Pa.
Definition at line 500 of file AbstractState.cpp.
|
inline |
Return the Prandtl number (dimensionless)
Definition at line 1462 of file AbstractState.h.
|
inline |
Return the vapor quality (mol/mol); Q = 0 for saturated liquid.
Definition at line 1082 of file AbstractState.h.
|
inline |
Return the mass density in kg/m^3.
Definition at line 1074 of file AbstractState.h.
| double CoolProp::AbstractState::rhomass_critical | ( | void | ) |
Return the critical mass density in kg/m^3.
Definition at line 521 of file AbstractState.cpp.
| double CoolProp::AbstractState::rhomass_reducing | ( | void | ) |
Return the mass density at the reducing point in kg/m^3.
Definition at line 530 of file AbstractState.cpp.
|
inline |
Return the molar density in mol/m^3.
Definition at line 1070 of file AbstractState.h.
| double CoolProp::AbstractState::rhomolar_critical | ( | void | ) |
Return the critical molar density in mol/m^3.
Definition at line 518 of file AbstractState.cpp.
| double CoolProp::AbstractState::rhomolar_reducing | ( | void | ) |
Return the molar density at the reducing point in mol/m^3.
Definition at line 524 of file AbstractState.cpp.
|
inline |
Get an output from the saturated liquid state by key.
Definition at line 1057 of file AbstractState.h.
|
inline |
Get an output from the saturated vapor state by key.
Definition at line 1061 of file AbstractState.h.
| double CoolProp::AbstractState::saturation_ancillary | ( | parameters | param, |
| int | Q, | ||
| parameters | given, | ||
| double | value | ||
| ) |
Return the value from a saturation ancillary curve (if the backend implements it)
| param | The key for the parameter to be returned |
| Q | The quality for the parameter that is given (0 = saturated liquid, 1 = saturated vapor) |
| given | The key for the parameter that is given |
| value | The value for the parameter that is given |
Definition at line 644 of file AbstractState.cpp.
|
inline |
The second partial derivative in homogeneous phases.
The first partial derivative (CoolProp::AbstractState::first_partial_deriv) can be expressed as
\[ \left(\frac{\partial A}{\partial B}\right)_C = \frac{\left(\frac{\partial A}{\partial T}\right)_\rho\left(\frac{\partial C}{\partial \rho}\right)_T-\left(\frac{\partial A}{\partial \rho}\right)_T\left(\frac{\partial C}{\partial T}\right)_\rho}{\left(\frac{\partial B}{\partial T}\right)_\rho\left(\frac{\partial C}{\partial \rho}\right)_T-\left(\frac{\partial B}{\partial \rho}\right)_T\left(\frac{\partial C}{\partial T}\right)_\rho} = \frac{N}{D}\]
and the second derivative can be expressed as
\[ \frac{\partial}{\partial D}\left(\left(\frac{\partial A}{\partial B}\right)_C\right)_E = \frac{\frac{\partial}{\partial T}\left( \left(\frac{\partial A}{\partial B}\right)_C \right)_\rho\left(\frac{\partial E}{\partial \rho}\right)_T-\frac{\partial}{\partial \rho}\left(\left(\frac{\partial A}{\partial B}\right)_C\right)_T\left(\frac{\partial E}{\partial T}\right)_\rho}{\left(\frac{\partial D}{\partial T}\right)_\rho\left(\frac{\partial E}{\partial \rho}\right)_T-\left(\frac{\partial D}{\partial \rho}\right)_T\left(\frac{\partial E}{\partial T}\right)_\rho} \]
which can be expressed in parts as
\[\left(\frac{\partial N}{\partial \rho}\right)_{T} = \left(\frac{\partial A}{\partial T}\right)_\rho\left(\frac{\partial^2 C}{\partial \rho^2}\right)_{T}+\left(\frac{\partial^2 A}{\partial T\partial\rho}\right)\left(\frac{\partial C}{\partial \rho}\right)_{T}-\left(\frac{\partial A}{\partial \rho}\right)_T\left(\frac{\partial^2 C}{\partial T\partial\rho}\right)-\left(\frac{\partial^2 A}{\partial \rho^2}\right)_{T}\left(\frac{\partial C}{\partial T}\right)_\rho\]
\[\left(\frac{\partial D}{\partial \rho}\right)_{T} = \left(\frac{\partial B}{\partial T}\right)_\rho\left(\frac{\partial^2 C}{\partial \rho^2}\right)_{T}+\left(\frac{\partial^2 B}{\partial T\partial\rho}\right)\left(\frac{\partial C}{\partial \rho}\right)_{T}-\left(\frac{\partial B}{\partial \rho}\right)_T\left(\frac{\partial^2 C}{\partial T\partial\rho}\right)-\left(\frac{\partial^2 B}{\partial \rho^2}\right)_{T}\left(\frac{\partial C}{\partial T}\right)_\rho\]
\[\left(\frac{\partial N}{\partial T}\right)_{\rho} = \left(\frac{\partial A}{\partial T}\right)_\rho\left(\frac{\partial^2 C}{\partial \rho\partial T}\right)+\left(\frac{\partial^2 A}{\partial T^2}\right)_\rho\left(\frac{\partial C}{\partial \rho}\right)_{T}-\left(\frac{\partial A}{\partial \rho}\right)_T\left(\frac{\partial^2 C}{\partial T^2}\right)_\rho-\left(\frac{\partial^2 A}{\partial \rho\partial T}\right)\left(\frac{\partial C}{\partial T}\right)_\rho\]
\[\left(\frac{\partial D}{\partial T}\right)_{\rho} = \left(\frac{\partial B}{\partial T}\right)_\rho\left(\frac{\partial^2 C}{\partial \rho\partial T}\right)+\left(\frac{\partial^2 B}{\partial T^2}\right)_\rho\left(\frac{\partial C}{\partial \rho}\right)_{T}-\left(\frac{\partial B}{\partial \rho}\right)_T\left(\frac{\partial^2 C}{\partial T^2}\right)_\rho-\left(\frac{\partial^2 B}{\partial \rho\partial T}\right)\left(\frac{\partial C}{\partial T}\right)_\rho\]
\[\frac{\partial}{\partial \rho}\left( \left(\frac{\partial A}{\partial B}\right)_C \right)_T = \frac{D\left(\frac{\partial N}{\partial \rho}\right)_{T}-N\left(\frac{\partial D}{\partial \rho}\right)_{\tau}}{D^2}\]
\[\frac{\partial}{\partial T}\left( \left(\frac{\partial A}{\partial B}\right)_C \right)_\rho = \frac{D\left(\frac{\partial N}{\partial T}\right)_{\rho}-N\left(\frac{\partial D}{\partial T}\right)_{\rho}}{D^2}\]
The terms \( N \) and \( D \) are the numerator and denominator from CoolProp::AbstractState::first_partial_deriv respectively
Definition at line 1290 of file AbstractState.h.
|
inline |
The second partial derivative along the saturation curve.
Implementing the algorithms and ideas of: Matthis Thorade, Ali Saadat, "Partial derivatives of thermodynamic state properties for dynamic simulation", Environmental Earth Sciences, December 2013, Volume 70, Issue 8, pp 3497-3503
Like with first_saturation_deriv, we can express the derivative as
\[ \left(\frac{\partial y}{\partial T}\right)_{\sigma} = \left(\frac{\partial y}{\partial T}\right)+\left(\frac{\partial y}{\partial p}\right)\left(\frac{\partial p}{\partial T}\right)_{\sigma} \]
where \( y \) is already a saturation derivative. So you might end up with something like
\[ \left(\frac{\partial \left(\frac{\partial T}{\partial p}\right)_{\sigma}}{\partial T}\right)_{\sigma} = \left(\frac{\partial \left(\frac{\partial T}{\partial p}\right)_{\sigma}}{\partial T}\right)+\left(\frac{\partial \left(\frac{\partial T}{\partial p}\right)_{\sigma}}{\partial p}\right)\left(\frac{\partial p}{\partial T}\right)_{\sigma} \]
| Of1 | The parameter that the first derivative is taken of |
| Wrt1 | The parameter that the first derivative is taken with respect to |
| Wrt2 | The parameter that the second derivative is taken with respect to |
Definition at line 1336 of file AbstractState.h.
|
inline |
Calculate the second "two-phase" derivative as described by Thorade and Sadaat, EAS, 2013.
Implementing the algorithms and ideas of: Matthis Thorade, Ali Saadat, "Partial derivatives of thermodynamic state properties for dynamic simulation", Environmental Earth Sciences, December 2013, Volume 70, Issue 8, pp 3497-3503
| Of | The parameter to be derived |
| Wrt1 | The parameter that the derivative is taken with respect to in the first derivative |
| Constant1 | The parameter that is held constant in the first derivative |
| Wrt2 | The parameter that the derivative is taken with respect to in the second derivative |
| Constant2 | The parameter that is held constant in the second derivative |
Definition at line 1379 of file AbstractState.h.
|
inlinevirtual |
Set binary mixture floating point parameter (EXPERT USE ONLY!!!)
Reimplemented in CoolProp::VTPRBackend, CoolProp::HelmholtzEOSMixtureBackend, CoolProp::REFPROPMixtureBackend, and CoolProp::AbstractCubicBackend.
Definition at line 888 of file AbstractState.h.
|
inlinevirtual |
Set binary mixture floating point parameter (EXPERT USE ONLY!!!)
Reimplemented in CoolProp::AbstractCubicBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 884 of file AbstractState.h.
|
inlinevirtual |
Set binary mixture string parameter (EXPERT USE ONLY!!!)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 897 of file AbstractState.h.
|
inlinevirtual |
Set binary mixture string parameter (EXPERT USE ONLY!!!)
Definition at line 892 of file AbstractState.h.
|
inlinevirtual |
Set the cubic alpha function's constants:
Reimplemented in CoolProp::AbstractCubicBackend, and CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 917 of file AbstractState.h.
|
inlinevirtual |
Set fluid parameter (currently the volume translation parameter for cubic)
Reimplemented in CoolProp::AbstractCubicBackend, and CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 921 of file AbstractState.h.
|
pure virtual |
|
pure virtual |
|
inlinevirtual |
Set the reference state based on a thermodynamic state point specified by temperature and molar density
| T | Temperature at reference state [K] |
| rhomolar | Molar density at reference state [mol/m^3] |
| hmolar0 | Molar enthalpy at reference state [J/mol] |
| smolar0 | Molar entropy at reference state [J/mol/K] |
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 796 of file AbstractState.h.
|
inlinevirtual |
Set the reference state based on a string representation.
| reference_state | The reference state to use, one of |
| Reference State | Description |
|---|---|
| "IIR" | h = 200 kJ/kg, s=1 kJ/kg/K at 0C saturated liquid |
| "ASHRAE" | h = 0, s = 0 @ -40C saturated liquid |
| "NBP" | h = 0, s = 0 @ 1.0 bar saturated liquid |
| "DEF" | Reset to the default reference state for the fluid |
| "RESET" | Remove the offset |
The offset in the ideal gas Helmholtz energy can be obtained from
\[ \displaystyle\frac{\Delta s}{R_u/M}+\frac{\Delta h}{(R_u/M)T}\tau \]
where \( \Delta s = s-s_{spec} \) and \( \Delta h = h-h_{spec} \)
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 786 of file AbstractState.h.
|
inline |
Set the internal variable T without a flash call (expert use only!)
Definition at line 746 of file AbstractState.h.
|
inlinevirtual |
Reimplemented in CoolProp::AbstractCubicBackend, CoolProp::IF97Backend, CoolProp::IncompressibleBackend, and CoolProp::PCSAFTBackend.
Definition at line 763 of file AbstractState.h.
|
inline |
Return the molar entropy in J/kg/K.
Definition at line 1134 of file AbstractState.h.
|
inline |
Return the molar entropy in J/kg/K.
Definition at line 1140 of file AbstractState.h.
| double CoolProp::AbstractState::smass_idealgas | ( | void | ) |
Return the ideal gas specific entropy in J/kg/K.
Definition at line 562 of file AbstractState.cpp.
| double CoolProp::AbstractState::smolar | ( | void | ) |
Return the molar entropy in J/mol/K.
Definition at line 551 of file AbstractState.cpp.
| double CoolProp::AbstractState::smolar_excess | ( | void | ) |
Return the molar entropy in J/mol/K.
Definition at line 573 of file AbstractState.cpp.
| double CoolProp::AbstractState::smolar_idealgas | ( | void | ) |
Return the ideal gas molar entropy in J/mol/K.
Definition at line 559 of file AbstractState.cpp.
| double CoolProp::AbstractState::smolar_residual | ( | void | ) |
Return the residual molar entropy (as a function of temperature and density) in J/mol/K.
Definition at line 555 of file AbstractState.cpp.
|
inline |
Specify the phase for all further calculations with this state class.
Definition at line 960 of file AbstractState.h.
| double CoolProp::AbstractState::speed_sound | ( | void | ) |
Return the speed of sound in m/s.
Definition at line 626 of file AbstractState.cpp.
| double CoolProp::AbstractState::surface_tension | ( | void | ) |
Return the surface tension in N/m.
Definition at line 647 of file AbstractState.cpp.
|
inline |
Return the temperature in K.
Definition at line 1066 of file AbstractState.h.
| double CoolProp::AbstractState::T_critical | ( | void | ) |
Return the critical temperature in K.
Definition at line 503 of file AbstractState.cpp.
| double CoolProp::AbstractState::T_reducing | ( | void | ) |
Return the reducing point temperature in K.
Definition at line 506 of file AbstractState.cpp.
|
inline |
Return the tangent plane distance for a given trial composition w
| T | Temperature (K) |
| p | Pressure (Pa) |
| w | The trial composition |
| rhomolar_guess | (mol/m^3) The molar density guess value (if <0 (default), not used; if >0, guess value will be used in flash evaluation) |
\[ tpd(w) = \sum_i w_i(\ln w_i + \ln \phi_i(w) - d_i) \]
with
\[ d_i = \ln z_i + \ln \phi_i(z) \]
Or you can express the \( tpd \) in terms of fugacity (See Table 7.3 from GERG 2004 monograph) since \( \ln \phi_i = \ln f_i - \ln p -\ln z_i\) thus
\[ d_i = \ln f_i(z) - \ln p\]
and
\[ tpd(w) = \sum_i w_i(\ln f_i(w) - \ln p - d_i) \]
and the \( \ln p \) cancel, leaving
\[ tpd(w) = \sum_i w_i(\ln f_i(w) - \ln f_i(z)) \]
Definition at line 1020 of file AbstractState.h.
| double CoolProp::AbstractState::tau | ( | void | ) |
Return the reciprocal of the reduced temperature ( \(\tau = T_c/T\))
Definition at line 483 of file AbstractState.cpp.
| double CoolProp::AbstractState::Tmax | ( | void | ) |
Get the maximum temperature in K.
Definition at line 494 of file AbstractState.cpp.
| double CoolProp::AbstractState::Tmin | ( | void | ) |
Get the minimum temperature in K.
Definition at line 491 of file AbstractState.cpp.
| double CoolProp::AbstractState::trivial_keyed_output | ( | parameters | key | ) |
A trivial keyed output like molar mass that does not depend on the state.
Definition at line 290 of file AbstractState.cpp.
|
inline |
Calculate the "true" critical point for pure fluids where dpdrho|T and d2p/drho2|T are equal to zero.
Definition at line 1240 of file AbstractState.h.
| double CoolProp::AbstractState::Ttriple | ( | void | ) |
Get the triple point temperature in K.
Definition at line 497 of file AbstractState.cpp.
|
inline |
Return the mass internal energy in J/kg.
Definition at line 1146 of file AbstractState.h.
|
inline |
Return the excess internal energy in J/kg.
Definition at line 1152 of file AbstractState.h.
| double CoolProp::AbstractState::umass_idealgas | ( | void | ) |
Return the ideal gas specific internal energy in J/kg.
Definition at line 588 of file AbstractState.cpp.
| double CoolProp::AbstractState::umolar | ( | void | ) |
Return the molar internal energy in J/mol.
Definition at line 577 of file AbstractState.cpp.
| double CoolProp::AbstractState::umolar_excess | ( | void | ) |
Return the excess internal energy in J/mol.
Definition at line 581 of file AbstractState.cpp.
| double CoolProp::AbstractState::umolar_idealgas | ( | void | ) |
Return the ideal gas molar internal energy in J/mol.
Definition at line 585 of file AbstractState.cpp.
|
inline |
Unspecify the phase and go back to calculating it based on the inputs.
Definition at line 964 of file AbstractState.h.
|
pure virtual |
Update the state using two state variables.
Implemented in CoolProp::AbstractCubicBackend, CoolProp::HelmholtzEOSMixtureBackend, CoolProp::IF97Backend, CoolProp::IncompressibleBackend, CoolProp::PCSAFTBackend, CoolProp::TabularBackend, and CoolProp::REFPROPMixtureBackend.
|
inlinevirtual |
Update the state for QT inputs for pure fluids when using the superancillary functions.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 850 of file AbstractState.h.
|
inlineprotectedvirtual |
Update the states after having changed the reference state for enthalpy and entropy.
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend.
Definition at line 563 of file AbstractState.h.
|
inlinevirtual |
Update the state using two state variables and providing guess values Some or all of the guesses will be used - this is backend dependent
Reimplemented in CoolProp::HelmholtzEOSMixtureBackend, and CoolProp::REFPROPMixtureBackend.
Definition at line 856 of file AbstractState.h.
|
pure virtual |
|
pure virtual |
|
pure virtual |
| double CoolProp::AbstractState::viscosity | ( | void | ) |
Return the viscosity in Pa-s.
Definition at line 630 of file AbstractState.cpp.
|
inline |
Return the viscosity contributions, each in Pa-s.
Definition at line 1450 of file AbstractState.h.
|
inline |
Return the excess volume in m^3/kg.
Definition at line 1206 of file AbstractState.h.
| double CoolProp::AbstractState::volumemolar_excess | ( | void | ) |
Return the excess volume in m^3/mol.
Definition at line 611 of file AbstractState.cpp.
|
protected |
Cached low-level elements for in-place calculation of other properties.
Definition at line 138 of file AbstractState.h.
|
protected |
Definition at line 140 of file AbstractState.h.
|
protected |
Definition at line 116 of file AbstractState.h.
|
protected |
Definition at line 119 of file AbstractState.h.
|
protected |
Definition at line 119 of file AbstractState.h.
|
protected |
Two important points.
Definition at line 102 of file AbstractState.h.
|
protected |
Definition at line 119 of file AbstractState.h.
|
protected |
Definition at line 139 of file AbstractState.h.
|
protected |
Definition at line 139 of file AbstractState.h.
|
protected |
Definition at line 138 of file AbstractState.h.
|
protected |
Definition at line 141 of file AbstractState.h.
|
protected |
Definition at line 146 of file AbstractState.h.
|
protected |
Definition at line 141 of file AbstractState.h.
|
protected |
Definition at line 146 of file AbstractState.h.
|
protected |
Definition at line 141 of file AbstractState.h.
|
protected |
Definition at line 140 of file AbstractState.h.
|
protected |
Definition at line 140 of file AbstractState.h.
|
protected |
Definition at line 139 of file AbstractState.h.
|
protected |
Definition at line 139 of file AbstractState.h.
|
protected |
Definition at line 142 of file AbstractState.h.
|
protected |
Definition at line 147 of file AbstractState.h.
|
protected |
Definition at line 143 of file AbstractState.h.
|
protected |
Definition at line 142 of file AbstractState.h.
|
protected |
Definition at line 142 of file AbstractState.h.
|
protected |
Definition at line 144 of file AbstractState.h.
|
protected |
Definition at line 144 of file AbstractState.h.
|
protected |
Definition at line 144 of file AbstractState.h.
|
protected |
Definition at line 143 of file AbstractState.h.
|
protected |
Definition at line 143 of file AbstractState.h.
|
protected |
Definition at line 138 of file AbstractState.h.
|
protected |
Definition at line 138 of file AbstractState.h.
|
protected |
Definition at line 141 of file AbstractState.h.
|
protected |
Definition at line 146 of file AbstractState.h.
|
protected |
Definition at line 140 of file AbstractState.h.
|
protected |
Definition at line 113 of file AbstractState.h.
|
protected |
Definition at line 135 of file AbstractState.h.
|
protected |
Definition at line 135 of file AbstractState.h.
|
protected |
Some administrative variables.
Definition at line 82 of file AbstractState.h.
|
protected |
Definition at line 132 of file AbstractState.h.
|
protected |
Universal gas constant [J/mol/K].
Definition at line 108 of file AbstractState.h.
|
protected |
Definition at line 119 of file AbstractState.h.
|
protected |
Definition at line 126 of file AbstractState.h.
|
protected |
Definition at line 123 of file AbstractState.h.
|
protected |
Definition at line 120 of file AbstractState.h.
|
protected |
Definition at line 127 of file AbstractState.h.
|
protected |
Definition at line 118 of file AbstractState.h.
|
protected |
Excess properties.
Definition at line 126 of file AbstractState.h.
|
protected |
Residual properties.
Definition at line 123 of file AbstractState.h.
|
protected |
Definition at line 118 of file AbstractState.h.
|
protected |
Definition at line 118 of file AbstractState.h.
|
protected |
Molar mass [mol/kg].
Definition at line 105 of file AbstractState.h.
|
protected |
Definition at line 111 of file AbstractState.h.
|
protected |
The key for the phase from CoolProp::phases enum.
Definition at line 83 of file AbstractState.h.
|
protected |
Definition at line 130 of file AbstractState.h.
|
protected |
Definition at line 130 of file AbstractState.h.
|
protected |
Definition at line 111 of file AbstractState.h.
|
protected |
Definition at line 102 of file AbstractState.h.
|
protected |
Smoothing values.
Definition at line 135 of file AbstractState.h.
|
protected |
Ancillary values.
Definition at line 130 of file AbstractState.h.
|
protected |
Two-Phase variables.
Definition at line 150 of file AbstractState.h.
|
protected |
Bulk values.
Definition at line 111 of file AbstractState.h.
|
protected |
Definition at line 130 of file AbstractState.h.
|
protected |
Definition at line 150 of file AbstractState.h.
|
protected |
Definition at line 118 of file AbstractState.h.
|
protected |
Definition at line 126 of file AbstractState.h.
|
protected |
Definition at line 123 of file AbstractState.h.
|
protected |
Definition at line 119 of file AbstractState.h.
|
protected |
Definition at line 116 of file AbstractState.h.
|
protected |
Definition at line 111 of file AbstractState.h.
|
protected |
Definition at line 113 of file AbstractState.h.
|
protected |
Definition at line 130 of file AbstractState.h.
|
protected |
Definition at line 130 of file AbstractState.h.
|
protected |
Definition at line 118 of file AbstractState.h.
|
protected |
Definition at line 126 of file AbstractState.h.
|
protected |
Transport properties.
Definition at line 116 of file AbstractState.h.
|
protected |
Definition at line 127 of file AbstractState.h.
|
protected |
Definition at line 98 of file AbstractState.h.
|
protected |
If the phase is imposed, the imposed phase index.
Definition at line 84 of file AbstractState.h.