|
CoolProp 7.2.0
An open-source fluid property and humid air property database
|
Definition at line 62 of file HelmholtzEOSMixtureBackend.h.
#include <HelmholtzEOSMixtureBackend.h>
Public Types | |
| enum | StationaryPointReturnFlag { ZERO_STATIONARY_POINTS , ONE_STATIONARY_POINT_FOUND , TWO_STATIONARY_POINTS_FOUND } |
Public Member Functions | |
| HelmholtzEOSMixtureBackend () | |
| HelmholtzEOSMixtureBackend (const std::vector< CoolPropFluid > &components, bool generate_SatL_and_SatV=true) | |
| HelmholtzEOSMixtureBackend (const std::vector< std::string > &component_names, bool generate_SatL_and_SatV=true) | |
| virtual HelmholtzEOSMixtureBackend * | get_copy (bool generate_SatL_and_SatV=true) |
| void | sync_linked_states (const HelmholtzEOSMixtureBackend *const) |
| virtual | ~HelmholtzEOSMixtureBackend () |
| std::string | backend_name (void) |
| bool | clear () |
| Clear all the cached values. More... | |
| bool | using_mole_fractions () |
| bool | using_mass_fractions () |
| bool | using_volu_fractions () |
| bool | is_pure () |
| bool | has_melting_line () |
| Return true if the fluid has a melting line - default is false, but can be re-implemented by derived class. More... | |
| CoolPropDbl | calc_melting_line (int param, int given, CoolPropDbl value) |
| std::string | fluid_param_string (const std::string &) |
| Return a string from the backend for the mixture/fluid. More... | |
| virtual void | set_reference_stateS (const std::string &reference_state) |
| brief Set the reference state based on a string representation More... | |
| virtual void | set_reference_stateD (double T, double rhomolar, double hmolar0, double smolar0) |
| Set the reference state based on a thermodynamic state point specified by temperature and molar density. 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. 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. More... | |
| void | set_binary_interaction_string (const std::size_t i, const std::size_t j, const std::string ¶meter, const std::string &value) |
| Set a binary interaction string. More... | |
| void | apply_simple_mixing_rule (std::size_t i, std::size_t j, const std::string &model) |
| Apply a simple mixing rule. 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... | |
| phases | calc_phase (void) |
| Using this backend, calculate the phase. More... | |
| void | calc_specify_phase (phases phase_index) |
| Specify the phase - this phase will always be used in calculations. More... | |
| void | calc_unspecify_phase () |
| Unspecify the phase - the phase is no longer imposed, different solvers can do as they like. More... | |
| CoolPropDbl | calc_saturation_ancillary (parameters param, int Q, parameters given, double value) |
| void | calc_ssat_max (void) |
| void | calc_hsat_max (void) |
| CoolPropDbl | calc_GWP20 () |
| Using this backend, calculate the 20-year global warming potential (GWP) More... | |
| CoolPropDbl | calc_GWP500 () |
| Using this backend, calculate the 500-year global warming potential (GWP) More... | |
| CoolPropDbl | calc_GWP100 () |
| Using this backend, calculate the 100-year global warming potential (GWP) More... | |
| CoolPropDbl | calc_ODP () |
| Using this backend, calculate the ozone depletion potential (ODP) More... | |
| CoolPropDbl | calc_first_saturation_deriv (parameters Of1, parameters Wrt1) |
| CoolPropDbl | calc_first_saturation_deriv (parameters Of1, parameters Wrt1, HelmholtzEOSMixtureBackend &SatL, HelmholtzEOSMixtureBackend &SatV) |
| CoolPropDbl | calc_second_saturation_deriv (parameters Of1, parameters Wrt1, parameters Wrt2) |
| CoolPropDbl | calc_first_two_phase_deriv (parameters Of, parameters Wrt, parameters Constant) |
| CoolPropDbl | calc_second_two_phase_deriv (parameters Of, parameters Wrt1, parameters Constant1, parameters Wrt2, parameters Constant2) |
| CoolPropDbl | calc_first_two_phase_deriv_splined (parameters Of, parameters Wrt, parameters Constant, CoolPropDbl x_end) |
| CriticalState | calc_critical_point (double rho0, double T0) |
| std::vector< CoolProp::CriticalState > | calc_all_critical_points () |
| An overload to make the compiler (clang in this case) happy. More... | |
| virtual void | get_critical_point_starting_values (double &delta0, double &tau0) |
| virtual void | get_critical_point_search_radii (double &R_delta, double &R_tau) |
| Get the search radius in delta and tau for the tracer. More... | |
| virtual bool | get_critical_is_terminated (double &delta, double &tau) |
| Checking function to see if we should stop the tracing of the critical contour. More... | |
| virtual void | calc_build_spinodal () |
| Build the spinodal curve. More... | |
| virtual SpinodalData | calc_get_spinodal_data () |
| Get the data from the spinodal curve. More... | |
| void | calc_criticality_contour_values (double &L1star, double &M1star) |
| Calculate the values \(\mathcal{L}_1^*\) and \(\mathcal{M}_1^*\). More... | |
| double | calc_tangent_plane_distance (const double T, const double p, const std::vector< double > &w, const double rhomolar_guess) |
| Calculate tangent plane distance for given trial composition w. More... | |
| void | recalculate_singlephase_phase () |
| Calculate the phase once the state is fully calculated but you aren't sure if it is liquid or gas or ... More... | |
| void | calc_change_EOS (const std::size_t i, const std::string &EOS_name) |
| Change the equation of state for one component. More... | |
| const CoolProp::SimpleState & | calc_state (const std::string &state) |
| virtual const double | get_fluid_constant (std::size_t i, parameters param) const |
| const std::vector< CoolPropFluid > & | get_components () const |
| std::vector< CoolPropFluid > & | get_components () |
| std::vector< CoolPropDbl > & | get_K () |
| std::vector< CoolPropDbl > & | get_lnK () |
| HelmholtzEOSMixtureBackend & | get_SatL () |
| HelmholtzEOSMixtureBackend & | get_SatV () |
| std::vector< CoolPropDbl > | calc_mole_fractions_liquid (void) |
| std::vector< CoolPropDbl > | calc_mole_fractions_vapor (void) |
| const std::vector< CoolPropDbl > | calc_mass_fractions (void) |
| const CoolProp::PhaseEnvelopeData & | calc_phase_envelope_data () |
| void | calc_conformal_state (const std::string &reference_fluid, CoolPropDbl &T, CoolPropDbl &rhomolar) |
| Calculate the conformal state (unity shape factors starting point if T < 0 and rhomolar < 0) More... | |
| void | resize (std::size_t N) |
| virtual void | update (CoolProp::input_pairs input_pair, double value1, double value2) |
| The standard update function. More... | |
| void | update_with_guesses (CoolProp::input_pairs input_pair, double Value1, double Value2, const GuessesStructure &guesses) |
| Update the state using guess values. More... | |
| void | update_internal (HelmholtzEOSMixtureBackend &HEOS) |
| Update all the internal variables for a state by copying from another state. More... | |
| void | update_TP_guessrho (CoolPropDbl T, CoolPropDbl p, CoolPropDbl rho_guess) |
| Update with TP and a guess for rho. More... | |
| void | update_DmolarT_direct (CoolPropDbl rhomolar, CoolPropDbl T) |
| void | update_TDmolarP_unchecked (CoolPropDbl T, CoolPropDbl rhomolarL, CoolPropDbl p) |
| void | update_QT_pure_superanc (CoolPropDbl Q, CoolPropDbl T) |
| Update the state for QT inputs for pure fluids when using the superancillary functions. More... | |
| void | update_HmolarQ_with_guessT (CoolPropDbl hmolar, CoolPropDbl Q, CoolPropDbl Tguess) |
| virtual void | set_components (const std::vector< CoolPropFluid > &components, bool generate_SatL_and_SatV=true) |
| Set the components of the mixture. More... | |
| void | set_mixture_parameters () |
| Set the mixture parameters - binary pair reducing functions, departure functions, F_ij, etc. More... | |
| void | set_mole_fractions (const std::vector< CoolPropDbl > &mf) |
| Set the mole fractions. More... | |
| const std::vector< CoolPropDbl > & | get_mole_fractions () |
| Get the mole fractions of the fluid. More... | |
| std::vector< CoolPropDbl > & | get_mole_fractions_ref () |
| std::vector< double > & | get_mole_fractions_doubleref (void) |
| void | set_mass_fractions (const std::vector< CoolPropDbl > &mass_fractions) |
| Set the mass fractions. More... | |
| void | calc_ideal_curve (const std::string &type, std::vector< double > &T, std::vector< double > &p) |
| CoolPropDbl | calc_molar_mass (void) |
| Using this backend, calculate the molar mass in kg/mol. More... | |
| CoolPropDbl | calc_gas_constant (void) |
| Using this backend, calculate the universal gas constant \(R_u\) in J/mol/K. More... | |
| CoolPropDbl | calc_acentric_factor (void) |
| Using this backend, calculate the acentric factor. More... | |
| CoolPropDbl | calc_Bvirial (void) |
| Using this backend, calculate the second virial coefficient. More... | |
| CoolPropDbl | calc_Cvirial (void) |
| Using this backend, calculate the third virial coefficient. More... | |
| CoolPropDbl | calc_dBvirial_dT (void) |
| Using this backend, calculate the derivative dB/dT. More... | |
| CoolPropDbl | calc_dCvirial_dT (void) |
| Using this backend, calculate the derivative dC/dT. More... | |
| CoolPropDbl | calc_pressure (void) |
| Using this backend, calculate the pressure in Pa. More... | |
| CoolPropDbl | calc_cvmolar (void) |
| Using this backend, calculate the molar constant-volume specific heat in J/mol/K. More... | |
| CoolPropDbl | calc_cpmolar (void) |
| Using this backend, calculate the molar constant-pressure specific heat in J/mol/K. More... | |
| CoolPropDbl | calc_gibbsmolar (void) |
| Using this backend, calculate the molar Gibbs function in J/mol. More... | |
| CoolPropDbl | calc_gibbsmolar_residual (void) |
| Using this backend, calculate the residual molar Gibbs function in J/mol. More... | |
| CoolPropDbl | calc_gibbsmolar_nocache (CoolPropDbl T, CoolPropDbl rhomolar) |
| CoolPropDbl | calc_helmholtzmolar (void) |
| Using this backend, calculate the molar Helmholtz energy in J/mol. More... | |
| CoolPropDbl | calc_cpmolar_idealgas (void) |
| Using this backend, calculate the ideal gas molar constant-pressure specific heat in J/mol/K. More... | |
| CoolPropDbl | calc_pressure_nocache (CoolPropDbl T, CoolPropDbl rhomolar) |
| CoolPropDbl | calc_smolar (void) |
| Using this backend, calculate the molar entropy in J/mol/K. More... | |
| CoolPropDbl | calc_smolar_residual (void) |
| Using this backend, calculate the residual molar entropy in J/mol/K. More... | |
| CoolPropDbl | calc_smolar_nocache (CoolPropDbl T, CoolPropDbl rhomolar) |
| CoolPropDbl | calc_hmolar (void) |
| Using this backend, calculate the molar enthalpy in J/mol. More... | |
| CoolPropDbl | calc_hmolar_residual (void) |
| Using this backend, calculate the residual molar enthalpy in J/mol. More... | |
| CoolPropDbl | calc_hmolar_nocache (CoolPropDbl T, CoolPropDbl rhomolar) |
| CoolPropDbl | calc_umolar_nocache (CoolPropDbl T, CoolPropDbl rhomolar) |
| CoolPropDbl | calc_umolar (void) |
| Using this backend, calculate the molar internal energy in J/mol. More... | |
| CoolPropDbl | calc_speed_sound (void) |
| Using this backend, calculate the speed of sound in m/s. More... | |
| void | calc_excess_properties () |
| Using this backend, calculate and cache the excess properties. More... | |
| CoolPropDbl | calc_phase_identification_parameter (void) |
| The phase identification parameter of Venkatarathnam et al., FPE, 2011. More... | |
| CoolPropDbl | calc_fugacity (std::size_t i) |
| Using this backend, calculate the fugacity in Pa. More... | |
| virtual CoolPropDbl | calc_fugacity_coefficient (std::size_t i) |
| Using this backend, calculate the fugacity coefficient (dimensionless) More... | |
| CoolPropDbl | calc_chemical_potential (std::size_t i) |
| Using this backend, calculate the chemical potential in J/mol. More... | |
| CoolPropDbl | calc_flame_hazard (void) |
| Using this backend, calculate the flame hazard. More... | |
| CoolPropDbl | calc_health_hazard (void) |
| Using this backend, calculate the health hazard. More... | |
| CoolPropDbl | calc_physical_hazard (void) |
| Using this backend, calculate the physical hazard. More... | |
| CoolPropDbl | calc_alphar (void) |
| Using this backend, calculate the residual Helmholtz energy term \(\alpha^r\) (dimensionless) More... | |
| CoolPropDbl | calc_dalphar_dDelta (void) |
| Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta}\) (dimensionless) More... | |
| CoolPropDbl | calc_dalphar_dTau (void) |
| Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\tau}\) (dimensionless) More... | |
| CoolPropDbl | calc_d2alphar_dDelta2 (void) |
| Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\delta}\) (dimensionless) More... | |
| CoolPropDbl | calc_d2alphar_dDelta_dTau (void) |
| Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\tau}\) (dimensionless) More... | |
| CoolPropDbl | calc_d2alphar_dTau2 (void) |
| Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\tau\tau}\) (dimensionless) More... | |
| CoolPropDbl | calc_d3alphar_dDelta3 (void) |
| Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\delta\delta}\) (dimensionless) More... | |
| CoolPropDbl | calc_d3alphar_dDelta2_dTau (void) |
| Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\delta\tau}\) (dimensionless) More... | |
| CoolPropDbl | calc_d3alphar_dDelta_dTau2 (void) |
| Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\tau\tau}\) (dimensionless) More... | |
| CoolPropDbl | calc_d3alphar_dTau3 (void) |
| Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\tau\tau\tau}\) (dimensionless) More... | |
| CoolPropDbl | calc_d4alphar_dDelta4 (void) |
| Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\delta\delta\delta}\) (dimensionless) More... | |
| CoolPropDbl | calc_d4alphar_dDelta3_dTau (void) |
| Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\delta\delta\tau}\) (dimensionless) More... | |
| CoolPropDbl | calc_d4alphar_dDelta2_dTau2 (void) |
| Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\delta\tau\tau}\) (dimensionless) More... | |
| CoolPropDbl | calc_d4alphar_dDelta_dTau3 (void) |
| Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\tau\tau\tau}\) (dimensionless) More... | |
| CoolPropDbl | calc_d4alphar_dTau4 (void) |
| Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\tau\tau\tau\tau}\) (dimensionless) More... | |
| CoolPropDbl | calc_alpha0 (void) |
| Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0\) (dimensionless) More... | |
| CoolPropDbl | calc_dalpha0_dDelta (void) |
| Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\delta}\) (dimensionless) More... | |
| CoolPropDbl | calc_dalpha0_dTau (void) |
| Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\tau}\) (dimensionless) More... | |
| CoolPropDbl | calc_d2alpha0_dDelta2 (void) |
| Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\delta\tau}\) (dimensionless) More... | |
| CoolPropDbl | calc_d2alpha0_dDelta_dTau (void) |
| Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\delta\delta}\) (dimensionless) More... | |
| CoolPropDbl | calc_d2alpha0_dTau2 (void) |
| Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\tau\tau}\) (dimensionless) More... | |
| CoolPropDbl | calc_d3alpha0_dDelta3 (void) |
| Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\delta\delta\delta}\) (dimensionless) More... | |
| CoolPropDbl | calc_d3alpha0_dDelta2_dTau (void) |
| Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\delta\delta\tau}\) (dimensionless) More... | |
| CoolPropDbl | calc_d3alpha0_dDelta_dTau2 (void) |
| Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\delta\tau\tau}\) (dimensionless) More... | |
| CoolPropDbl | calc_d3alpha0_dTau3 (void) |
| Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\tau\tau\tau}\) (dimensionless) More... | |
| CoolPropDbl | calc_surface_tension (void) |
| Using this backend, calculate the surface tension in N/m. More... | |
| CoolPropDbl | calc_viscosity (void) |
| Using this backend, calculate the viscosity in Pa-s. More... | |
| CoolPropDbl | calc_viscosity_dilute (void) |
| CoolPropDbl | calc_viscosity_background (void) |
| CoolPropDbl | calc_viscosity_background (CoolPropDbl eta_dilute, CoolPropDbl &initial_density, CoolPropDbl &residual) |
| CoolPropDbl | calc_conductivity (void) |
| Using this backend, calculate the thermal conductivity in W/m/K. More... | |
| CoolPropDbl | calc_conductivity_background (void) |
| void | calc_viscosity_contributions (CoolPropDbl &dilute, CoolPropDbl &initial_density, CoolPropDbl &residual, CoolPropDbl &critical) |
| Calculate each of the contributions to the viscosity. More... | |
| void | calc_conductivity_contributions (CoolPropDbl &dilute, CoolPropDbl &initial_density, CoolPropDbl &residual, CoolPropDbl &critical) |
| Calculate each of the contributions to the conductivity. More... | |
| CoolPropDbl | calc_saturated_liquid_keyed_output (parameters key) |
| CoolPropDbl | calc_saturated_vapor_keyed_output (parameters key) |
| CoolPropDbl | calc_Tmin (void) |
| Using this backend, calculate the minimum temperature in K. More... | |
| CoolPropDbl | calc_Tmax (void) |
| Using this backend, calculate the maximum temperature in K. More... | |
| CoolPropDbl | calc_pmax (void) |
| Using this backend, calculate the maximum pressure in Pa. More... | |
| CoolPropDbl | calc_Ttriple (void) |
| Using this backend, get the triple point temperature in K. More... | |
| CoolPropDbl | calc_p_triple (void) |
| Using this backend, get the triple point pressure in Pa. More... | |
| CoolPropDbl | calc_pmax_sat (void) |
| CoolPropDbl | calc_Tmax_sat (void) |
| void | calc_Tmin_sat (CoolPropDbl &Tmin_satL, CoolPropDbl &Tmin_satV) |
| void | calc_pmin_sat (CoolPropDbl &pmin_satL, CoolPropDbl &pmin_satV) |
| virtual CoolPropDbl | calc_T_critical (void) |
| Using this backend, get the critical 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_rhomolar_critical (void) |
| Using this backend, get the critical point molar density in mol/m^3. More... | |
| CoolPropDbl | calc_T_reducing (void) |
| Using this backend, get the reducing point temperature in K. More... | |
| CoolPropDbl | calc_rhomolar_reducing (void) |
| Using this backend, get the reducing point molar density in mol/m^3. More... | |
| CoolPropDbl | calc_p_reducing (void) |
| Using this backend, get the reducing point pressure in Pa. More... | |
| CoolPropDbl | calc_PIP (void) |
| Using this backend, calculate the phase identification parameter (PIP) More... | |
| std::string | calc_name (void) |
| Using this backend, get the name of the fluid. More... | |
| std::vector< std::string > | calc_fluid_names (void) |
| Using this backend, get a vector of fluid names. More... | |
| void | calc_all_alphar_deriv_cache (const std::vector< CoolPropDbl > &mole_fractions, const CoolPropDbl &tau, const CoolPropDbl &delta) |
| virtual CoolPropDbl | calc_alphar_deriv_nocache (const int nTau, const int nDelta, const std::vector< CoolPropDbl > &mole_fractions, const CoolPropDbl &tau, const CoolPropDbl &delta) |
| CoolPropDbl | calc_alpha0_deriv_nocache (const int nTau, const int nDelta, const std::vector< CoolPropDbl > &mole_fractions, const CoolPropDbl &tau, const CoolPropDbl &delta, const CoolPropDbl &Tr, const CoolPropDbl &rhor) |
| Take derivatives of the ideal-gas part of the Helmholtz energy, don't use any cached values, or store any cached values. More... | |
| HelmholtzDerivatives | calc_all_alpha0_derivs_nocache (const std::vector< CoolPropDbl > &mole_fractions, const CoolPropDbl &tau, const CoolPropDbl &delta, const CoolPropDbl &Tr, const CoolPropDbl &rhor) |
| virtual void | calc_reducing_state (void) |
| virtual SimpleState | calc_reducing_state_nocache (const std::vector< CoolPropDbl > &mole_fractions) |
| const CoolProp::SimpleState & | get_reducing_state () |
| void | update_states () |
| Update the states after having changed the reference state for enthalpy and entropy. More... | |
| CoolPropDbl | calc_compressibility_factor (void) |
| Using this backend, calculate the compressibility factor Z \( Z = p/(\rho R T) \). More... | |
| 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... | |
| CoolPropDbl | SRK_covolume () |
| void | T_phase_determination_pure_or_pseudopure (int other, CoolPropDbl value) |
| void | p_phase_determination_pure_or_pseudopure (int other, CoolPropDbl value, bool &saturation_called) |
| void | DmolarP_phase_determination () |
| virtual CoolPropDbl | solver_rho_Tp (CoolPropDbl T, CoolPropDbl p, CoolPropDbl rho_guess=-1) |
| virtual CoolPropDbl | solver_rho_Tp_SRK (CoolPropDbl T, CoolPropDbl p, phases phase) |
| virtual StationaryPointReturnFlag | solver_dpdrho0_Tp (CoolPropDbl T, CoolPropDbl p, CoolPropDbl rhomax, CoolPropDbl &light, CoolPropDbl &heavy) |
| virtual CoolPropDbl | solver_rho_Tp_global (CoolPropDbl T, CoolPropDbl p, CoolPropDbl rhomax) |
Public Member Functions inherited from CoolProp::AbstractState | |
| 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... | |
Public Attributes | |
| shared_ptr< ReducingFunction > | Reducing |
| shared_ptr< ResidualHelmholtz > | residual_helmholtz |
| PhaseEnvelopeData | PhaseEnvelope |
| SimpleState | hsat_max |
| SsatSimpleState | ssat_max |
| SpinodalData | spinodal_values |
| shared_ptr< HelmholtzEOSMixtureBackend > | SatL |
| shared_ptr< HelmholtzEOSMixtureBackend > | SatV |
Protected Member Functions | |
| void | pre_update (CoolProp::input_pairs &input_pair, CoolPropDbl &value1, CoolPropDbl &value2) |
| void | post_update (bool optional_checks=true) |
| virtual void | add_TPD_state () |
| Update the state class used to calculate the tangent-plane-distance. More... | |
| virtual void | add_critical_state () |
| Update the state class used to calculate the critical point(s) More... | |
| virtual void | add_transient_pure_state () |
| Update the state class used to calculate the critical point(s) More... | |
| std::vector< CoolProp::CriticalState > | _calc_all_critical_points (bool find_critical_points=true) |
| This overload is protected because it doesn't follow the base class definition, since this function is needed for constructing spinodals. More... | |
| std::optional< EquationOfState::SuperAncillary_t > & | get_superanc_optional () |
Protected Member Functions inherited from CoolProp::AbstractState | |
| 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... | |
Static Protected Member Functions | |
| static void | set_fluid_enthalpy_entropy_offset (CoolPropFluid &component, double delta_a1, double delta_a2, const std::string &ref) |
Protected Attributes | |
| std::vector< shared_ptr< HelmholtzEOSMixtureBackend > > | linked_states |
| States that are linked to this one, and should be updated (BIP, reference state, etc.) More... | |
| shared_ptr< HelmholtzEOSMixtureBackend > | transient_pure_state |
| A temporary state used for calculations of pure fluid properties. More... | |
| shared_ptr< HelmholtzEOSMixtureBackend > | TPD_state |
| A temporary state used for calculations of the tangent-plane-distance. More... | |
| shared_ptr< HelmholtzEOSMixtureBackend > | critical_state |
| std::vector< CoolPropFluid > | components |
| The components that are in use. More... | |
| bool | is_pure_or_pseudopure |
| A flag for whether the substance is a pure or pseudo-pure fluid (true) or a mixture (false) More... | |
| MoleFractions | mole_fractions |
| The bulk mole fractions of the mixture. More... | |
| std::vector< CoolPropDbl > | K |
| The K factors for the components. More... | |
| std::vector< CoolPropDbl > | lnK |
| The natural logarithms of the K factors of the components. More... | |
| SimpleState | _crit |
| std::size_t | N |
| Number of components. More... | |
Protected Attributes inherited from CoolProp::AbstractState | |
| 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() |
Friends | |
| class | FlashRoutines |
| class | TransportRoutines |
| class | MixtureDerivatives |
| class | PhaseEnvelopeRoutines |
| class | MixtureParameters |
| class | CorrespondingStatesTerm |
Additional Inherited Members | |
Static Public Member Functions inherited from CoolProp::AbstractState | |
| 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 inherited from CoolProp::AbstractState | |
| using | CAE = CacheArrayElement< double > |
| Enumerator | |
|---|---|
| ZERO_STATIONARY_POINTS | |
| ONE_STATIONARY_POINT_FOUND | |
| TWO_STATIONARY_POINTS_FOUND | |
Definition at line 637 of file HelmholtzEOSMixtureBackend.h.
| CoolProp::HelmholtzEOSMixtureBackend::HelmholtzEOSMixtureBackend | ( | ) |
Definition at line 58 of file HelmholtzEOSMixtureBackend.cpp.
| CoolProp::HelmholtzEOSMixtureBackend::HelmholtzEOSMixtureBackend | ( | const std::vector< CoolPropFluid > & | components, |
| bool | generate_SatL_and_SatV = true |
||
| ) |
Definition at line 81 of file HelmholtzEOSMixtureBackend.cpp.
| CoolProp::HelmholtzEOSMixtureBackend::HelmholtzEOSMixtureBackend | ( | const std::vector< std::string > & | component_names, |
| bool | generate_SatL_and_SatV = true |
||
| ) |
Definition at line 66 of file HelmholtzEOSMixtureBackend.cpp.
|
inlinevirtual |
Definition at line 123 of file HelmholtzEOSMixtureBackend.h.
|
protected |
This overload is protected because it doesn't follow the base class definition, since this function is needed for constructing spinodals.
Definition at line 4076 of file HelmholtzEOSMixtureBackend.cpp.
|
inlineprotectedvirtual |
Update the state class used to calculate the critical point(s)
Definition at line 82 of file HelmholtzEOSMixtureBackend.h.
|
inlineprotectedvirtual |
Update the state class used to calculate the tangent-plane-distance.
Definition at line 74 of file HelmholtzEOSMixtureBackend.h.
|
inlineprotectedvirtual |
Update the state class used to calculate the critical point(s)
Definition at line 90 of file HelmholtzEOSMixtureBackend.h.
|
virtual |
Apply a simple mixing rule.
Reimplemented from CoolProp::AbstractState.
Definition at line 301 of file HelmholtzEOSMixtureBackend.cpp.
|
inlinevirtual |
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
Implements CoolProp::AbstractState.
Definition at line 124 of file HelmholtzEOSMixtureBackend.h.
|
virtual |
Using this backend, calculate the acentric factor.
Reimplemented from CoolProp::AbstractState.
Definition at line 508 of file HelmholtzEOSMixtureBackend.cpp.
| HelmholtzDerivatives CoolProp::HelmholtzEOSMixtureBackend::calc_all_alpha0_derivs_nocache | ( | const std::vector< CoolPropDbl > & | mole_fractions, |
| const CoolPropDbl & | tau, | ||
| const CoolPropDbl & | delta, | ||
| const CoolPropDbl & | Tr, | ||
| const CoolPropDbl & | rhor | ||
| ) |
Definition at line 3369 of file HelmholtzEOSMixtureBackend.cpp.
| void CoolProp::HelmholtzEOSMixtureBackend::calc_all_alphar_deriv_cache | ( | const std::vector< CoolPropDbl > & | mole_fractions, |
| const CoolPropDbl & | tau, | ||
| const CoolPropDbl & | delta | ||
| ) |
Definition at line 3253 of file HelmholtzEOSMixtureBackend.cpp.
|
inlinevirtual |
An overload to make the compiler (clang in this case) happy.
Reimplemented from CoolProp::AbstractState.
Definition at line 239 of file HelmholtzEOSMixtureBackend.h.
|
virtual |
Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0\) (dimensionless)
Reimplemented from CoolProp::AbstractState.
Definition at line 3493 of file HelmholtzEOSMixtureBackend.cpp.
| CoolPropDbl CoolProp::HelmholtzEOSMixtureBackend::calc_alpha0_deriv_nocache | ( | const int | nTau, |
| const int | nDelta, | ||
| const std::vector< CoolPropDbl > & | mole_fractions, | ||
| const CoolPropDbl & | tau, | ||
| const CoolPropDbl & | delta, | ||
| const CoolPropDbl & | Tr, | ||
| const CoolPropDbl & | rhor | ||
| ) |
Take derivatives of the ideal-gas part of the Helmholtz energy, don't use any cached values, or store any cached values.
| nTau | How many derivatives with respect to \(\tau\) to take |
| nDelta | How many derivatives with respect to \(\delta\) to take |
| mole_fractions | Mole fractions |
| tau | Reciprocal reduced temperature where \(\tau=T_r / T\) |
| delta | Reduced density where \(\delta = \rho / \rho_r \) |
| Tr | Reducing temperature of the mixture [K] |
| rhor | Reducing molar density of the mixture [mol/m^3] |
\[ \alpha^0 = \displaystyle\sum_{i=1}^{N}x_i[\alpha^0_{oi}(\rho,T) + \ln x_i] \]
where in this case, we use the \(\alpha^0\) for the given fluid, which uses the inputs \(\tau_i\) and \(\delta_i\), so we do the conversion between mixture and component reduced states with
\[ \tau_i = \frac{T_{c,i}}{T} = \frac{\tau T_{c,i}}{T_r} \]
\[ \delta_i = \frac{\rho}{\rho_{c,i}} = \frac{\delta\rho_r}{\rho_{c,i}} \]
Definition at line 3281 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r\) (dimensionless)
Reimplemented from CoolProp::AbstractState.
Definition at line 3431 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Reimplemented in CoolProp::AbstractCubicBackend.
Definition at line 3275 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Build the spinodal curve.
Reimplemented from CoolProp::AbstractState.
Definition at line 4147 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the second virial coefficient.
Reimplemented from CoolProp::AbstractState.
Definition at line 1523 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Change the equation of state for one component.
Reimplemented from CoolProp::AbstractState.
Definition at line 401 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the chemical potential in J/mol.
Reimplemented from CoolProp::AbstractState.
Definition at line 3219 of file HelmholtzEOSMixtureBackend.cpp.
|
inlinevirtual |
Using this backend, calculate the compressibility factor Z \( Z = p/(\rho R T) \).
Reimplemented from CoolProp::AbstractState.
Definition at line 612 of file HelmholtzEOSMixtureBackend.h.
|
virtual |
Using this backend, calculate the thermal conductivity in W/m/K.
Reimplemented from CoolProp::AbstractState.
Definition at line 935 of file HelmholtzEOSMixtureBackend.cpp.
| CoolPropDbl CoolProp::HelmholtzEOSMixtureBackend::calc_conductivity_background | ( | void | ) |
Definition at line 919 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Calculate each of the contributions to the conductivity.
If the conductivity model is hardcoded or ECS is being used, there will only be one entry in initial_density and all others will be zero
Reimplemented from CoolProp::AbstractState.
Definition at line 807 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Calculate the conformal state (unity shape factors starting point if T < 0 and rhomolar < 0)
Reimplemented from CoolProp::AbstractState.
Definition at line 951 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the molar constant-pressure specific heat in J/mol/K.
Reimplemented from CoolProp::AbstractState.
Definition at line 3065 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the ideal gas molar constant-pressure specific heat in J/mol/K.
Reimplemented from CoolProp::AbstractState.
Definition at line 3086 of file HelmholtzEOSMixtureBackend.cpp.
| CoolProp::CriticalState CoolProp::HelmholtzEOSMixtureBackend::calc_critical_point | ( | double | rho0, |
| double | T0 | ||
| ) |
Not used, for testing purposes
Definition at line 3801 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Calculate the values \(\mathcal{L}_1^*\) and \(\mathcal{M}_1^*\).
Reimplemented from CoolProp::AbstractState.
Definition at line 4065 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the third virial coefficient.
Reimplemented from CoolProp::AbstractState.
Definition at line 1531 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the molar constant-volume specific heat in J/mol/K.
Reimplemented from CoolProp::AbstractState.
Definition at line 3050 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\delta\tau}\) (dimensionless)
Reimplemented from CoolProp::AbstractState.
Definition at line 3505 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\delta\delta}\) (dimensionless)
Reimplemented from CoolProp::AbstractState.
Definition at line 3509 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\tau\tau}\) (dimensionless)
Reimplemented from CoolProp::AbstractState.
Definition at line 3513 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\delta}\) (dimensionless)
Reimplemented from CoolProp::AbstractState.
Definition at line 3451 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\tau}\) (dimensionless)
Reimplemented from CoolProp::AbstractState.
Definition at line 3447 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\tau\tau}\) (dimensionless)
Reimplemented from CoolProp::AbstractState.
Definition at line 3443 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\delta\delta\tau}\) (dimensionless)
Reimplemented from CoolProp::AbstractState.
Definition at line 3521 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\delta\delta\delta}\) (dimensionless)
Reimplemented from CoolProp::AbstractState.
Definition at line 3517 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\delta\tau\tau}\) (dimensionless)
Reimplemented from CoolProp::AbstractState.
Definition at line 3525 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\tau\tau\tau}\) (dimensionless)
Reimplemented from CoolProp::AbstractState.
Definition at line 3529 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\delta\tau}\) (dimensionless)
Reimplemented from CoolProp::AbstractState.
Definition at line 3459 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\delta\delta}\) (dimensionless)
Reimplemented from CoolProp::AbstractState.
Definition at line 3455 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\tau\tau}\) (dimensionless)
Reimplemented from CoolProp::AbstractState.
Definition at line 3463 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\tau\tau\tau}\) (dimensionless)
Reimplemented from CoolProp::AbstractState.
Definition at line 3467 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\delta\tau\tau}\) (dimensionless)
Reimplemented from CoolProp::AbstractState.
Definition at line 3480 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\delta\delta\tau}\) (dimensionless)
Reimplemented from CoolProp::AbstractState.
Definition at line 3476 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\delta\delta\delta}\) (dimensionless)
Reimplemented from CoolProp::AbstractState.
Definition at line 3472 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\tau\tau\tau}\) (dimensionless)
Reimplemented from CoolProp::AbstractState.
Definition at line 3484 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\tau\tau\tau\tau}\) (dimensionless)
Reimplemented from CoolProp::AbstractState.
Definition at line 3488 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\delta}\) (dimensionless)
Reimplemented from CoolProp::AbstractState.
Definition at line 3497 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\tau}\) (dimensionless)
Reimplemented from CoolProp::AbstractState.
Definition at line 3501 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta}\) (dimensionless)
Reimplemented from CoolProp::AbstractState.
Definition at line 3435 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\tau}\) (dimensionless)
Reimplemented from CoolProp::AbstractState.
Definition at line 3439 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the derivative dB/dT.
Reimplemented from CoolProp::AbstractState.
Definition at line 1526 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the derivative dC/dT.
Reimplemented from CoolProp::AbstractState.
Definition at line 1534 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate and cache the excess properties.
Reimplemented from CoolProp::AbstractState.
Definition at line 3170 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Reimplemented from CoolProp::AbstractState.
Definition at line 3556 of file HelmholtzEOSMixtureBackend.cpp.
| CoolPropDbl CoolProp::HelmholtzEOSMixtureBackend::calc_first_saturation_deriv | ( | parameters | Of1, |
| parameters | Wrt1, | ||
| HelmholtzEOSMixtureBackend & | SatL, | ||
| HelmholtzEOSMixtureBackend & | SatV | ||
| ) |
Definition at line 3533 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Reimplemented from CoolProp::AbstractState.
Definition at line 3655 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Reimplemented from CoolProp::AbstractState.
Definition at line 3687 of file HelmholtzEOSMixtureBackend.cpp.
|
inlinevirtual |
Using this backend, calculate the flame hazard.
Reimplemented from CoolProp::AbstractState.
Definition at line 455 of file HelmholtzEOSMixtureBackend.h.
|
virtual |
Using this backend, get a vector of fluid names.
Reimplemented from CoolProp::AbstractState.
Definition at line 1023 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the fugacity in Pa.
Reimplemented from CoolProp::AbstractState.
Definition at line 3215 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the fugacity coefficient (dimensionless)
Reimplemented from CoolProp::AbstractState.
Reimplemented in CoolProp::VTPRBackend.
Definition at line 3211 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the universal gas constant \(R_u\) in J/mol/K.
Reimplemented from CoolProp::AbstractState.
Definition at line 515 of file HelmholtzEOSMixtureBackend.cpp.
|
inlinevirtual |
Get the data from the spinodal curve.
Reimplemented from CoolProp::AbstractState.
Definition at line 259 of file HelmholtzEOSMixtureBackend.h.
|
virtual |
Using this backend, calculate the molar Gibbs function in J/mol.
Reimplemented from CoolProp::AbstractState.
Definition at line 3146 of file HelmholtzEOSMixtureBackend.cpp.
| CoolPropDbl CoolProp::HelmholtzEOSMixtureBackend::calc_gibbsmolar_nocache | ( | CoolPropDbl | T, |
| CoolPropDbl | rhomolar | ||
| ) |
Definition at line 3133 of file HelmholtzEOSMixtureBackend.cpp.
|
inlinevirtual |
Using this backend, calculate the residual molar Gibbs function in J/mol.
Reimplemented from CoolProp::AbstractState.
Definition at line 418 of file HelmholtzEOSMixtureBackend.h.
|
virtual |
Using this backend, calculate the 100-year global warming potential (GWP)
Reimplemented from CoolProp::AbstractState.
Definition at line 1052 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the 20-year global warming potential (GWP)
Reimplemented from CoolProp::AbstractState.
Definition at line 1041 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the 500-year global warming potential (GWP)
Reimplemented from CoolProp::AbstractState.
Definition at line 1063 of file HelmholtzEOSMixtureBackend.cpp.
|
inlinevirtual |
Using this backend, calculate the health hazard.
Reimplemented from CoolProp::AbstractState.
Definition at line 459 of file HelmholtzEOSMixtureBackend.h.
|
virtual |
Using this backend, calculate the molar Helmholtz energy in J/mol.
Reimplemented from CoolProp::AbstractState.
Definition at line 3188 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the molar enthalpy in J/mol.
Reimplemented from CoolProp::AbstractState.
Definition at line 2926 of file HelmholtzEOSMixtureBackend.cpp.
| CoolPropDbl CoolProp::HelmholtzEOSMixtureBackend::calc_hmolar_nocache | ( | CoolPropDbl | T, |
| CoolPropDbl | rhomolar | ||
| ) |
Definition at line 2911 of file HelmholtzEOSMixtureBackend.cpp.
|
inlinevirtual |
Using this backend, calculate the residual molar enthalpy in J/mol.
Reimplemented from CoolProp::AbstractState.
Definition at line 433 of file HelmholtzEOSMixtureBackend.h.
| void CoolProp::HelmholtzEOSMixtureBackend::calc_hsat_max | ( | void | ) |
Definition at line 1972 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Reimplemented from CoolProp::AbstractState.
Definition at line 1006 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Reimplemented from CoolProp::AbstractState.
Definition at line 1440 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Reimplemented from CoolProp::AbstractState.
Definition at line 585 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the molar mass in kg/mol.
Reimplemented from CoolProp::AbstractState.
Definition at line 531 of file HelmholtzEOSMixtureBackend.cpp.
|
inlinevirtual |
Reimplemented from CoolProp::AbstractState.
Definition at line 324 of file HelmholtzEOSMixtureBackend.h.
|
inlinevirtual |
Reimplemented from CoolProp::AbstractState.
Definition at line 327 of file HelmholtzEOSMixtureBackend.h.
|
virtual |
Using this backend, get the name of the fluid.
Reimplemented from CoolProp::AbstractState.
Definition at line 987 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the ozone depletion potential (ODP)
Reimplemented from CoolProp::AbstractState.
Definition at line 1030 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, get the critical point pressure in Pa.
Reimplemented from CoolProp::AbstractState.
Reimplemented in CoolProp::AbstractCubicBackend.
Definition at line 1093 of file HelmholtzEOSMixtureBackend.cpp.
|
inlinevirtual |
Using this backend, get the reducing point pressure in Pa.
Reimplemented from CoolProp::AbstractState.
Definition at line 553 of file HelmholtzEOSMixtureBackend.h.
|
virtual |
Using this backend, get the triple point pressure in Pa.
Reimplemented from CoolProp::AbstractState.
Definition at line 980 of file HelmholtzEOSMixtureBackend.cpp.
|
inlinevirtual |
Using this backend, calculate the phase.
Reimplemented from CoolProp::AbstractState.
Definition at line 205 of file HelmholtzEOSMixtureBackend.h.
|
virtual |
Using this backend, construct the phase envelope, the variable type describes the type of phase envelope to be built.
Reimplemented from CoolProp::AbstractState.
Definition at line 449 of file HelmholtzEOSMixtureBackend.cpp.
|
inlinevirtual |
Reimplemented from CoolProp::AbstractState.
Definition at line 333 of file HelmholtzEOSMixtureBackend.h.
| CoolPropDbl CoolProp::HelmholtzEOSMixtureBackend::calc_phase_identification_parameter | ( | void | ) |
The phase identification parameter of Venkatarathnam et al., FPE, 2011.
Definition at line 3228 of file HelmholtzEOSMixtureBackend.cpp.
|
inlinevirtual |
Using this backend, calculate the physical hazard.
Reimplemented from CoolProp::AbstractState.
Definition at line 463 of file HelmholtzEOSMixtureBackend.h.
|
inlinevirtual |
Using this backend, calculate the phase identification parameter (PIP)
Reimplemented from CoolProp::AbstractState.
Definition at line 558 of file HelmholtzEOSMixtureBackend.h.
|
virtual |
Using this backend, calculate the maximum pressure in Pa.
Reimplemented from CoolProp::AbstractState.
Definition at line 1196 of file HelmholtzEOSMixtureBackend.cpp.
| CoolPropDbl CoolProp::HelmholtzEOSMixtureBackend::calc_pmax_sat | ( | void | ) |
Definition at line 1131 of file HelmholtzEOSMixtureBackend.cpp.
| void CoolProp::HelmholtzEOSMixtureBackend::calc_pmin_sat | ( | CoolPropDbl & | pmin_satL, |
| CoolPropDbl & | pmin_satV | ||
| ) |
Definition at line 1169 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the pressure in Pa.
Reimplemented from CoolProp::AbstractState.
Definition at line 2892 of file HelmholtzEOSMixtureBackend.cpp.
| CoolPropDbl CoolProp::HelmholtzEOSMixtureBackend::calc_pressure_nocache | ( | CoolPropDbl | T, |
| CoolPropDbl | rhomolar | ||
| ) |
Definition at line 2456 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Reimplemented from CoolProp::AbstractState.
Definition at line 3245 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Reimplemented in CoolProp::AbstractCubicBackend.
Definition at line 3235 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, get the critical point molar density in mol/m^3.
Reimplemented from CoolProp::AbstractState.
Reimplemented in CoolProp::AbstractCubicBackend.
Definition at line 1112 of file HelmholtzEOSMixtureBackend.cpp.
|
inlinevirtual |
Using this backend, get the reducing point molar density in mol/m^3.
Reimplemented from CoolProp::AbstractState.
Definition at line 550 of file HelmholtzEOSMixtureBackend.h.
|
virtual |
Reimplemented from CoolProp::AbstractState.
Definition at line 995 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Reimplemented from CoolProp::AbstractState.
Definition at line 1000 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
| 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 from CoolProp::AbstractState.
Definition at line 538 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Reimplemented from CoolProp::AbstractState.
Definition at line 3580 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Reimplemented from CoolProp::AbstractState.
Definition at line 3608 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the molar entropy in J/mol/K.
Reimplemented from CoolProp::AbstractState.
Definition at line 2974 of file HelmholtzEOSMixtureBackend.cpp.
| CoolPropDbl CoolProp::HelmholtzEOSMixtureBackend::calc_smolar_nocache | ( | CoolPropDbl | T, |
| CoolPropDbl | rhomolar | ||
| ) |
Definition at line 2958 of file HelmholtzEOSMixtureBackend.cpp.
|
inlinevirtual |
Using this backend, calculate the residual molar entropy in J/mol/K.
Reimplemented from CoolProp::AbstractState.
Definition at line 427 of file HelmholtzEOSMixtureBackend.h.
|
inlinevirtual |
Specify the phase - this phase will always be used in calculations.
| phase_index | The index from CoolProp::phases |
Reimplemented from CoolProp::AbstractState.
Definition at line 213 of file HelmholtzEOSMixtureBackend.h.
|
virtual |
Using this backend, calculate the speed of sound in m/s.
Reimplemented from CoolProp::AbstractState.
Definition at line 3098 of file HelmholtzEOSMixtureBackend.cpp.
| void CoolProp::HelmholtzEOSMixtureBackend::calc_ssat_max | ( | void | ) |
Definition at line 1937 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
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 from CoolProp::AbstractState.
Definition at line 481 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the surface tension in N/m.
Reimplemented from CoolProp::AbstractState.
Definition at line 592 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, get the critical point temperature in K.
Reimplemented from CoolProp::AbstractState.
Reimplemented in CoolProp::AbstractCubicBackend.
Definition at line 1074 of file HelmholtzEOSMixtureBackend.cpp.
|
inlinevirtual |
Using this backend, get the reducing point temperature in K.
Reimplemented from CoolProp::AbstractState.
Definition at line 547 of file HelmholtzEOSMixtureBackend.h.
|
virtual |
Calculate tangent plane distance for given trial composition w.
Reimplemented from CoolProp::AbstractState.
Definition at line 4126 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the maximum temperature in K.
Reimplemented from CoolProp::AbstractState.
Definition at line 1182 of file HelmholtzEOSMixtureBackend.cpp.
| CoolPropDbl CoolProp::HelmholtzEOSMixtureBackend::calc_Tmax_sat | ( | void | ) |
Definition at line 1142 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the minimum temperature in K.
Reimplemented from CoolProp::AbstractState.
Definition at line 1189 of file HelmholtzEOSMixtureBackend.cpp.
| void CoolProp::HelmholtzEOSMixtureBackend::calc_Tmin_sat | ( | CoolPropDbl & | Tmin_satL, |
| CoolPropDbl & | Tmin_satV | ||
| ) |
Definition at line 1159 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, get the triple point temperature in K.
Reimplemented from CoolProp::AbstractState.
Definition at line 973 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Using this backend, calculate the molar internal energy in J/mol.
Reimplemented from CoolProp::AbstractState.
Definition at line 3021 of file HelmholtzEOSMixtureBackend.cpp.
| CoolPropDbl CoolProp::HelmholtzEOSMixtureBackend::calc_umolar_nocache | ( | CoolPropDbl | T, |
| CoolPropDbl | rhomolar | ||
| ) |
Definition at line 3008 of file HelmholtzEOSMixtureBackend.cpp.
|
inlinevirtual |
Unspecify the phase - the phase is no longer imposed, different solvers can do as they like.
Reimplemented from CoolProp::AbstractState.
Definition at line 219 of file HelmholtzEOSMixtureBackend.h.
|
virtual |
Using this backend, calculate the viscosity in Pa-s.
Reimplemented from CoolProp::AbstractState.
Definition at line 699 of file HelmholtzEOSMixtureBackend.cpp.
| CoolPropDbl CoolProp::HelmholtzEOSMixtureBackend::calc_viscosity_background | ( | CoolPropDbl | eta_dilute, |
| CoolPropDbl & | initial_density, | ||
| CoolPropDbl & | residual | ||
| ) |
Definition at line 644 of file HelmholtzEOSMixtureBackend.cpp.
| CoolPropDbl CoolProp::HelmholtzEOSMixtureBackend::calc_viscosity_background | ( | void | ) |
Definition at line 640 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Calculate each of the contributions to the viscosity.
If the viscosity model is hardcoded or ECS is being used, there will only be one entry in critical and all others will be zero
Reimplemented from CoolProp::AbstractState.
Definition at line 715 of file HelmholtzEOSMixtureBackend.cpp.
| CoolPropDbl CoolProp::HelmholtzEOSMixtureBackend::calc_viscosity_dilute | ( | void | ) |
Definition at line 603 of file HelmholtzEOSMixtureBackend.cpp.
|
inlinevirtual |
Clear all the cached values.
Bulk values
Reimplemented from CoolProp::AbstractState.
Definition at line 134 of file HelmholtzEOSMixtureBackend.h.
| void CoolProp::HelmholtzEOSMixtureBackend::DmolarP_phase_determination | ( | ) |
|
virtual |
Return a string from the backend for the mixture/fluid.
Reimplemented from CoolProp::AbstractState.
Definition at line 199 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Get binary mixture double value.
Get binary mixture floating point parameter for this instance.
Reimplemented from CoolProp::AbstractState.
Reimplemented in CoolProp::VTPRBackend, and CoolProp::AbstractCubicBackend.
Definition at line 355 of file HelmholtzEOSMixtureBackend.cpp.
|
inline |
Definition at line 308 of file HelmholtzEOSMixtureBackend.h.
|
inline |
Definition at line 305 of file HelmholtzEOSMixtureBackend.h.
|
virtual |
Reimplemented in CoolProp::SRKBackend, CoolProp::PengRobinsonBackend, CoolProp::VTPRBackend, and CoolProp::AbstractCubicBackend.
Definition at line 144 of file HelmholtzEOSMixtureBackend.cpp.
|
inlinevirtual |
Checking function to see if we should stop the tracing of the critical contour.
Reimplemented in CoolProp::AbstractCubicBackend.
Definition at line 251 of file HelmholtzEOSMixtureBackend.h.
|
virtual |
Get the search radius in delta and tau for the tracer.
Reimplemented in CoolProp::AbstractCubicBackend.
Definition at line 4072 of file HelmholtzEOSMixtureBackend.cpp.
|
inlinevirtual |
Reimplemented in CoolProp::AbstractCubicBackend.
Definition at line 244 of file HelmholtzEOSMixtureBackend.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 from CoolProp::AbstractState.
Reimplemented in CoolProp::AbstractCubicBackend.
Definition at line 277 of file HelmholtzEOSMixtureBackend.h.
|
virtual |
Double fluid parameter (currently the volume translation parameter for cubic)
Reimplemented from CoolProp::AbstractState.
Reimplemented in CoolProp::AbstractCubicBackend.
Definition at line 273 of file HelmholtzEOSMixtureBackend.cpp.
|
inline |
Definition at line 311 of file HelmholtzEOSMixtureBackend.h.
|
inline |
Definition at line 314 of file HelmholtzEOSMixtureBackend.h.
|
inlinevirtual |
Get the mole fractions of the fluid.
Implements CoolProp::AbstractState.
Definition at line 387 of file HelmholtzEOSMixtureBackend.h.
|
inline |
Definition at line 393 of file HelmholtzEOSMixtureBackend.h.
|
inline |
Definition at line 390 of file HelmholtzEOSMixtureBackend.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 from CoolProp::AbstractState.
Definition at line 605 of file HelmholtzEOSMixtureBackend.h.
|
inline |
Definition at line 317 of file HelmholtzEOSMixtureBackend.h.
|
inline |
Definition at line 320 of file HelmholtzEOSMixtureBackend.h.
|
protected |
Definition at line 266 of file HelmholtzEOSMixtureBackend.cpp.
|
inlinevirtual |
Return true if the fluid has a melting line - default is false, but can be re-implemented by derived class.
Reimplemented from CoolProp::AbstractState.
Definition at line 170 of file HelmholtzEOSMixtureBackend.h.
|
inline |
Definition at line 167 of file HelmholtzEOSMixtureBackend.h.
| void CoolProp::HelmholtzEOSMixtureBackend::p_phase_determination_pure_or_pseudopure | ( | int | other, |
| CoolPropDbl | value, | ||
| bool & | saturation_called | ||
| ) |
Definition at line 1539 of file HelmholtzEOSMixtureBackend.cpp.
|
protected |
Definition at line 1487 of file HelmholtzEOSMixtureBackend.cpp.
|
protected |
Definition at line 1313 of file HelmholtzEOSMixtureBackend.cpp.
| void CoolProp::HelmholtzEOSMixtureBackend::recalculate_singlephase_phase | ( | ) |
Calculate the phase once the state is fully calculated but you aren't sure if it is liquid or gas or ...
Definition at line 179 of file HelmholtzEOSMixtureBackend.cpp.
| void CoolProp::HelmholtzEOSMixtureBackend::resize | ( | std::size_t | N | ) |
Definition at line 170 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Set binary mixture floating point parameter.
Set binary mixture floating point parameter for this instance.
Also set the parameters in the managed pointers for other states
Reimplemented from CoolProp::AbstractState.
Reimplemented in CoolProp::VTPRBackend, and CoolProp::AbstractCubicBackend.
Definition at line 331 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Set a binary interaction string.
Set binary mixture floating point parameter for this instance.
Also set the parameters in the managed pointers for other states
Reimplemented from CoolProp::AbstractState.
Definition at line 377 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Set the components of the mixture.
| components | The components that are to be used in this mixture |
| generate_SatL_and_SatV | true if SatL and SatV classes should be added, false otherwise. Added so that saturation classes can be added without infinite recursion of adding saturation classes |
Definition at line 92 of file HelmholtzEOSMixtureBackend.cpp.
|
inlinevirtual |
Set the cubic alpha function's constants:
Reimplemented from CoolProp::AbstractState.
Reimplemented in CoolProp::AbstractCubicBackend.
Definition at line 195 of file HelmholtzEOSMixtureBackend.h.
|
staticprotected |
Definition at line 4237 of file HelmholtzEOSMixtureBackend.cpp.
|
inlinevirtual |
Set fluid parameter (currently the volume translation parameter for cubic)
Reimplemented from CoolProp::AbstractState.
Reimplemented in CoolProp::AbstractCubicBackend.
Definition at line 200 of file HelmholtzEOSMixtureBackend.h.
|
virtual |
Set the mass fractions.
| mass_fractions | The vector of mass fractions of the components |
Implements CoolProp::AbstractState.
Definition at line 152 of file HelmholtzEOSMixtureBackend.cpp.
| void CoolProp::HelmholtzEOSMixtureBackend::set_mixture_parameters | ( | ) |
Set the mixture parameters - binary pair reducing functions, departure functions, F_ij, etc.
Definition at line 457 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Set the mole fractions.
| mole_fractions | The vector of mole fractions of the components |
Implements CoolProp::AbstractState.
Definition at line 125 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Set the reference state based on a thermodynamic state point specified by temperature and molar density.
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 from CoolProp::AbstractState.
Definition at line 4219 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
brief Set the reference state based on a string representation
Reimplemented from CoolProp::AbstractState.
Definition at line 4153 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
The residual to be used to find the location where dpdrho=0 for given T
Definition at line 2468 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Reimplemented in CoolProp::AbstractCubicBackend.
Definition at line 2707 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Reimplemented in CoolProp::AbstractCubicBackend.
Definition at line 2636 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Reimplemented in CoolProp::AbstractCubicBackend.
Definition at line 2819 of file HelmholtzEOSMixtureBackend.cpp.
| CoolPropDbl CoolProp::HelmholtzEOSMixtureBackend::SRK_covolume | ( | ) |
Definition at line 2626 of file HelmholtzEOSMixtureBackend.cpp.
| void CoolProp::HelmholtzEOSMixtureBackend::sync_linked_states | ( | const HelmholtzEOSMixtureBackend * const | source | ) |
Definition at line 134 of file HelmholtzEOSMixtureBackend.cpp.
| void CoolProp::HelmholtzEOSMixtureBackend::T_phase_determination_pure_or_pseudopure | ( | int | other, |
| CoolPropDbl | value | ||
| ) |
Definition at line 1998 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
The standard update function.
| input_pair | The pair of inputs that will be provided |
| value1 | The first input value |
| value2 | The second input value |
Implements CoolProp::AbstractState.
Reimplemented in CoolProp::AbstractCubicBackend.
Definition at line 1332 of file HelmholtzEOSMixtureBackend.cpp.
| void CoolProp::HelmholtzEOSMixtureBackend::update_DmolarT_direct | ( | CoolPropDbl | rhomolar, |
| CoolPropDbl | T | ||
| ) |
Definition at line 1249 of file HelmholtzEOSMixtureBackend.cpp.
| void CoolProp::HelmholtzEOSMixtureBackend::update_HmolarQ_with_guessT | ( | CoolPropDbl | hmolar, |
| CoolPropDbl | Q, | ||
| CoolPropDbl | Tguess | ||
| ) |
Definition at line 1275 of file HelmholtzEOSMixtureBackend.cpp.
| void CoolProp::HelmholtzEOSMixtureBackend::update_internal | ( | HelmholtzEOSMixtureBackend & | HEOS | ) |
Update all the internal variables for a state by copying from another state.
Definition at line 1287 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Update the state for QT inputs for pure fluids when using the superancillary functions.
Reimplemented from CoolProp::AbstractState.
Definition at line 1204 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Update the states after having changed the reference state for enthalpy and entropy.
Reimplemented from CoolProp::AbstractState.
Definition at line 461 of file HelmholtzEOSMixtureBackend.cpp.
| void CoolProp::HelmholtzEOSMixtureBackend::update_TDmolarP_unchecked | ( | CoolPropDbl | T, |
| CoolPropDbl | rhomolarL, | ||
| CoolPropDbl | p | ||
| ) |
Definition at line 1236 of file HelmholtzEOSMixtureBackend.cpp.
| void CoolProp::HelmholtzEOSMixtureBackend::update_TP_guessrho | ( | CoolPropDbl | T, |
| CoolPropDbl | p, | ||
| CoolPropDbl | rho_guess | ||
| ) |
Update with TP and a guess for rho.
| T | Temperature in K |
| p | Pressure in Pa |
| rho_guess | Density in mol/m^3 guessed |
Definition at line 1299 of file HelmholtzEOSMixtureBackend.cpp.
|
virtual |
Update the state using guess values.
Reimplemented from CoolProp::AbstractState.
Definition at line 1452 of file HelmholtzEOSMixtureBackend.cpp.
|
inlinevirtual |
Implements CoolProp::AbstractState.
Definition at line 161 of file HelmholtzEOSMixtureBackend.h.
|
inlinevirtual |
Implements CoolProp::AbstractState.
Definition at line 158 of file HelmholtzEOSMixtureBackend.h.
|
inlinevirtual |
Implements CoolProp::AbstractState.
Definition at line 164 of file HelmholtzEOSMixtureBackend.h.
|
friend |
Definition at line 155 of file HelmholtzEOSMixtureBackend.h.
|
friend |
Definition at line 145 of file HelmholtzEOSMixtureBackend.h.
|
friend |
Definition at line 149 of file HelmholtzEOSMixtureBackend.h.
|
friend |
Definition at line 153 of file HelmholtzEOSMixtureBackend.h.
|
friend |
Definition at line 151 of file HelmholtzEOSMixtureBackend.h.
|
friend |
Definition at line 147 of file HelmholtzEOSMixtureBackend.h.
|
protected |
Definition at line 104 of file HelmholtzEOSMixtureBackend.h.
|
protected |
The components that are in use.
Definition at line 98 of file HelmholtzEOSMixtureBackend.h.
|
protected |
A temporary state used for calculations of the critical point(s)
Definition at line 72 of file HelmholtzEOSMixtureBackend.h.
| SimpleState CoolProp::HelmholtzEOSMixtureBackend::hsat_max |
Definition at line 130 of file HelmholtzEOSMixtureBackend.h.
|
protected |
A flag for whether the substance is a pure or pseudo-pure fluid (true) or a mixture (false)
Definition at line 99 of file HelmholtzEOSMixtureBackend.h.
|
protected |
The K factors for the components.
Definition at line 101 of file HelmholtzEOSMixtureBackend.h.
|
protected |
States that are linked to this one, and should be updated (BIP, reference state, etc.)
Definition at line 69 of file HelmholtzEOSMixtureBackend.h.
|
protected |
The natural logarithms of the K factors of the components.
Definition at line 102 of file HelmholtzEOSMixtureBackend.h.
|
protected |
The bulk mole fractions of the mixture.
Definition at line 100 of file HelmholtzEOSMixtureBackend.h.
|
protected |
Number of components.
Definition at line 105 of file HelmholtzEOSMixtureBackend.h.
| PhaseEnvelopeData CoolProp::HelmholtzEOSMixtureBackend::PhaseEnvelope |
Definition at line 129 of file HelmholtzEOSMixtureBackend.h.
| shared_ptr<ReducingFunction> CoolProp::HelmholtzEOSMixtureBackend::Reducing |
Definition at line 127 of file HelmholtzEOSMixtureBackend.h.
| shared_ptr<ResidualHelmholtz> CoolProp::HelmholtzEOSMixtureBackend::residual_helmholtz |
Definition at line 128 of file HelmholtzEOSMixtureBackend.h.
| shared_ptr<HelmholtzEOSMixtureBackend> CoolProp::HelmholtzEOSMixtureBackend::SatL |
Definition at line 341 of file HelmholtzEOSMixtureBackend.h.
| shared_ptr<HelmholtzEOSMixtureBackend> CoolProp::HelmholtzEOSMixtureBackend::SatV |
Definition at line 341 of file HelmholtzEOSMixtureBackend.h.
| SpinodalData CoolProp::HelmholtzEOSMixtureBackend::spinodal_values |
Definition at line 132 of file HelmholtzEOSMixtureBackend.h.
| SsatSimpleState CoolProp::HelmholtzEOSMixtureBackend::ssat_max |
Definition at line 131 of file HelmholtzEOSMixtureBackend.h.
|
protected |
A temporary state used for calculations of the tangent-plane-distance.
Definition at line 71 of file HelmholtzEOSMixtureBackend.h.
|
protected |
A temporary state used for calculations of pure fluid properties.
Definition at line 70 of file HelmholtzEOSMixtureBackend.h.