8 #ifndef ABSTRACTSTATE_H_ 9 #define ABSTRACTSTATE_H_ 11 #include "CachedElement.h" 12 #include "Exceptions.h" 13 #include "DataStructures.h" 14 #include "PhaseEnvelope.h" 15 #include "crossplatform_shared_ptr.h" 26 std::vector<double>
tau,
43 std::vector<double> x,
86 bool isSupercriticalPhase(
void) {
90 bool isHomogeneousPhase(
void) {
94 bool isTwoPhase(
void) {
115 CachedElement _hmolar, _smolar, _umolar, _logp, _logrhomolar, _cpmolar, _cp0molar, _cvmolar, _speed_sound, _gibbsmolar, _helmholtzmolar;
132 CachedElement _alpha0, _dalpha0_dTau, _dalpha0_dDelta, _d2alpha0_dTau2, _d2alpha0_dDelta_dTau, _d2alpha0_dDelta2, _d3alpha0_dTau3,
133 _d3alpha0_dDelta_dTau2, _d3alpha0_dDelta2_dTau, _d3alpha0_dDelta3, _alphar, _dalphar_dTau, _dalphar_dDelta, _d2alphar_dTau2,
134 _d2alphar_dDelta_dTau, _d2alphar_dDelta2, _d3alphar_dTau3, _d3alphar_dDelta_dTau2, _d3alphar_dDelta2_dTau, _d3alphar_dDelta3, _d4alphar_dTau4,
135 _d4alphar_dDelta_dTau3, _d4alphar_dDelta2_dTau2, _d4alphar_dDelta3_dTau, _d4alphar_dDelta4;
137 CachedElement _dalphar_dDelta_lim, _d2alphar_dDelta2_lim, _d2alphar_dDelta_dTau_lim, _d3alphar_dDelta2_dTau_lim;
196 throw NotImplementedError(
"calc_isothermal_compressibility is not implemented for this backend");
200 throw NotImplementedError(
"calc_isobaric_expansion_coefficient is not implemented for this backend");
204 throw NotImplementedError(
"calc_isentropic_expansion_coefficient is not implemented for this backend");
236 throw NotImplementedError(
"calc_fugacity_coefficient is not implemented for this backend");
240 throw NotImplementedError(
"calc_fugacity_coefficients is not implemented for this backend");
280 throw NotImplementedError(
"calc_d2alphar_dDelta_dTau is not implemented for this backend");
292 throw NotImplementedError(
"calc_d3alphar_dDelta2_dTau is not implemented for this backend");
296 throw NotImplementedError(
"calc_d3alphar_dDelta_dTau2 is not implemented for this backend");
309 throw NotImplementedError(
"calc_d4alphar_dDelta3_dTau is not implemented for this backend");
313 throw NotImplementedError(
"calc_d4alphar_dDelta2_dTau2 is not implemented for this backend");
317 throw NotImplementedError(
"calc_d4alphar_dDelta_dTau3 is not implemented for this backend");
339 throw NotImplementedError(
"calc_d2alpha0_dDelta_dTau is not implemented for this backend");
355 throw NotImplementedError(
"calc_d3alpha0_dDelta2_dTau is not implemented for this backend");
359 throw NotImplementedError(
"calc_d3alpha0_dDelta_dTau2 is not implemented for this backend");
366 virtual void calc_reducing_state(
void) {
427 throw NotImplementedError(
"calc_reciprocal_reduced_temperature is not implemented for this backend");
448 throw NotImplementedError(
"calc_compressibility_factor is not implemented for this backend");
502 virtual CoolPropDbl calc_rhomass(
void) {
503 return rhomolar() * molar_mass();
505 virtual CoolPropDbl calc_hmass(
void) {
506 return hmolar() / molar_mass();
508 virtual CoolPropDbl calc_hmass_excess(
void) {
509 return hmolar_excess() / molar_mass();
511 virtual CoolPropDbl calc_smass(
void) {
512 return smolar() / molar_mass();
514 virtual CoolPropDbl calc_smass_excess(
void) {
515 return smolar_excess() / molar_mass();
517 virtual CoolPropDbl calc_cpmass(
void) {
518 return cpmolar() / molar_mass();
520 virtual CoolPropDbl calc_cp0mass(
void) {
521 return cp0molar() / molar_mass();
523 virtual CoolPropDbl calc_cvmass(
void) {
524 return cvmolar() / molar_mass();
526 virtual CoolPropDbl calc_umass(
void) {
527 return umolar() / molar_mass();
529 virtual CoolPropDbl calc_umass_excess(
void) {
530 return umolar_excess() / molar_mass();
532 virtual CoolPropDbl calc_gibbsmass(
void) {
533 return gibbsmolar() / molar_mass();
535 virtual CoolPropDbl calc_gibbsmass_excess(
void) {
536 return gibbsmolar_excess() / molar_mass();
538 virtual CoolPropDbl calc_helmholtzmass(
void) {
539 return helmholtzmolar() / molar_mass();
541 virtual CoolPropDbl calc_helmholtzmass_excess(
void) {
542 return helmholtzmolar_excess() / molar_mass();
544 virtual CoolPropDbl calc_volumemass_excess(
void) {
545 return volumemolar_excess() / molar_mass();
553 virtual CoolPropDbl calc_melting_line(
int param,
int given, CoolPropDbl value) {
575 throw NotImplementedError(
"This backend does not implement calc_unspecify_phase function");
591 virtual std::vector<CoolPropDbl> calc_mole_fractions_liquid(
void) {
592 throw NotImplementedError(
"calc_mole_fractions_liquid is not implemented for this backend");
594 virtual std::vector<CoolPropDbl> calc_mole_fractions_vapor(
void) {
595 throw NotImplementedError(
"calc_mole_fractions_vapor is not implemented for this backend");
597 virtual const std::vector<CoolPropDbl> calc_mass_fractions(
void) {
609 virtual CoolPropDbl calc_T_freeze(
void) {
614 throw NotImplementedError(
"calc_first_saturation_deriv is not implemented for this backend");
617 throw NotImplementedError(
"calc_second_saturation_deriv is not implemented for this backend");
620 throw NotImplementedError(
"calc_first_two_phase_deriv is not implemented for this backend");
623 throw NotImplementedError(
"calc_second_two_phase_deriv is not implemented for this backend");
626 throw NotImplementedError(
"calc_first_two_phase_deriv_splined is not implemented for this backend");
629 virtual CoolPropDbl calc_saturated_liquid_keyed_output(
parameters key) {
630 throw NotImplementedError(
"calc_saturated_liquid_keyed_output is not implemented for this backend");
632 virtual CoolPropDbl calc_saturated_vapor_keyed_output(
parameters key) {
633 throw NotImplementedError(
"calc_saturated_vapor_keyed_output is not implemented for this backend");
635 virtual void calc_ideal_curve(
const std::string& type, std::vector<double>& T, std::vector<double>& p) {
650 throw NotImplementedError(
"calc_tangent_plane_distance is not implemented for this backend");
658 virtual void calc_conformal_state(
const std::string& reference_fluid, CoolPropDbl& T, CoolPropDbl& rhomolar) {
662 virtual void calc_viscosity_contributions(CoolPropDbl& dilute, CoolPropDbl& initial_density, CoolPropDbl& residual, CoolPropDbl& critical) {
663 throw NotImplementedError(
"calc_viscosity_contributions is not implemented for this backend");
665 virtual void calc_conductivity_contributions(CoolPropDbl& dilute, CoolPropDbl& initial_density, CoolPropDbl& residual, CoolPropDbl& critical) {
666 throw NotImplementedError(
"calc_conductivity_contributions is not implemented for this backend");
668 virtual std::vector<CriticalState> calc_all_critical_points(
void) {
671 virtual void calc_build_spinodal() {
677 virtual void calc_criticality_contour_values(
double& L1star,
double& M1star) {
678 throw NotImplementedError(
"calc_criticality_contour_values is not implemented for this backend");
682 virtual void mass_to_molar_inputs(
CoolProp::input_pairs& input_pair, CoolPropDbl& value1, CoolPropDbl& value2);
703 return factory(backend, strsplit(fluid_names,
'&'));
729 static AbstractState* factory(
const std::string& backend,
const std::vector<std::string>& fluid_names);
740 virtual std::string backend_name(
void) = 0;
743 virtual bool using_mole_fractions(
void) = 0;
744 virtual bool using_mass_fractions(
void) = 0;
745 virtual bool using_volu_fractions(
void) = 0;
747 virtual void set_mole_fractions(
const std::vector<CoolPropDbl>& mole_fractions) = 0;
748 virtual void set_mass_fractions(
const std::vector<CoolPropDbl>& mass_fractions) = 0;
749 virtual void set_volu_fractions(
const std::vector<CoolPropDbl>& mass_fractions) {
774 "Setting reference state has not been implemented for this backend. Try using CoolProp::set_reference_stateD instead.");
784 "Setting reference state has not been implemented for this backend. Try using CoolProp::set_reference_stateD instead.");
787 #ifndef COOLPROPDBL_MAPS_TO_DOUBLE 788 void set_mole_fractions(
const std::vector<double>& mole_fractions) {
789 set_mole_fractions(std::vector<CoolPropDbl>(mole_fractions.begin(), mole_fractions.end()));
791 void set_mass_fractions(
const std::vector<double>& mass_fractions) {
792 set_mass_fractions(std::vector<CoolPropDbl>(mass_fractions.begin(), mass_fractions.end()));
794 void set_volu_fractions(
const std::vector<double>& volu_fractions) {
795 set_volu_fractions(std::vector<CoolPropDbl>(volu_fractions.begin(), volu_fractions.end()));
800 void set_mole_fractions_double(
const std::vector<double>& mole_fractions) {
801 set_mole_fractions(std::vector<CoolPropDbl>(mole_fractions.begin(), mole_fractions.end()));
807 return calc_mole_fractions_liquid();
811 std::vector<CoolPropDbl> x = calc_mole_fractions_liquid();
812 return std::vector<double>(x.begin(), x.end());
817 return calc_mole_fractions_vapor();
821 std::vector<CoolPropDbl> y = calc_mole_fractions_vapor();
822 return std::vector<double>(y.begin(), y.end());
826 virtual const std::vector<CoolPropDbl>& get_mole_fractions(
void) = 0;
829 return this->calc_mass_fractions();
854 std::vector<std::string> fluid_names(
void);
867 throw NotImplementedError(
"set_binary_interaction_double is not implemented for this backend");
871 throw NotImplementedError(
"set_binary_interaction_double is not implemented for this backend");
875 const std::string& value) {
876 throw NotImplementedError(
"set_binary_interaction_string is not implemented for this backend");
880 throw NotImplementedError(
"set_binary_interaction_string is not implemented for this backend");
884 throw NotImplementedError(
"get_binary_interaction_double is not implemented for this backend");
888 throw NotImplementedError(
"get_binary_interaction_double is not implemented for this backend");
892 throw NotImplementedError(
"get_binary_interaction_string is not implemented for this backend");
899 virtual void set_cubic_alpha_C(
const size_t i,
const std::string& parameter,
const double c1,
const double c2,
const double c3) {
900 throw ValueError(
"set_cubic_alpha_C only defined for cubic backends");
904 throw ValueError(
"set_fluid_parameter_double only defined for cubic backends");
908 throw ValueError(
"get_fluid_parameter_double only defined for cubic backends");
912 virtual bool clear();
914 virtual bool clear_comp_change();
925 return calc_state(state);
935 double Ttriple(
void);
943 calc_specify_phase(phase);
947 calc_unspecify_phase();
951 double T_critical(
void);
953 double p_critical(
void);
955 double rhomolar_critical(
void);
957 double rhomass_critical(
void);
961 return calc_all_critical_points();
966 calc_build_spinodal();
971 return calc_get_spinodal_data();
976 return calc_criticality_contour_values(L1star, M1star);
1002 double tangent_plane_distance(
const double T,
const double p,
const std::vector<double>& w,
const double rhomolar_guess = -1) {
1003 return calc_tangent_plane_distance(T, p, w, rhomolar_guess);
1007 double T_reducing(
void);
1009 double rhomolar_reducing(
void);
1011 double rhomass_reducing(
void);
1014 double p_triple(
void);
1022 return calc_description();
1027 return calc_dipole_moment();
1040 return calc_saturated_liquid_keyed_output(key);
1044 return calc_saturated_vapor_keyed_output(key);
1053 return calc_rhomolar();
1057 return calc_rhomass();
1072 double molar_mass(
void);
1074 double acentric_factor(
void);
1076 double gas_constant(
void);
1078 double Bvirial(
void);
1080 double dBvirial_dT(
void);
1082 double Cvirial(
void);
1084 double dCvirial_dT(
void);
1086 double compressibility_factor(
void);
1088 double hmolar(
void);
1090 double hmolar_residual(
void);
1093 return calc_hmass();
1096 double hmolar_excess(
void);
1099 return calc_hmass_excess();
1102 double smolar(
void);
1104 double smolar_residual(
void);
1107 return calc_smass();
1110 double smolar_excess(
void);
1113 return calc_smass_excess();
1116 double umolar(
void);
1119 return calc_umass();
1122 double umolar_excess(
void);
1125 return calc_umass_excess();
1128 double cpmolar(
void);
1131 return calc_cpmass();
1134 double cp0molar(
void);
1137 return calc_cp0mass();
1140 double cvmolar(
void);
1143 return calc_cvmass();
1146 double gibbsmolar(
void);
1148 double gibbsmolar_residual(
void);
1151 return calc_gibbsmass();
1154 double gibbsmolar_excess(
void);
1157 return calc_gibbsmass_excess();
1160 double helmholtzmolar(
void);
1163 return calc_helmholtzmass();
1166 double helmholtzmolar_excess(
void);
1169 return calc_helmholtzmass_excess();
1172 double volumemolar_excess(
void);
1175 return calc_volumemass_excess();
1178 double speed_sound(
void);
1180 double isothermal_compressibility(
void);
1182 double isobaric_expansion_coefficient(
void);
1184 double isentropic_expansion_coefficient(
void);
1186 double fugacity_coefficient(std::size_t i);
1188 std::vector<double> fugacity_coefficients();
1190 double fugacity(std::size_t i);
1192 double chemical_potential(std::size_t i);
1201 double fundamental_derivative_of_gas_dynamics(
void);
1209 calc_true_critical_point(T, rho);
1219 void ideal_curve(
const std::string& type, std::vector<double>& T, std::vector<double>& p) {
1220 calc_ideal_curve(type, T, p);
1232 return calc_first_partial_deriv(Of, Wrt, Constant);
1259 return calc_second_partial_deriv(Of1, Wrt1, Constant1, Wrt2, Constant2);
1284 return calc_first_saturation_deriv(Of1, Wrt1);
1305 return calc_second_saturation_deriv(Of1, Wrt1, Wrt2);
1328 return calc_first_two_phase_deriv(Of, Wrt, Constant);
1348 return calc_second_two_phase_deriv(Of, Wrt1, Constant1, Wrt2, Constant2);
1372 return calc_first_two_phase_deriv_splined(Of, Wrt, Constant, x_end);
1384 void build_phase_envelope(
const std::string& type =
"");
1389 return calc_phase_envelope_data();
1404 double melting_line(
int param,
int given,
double value);
1416 double viscosity(
void);
1419 calc_viscosity_contributions(dilute, initial_density, residual, critical);
1422 double conductivity(
void);
1425 calc_conductivity_contributions(dilute, initial_density, residual, critical);
1428 double surface_tension(
void);
1431 return cpmass() * viscosity() / conductivity();
1439 void conformal_state(
const std::string& reference_fluid, CoolPropDbl& T, CoolPropDbl& rhomolar) {
1440 return calc_conformal_state(reference_fluid, T, rhomolar);
1447 void change_EOS(
const std::size_t i,
const std::string& EOS_name) {
1448 calc_change_EOS(i, EOS_name);
1456 if (!_alpha0) _alpha0 = calc_alpha0();
1461 if (!_dalpha0_dDelta) _dalpha0_dDelta = calc_dalpha0_dDelta();
1462 return _dalpha0_dDelta;
1466 if (!_dalpha0_dTau) _dalpha0_dTau = calc_dalpha0_dTau();
1467 return _dalpha0_dTau;
1471 if (!_d2alpha0_dDelta2) _d2alpha0_dDelta2 = calc_d2alpha0_dDelta2();
1472 return _d2alpha0_dDelta2;
1476 if (!_d2alpha0_dDelta_dTau) _d2alpha0_dDelta_dTau = calc_d2alpha0_dDelta_dTau();
1477 return _d2alpha0_dDelta_dTau;
1481 if (!_d2alpha0_dTau2) _d2alpha0_dTau2 = calc_d2alpha0_dTau2();
1482 return _d2alpha0_dTau2;
1486 if (!_d3alpha0_dTau3) _d3alpha0_dTau3 = calc_d3alpha0_dTau3();
1487 return _d3alpha0_dTau3;
1491 if (!_d3alpha0_dDelta_dTau2) _d3alpha0_dDelta_dTau2 = calc_d3alpha0_dDelta_dTau2();
1492 return _d3alpha0_dDelta_dTau2;
1496 if (!_d3alpha0_dDelta2_dTau) _d3alpha0_dDelta2_dTau = calc_d3alpha0_dDelta2_dTau();
1497 return _d3alpha0_dDelta2_dTau;
1501 if (!_d3alpha0_dDelta3) _d3alpha0_dDelta3 = calc_d3alpha0_dDelta3();
1502 return _d3alpha0_dDelta3;
1507 if (!_alphar) _alphar = calc_alphar();
1512 if (!_dalphar_dDelta) _dalphar_dDelta = calc_dalphar_dDelta();
1513 return _dalphar_dDelta;
1517 if (!_dalphar_dTau) _dalphar_dTau = calc_dalphar_dTau();
1518 return _dalphar_dTau;
1522 if (!_d2alphar_dDelta2) _d2alphar_dDelta2 = calc_d2alphar_dDelta2();
1523 return _d2alphar_dDelta2;
1527 if (!_d2alphar_dDelta_dTau) _d2alphar_dDelta_dTau = calc_d2alphar_dDelta_dTau();
1528 return _d2alphar_dDelta_dTau;
1532 if (!_d2alphar_dTau2) _d2alphar_dTau2 = calc_d2alphar_dTau2();
1533 return _d2alphar_dTau2;
1537 if (!_d3alphar_dDelta3) _d3alphar_dDelta3 = calc_d3alphar_dDelta3();
1538 return _d3alphar_dDelta3;
1542 if (!_d3alphar_dDelta2_dTau) _d3alphar_dDelta2_dTau = calc_d3alphar_dDelta2_dTau();
1543 return _d3alphar_dDelta2_dTau;
1547 if (!_d3alphar_dDelta_dTau2) _d3alphar_dDelta_dTau2 = calc_d3alphar_dDelta_dTau2();
1548 return _d3alphar_dDelta_dTau2;
1552 if (!_d3alphar_dTau3) _d3alphar_dTau3 = calc_d3alphar_dTau3();
1553 return _d3alphar_dTau3;
1557 if (!_d4alphar_dDelta4) _d4alphar_dDelta4 = calc_d4alphar_dDelta4();
1558 return _d4alphar_dDelta4;
1562 if (!_d4alphar_dDelta3_dTau) _d4alphar_dDelta3_dTau = calc_d4alphar_dDelta3_dTau();
1563 return _d4alphar_dDelta3_dTau;
1567 if (!_d4alphar_dDelta2_dTau2) _d4alphar_dDelta2_dTau2 = calc_d4alphar_dDelta2_dTau2();
1568 return _d4alphar_dDelta2_dTau2;
1572 if (!_d4alphar_dDelta_dTau3) _d4alphar_dDelta_dTau3 = calc_d4alphar_dDelta_dTau3();
1573 return _d4alphar_dDelta_dTau3;
1577 if (!_d4alphar_dTau4) _d4alphar_dTau4 = calc_d4alphar_dTau4();
1578 return _d4alphar_dTau4;
1592 virtual AbstractState* get_AbstractState(
const std::vector<std::string>& fluid_names) = 0;
double T(void)
Return the temperature in K.
Definition: AbstractState.h:1048
phases imposed_phase_index
If the phase is imposed, the imposed phase index.
Definition: AbstractState.h:84
virtual CoolPropDbl calc_dalphar_dDelta(void)
Using this backend, calculate the residual Helmholtz energy term (dimensionless) ...
Definition: AbstractState.h:267
std::vector< double > delta
The reduced density ( )
Definition: AbstractState.h:26
CoolPropDbl dalphar_dDelta(void)
Return the term .
Definition: AbstractState.h:1511
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!!!)
Definition: AbstractState.h:891
double cvmass(void)
Return the mass constant volume specific heat in J/kg/K.
Definition: AbstractState.h:1142
virtual CoolPropDbl calc_hmolar_residual(void)
Using this backend, calculate the residual molar enthalpy in J/mol.
Definition: AbstractState.h:151
CoolPropDbl d4alphar_dDelta2_dTau2(void)
Return the term .
Definition: AbstractState.h:1566
Supercritical liquid (p > pc, T < Tc)
Definition: DataStructures.h:181
virtual CoolPropDbl calc_physical_hazard(void)
Using this backend, calculate the physical hazard.
Definition: AbstractState.h:408
SimpleState _critical
Two important points.
Definition: AbstractState.h:99
double T
temperature in K
Definition: AbstractState.h:36
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.
Definition: AbstractState.h:702
double rhomolar(void)
Return the molar density in mol/m^3.
Definition: AbstractState.h:1052
virtual CoolPropDbl calc_flame_hazard(void)
Using this backend, calculate the flame hazard.
Definition: AbstractState.h:400
CoolPropDbl first_partial_deriv(parameters Of, parameters Wrt, parameters Constant)
The first partial derivative in homogeneous phases.
Definition: AbstractState.h:1231
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!!!)
Definition: AbstractState.h:895
virtual CoolPropDbl calc_ODP(void)
Using this backend, calculate the ozone depletion potential (ODP)
Definition: AbstractState.h:396
std::vector< double > tau
The reciprocal reduced temperature ( )
Definition: AbstractState.h:26
CoolPropDbl dalpha0_dTau(void)
Return the term .
Definition: AbstractState.h:1465
CoolPropDbl d2alpha0_dDelta2(void)
Return the term .
Definition: AbstractState.h:1470
std::string name()
Return the name - backend dependent.
Definition: AbstractState.h:1017
virtual CoolPropDbl calc_Bvirial(void)
Using this backend, calculate the second virial coefficient.
Definition: AbstractState.h:431
virtual CoolPropDbl calc_dalpha0_dDelta(void)
Using this backend, calculate the ideal-gas Helmholtz energy term (dimensionless) ...
Definition: AbstractState.h:330
double tangent_plane_distance(const double T, const double p, const std::vector< double > &w, const double rhomolar_guess=-1)
Return the tangent plane distance for a given trial composition w.
Definition: AbstractState.h:1002
virtual CoolPropDbl calc_saturation_ancillary(parameters param, int Q, parameters given, double value)
Definition: AbstractState.h:561
virtual CoolPropDbl calc_d2alpha0_dTau2(void)
Using this backend, calculate the ideal-gas Helmholtz energy term (dimensionless) ...
Definition: AbstractState.h:346
CoolPropDbl d3alpha0_dDelta3(void)
Return the term .
Definition: AbstractState.h:1500
virtual CoolPropDbl calc_d2alphar_dTau2(void)
Using this backend, calculate the residual Helmholtz energy term (dimensionless) ...
Definition: AbstractState.h:283
virtual CoolPropDbl calc_reduced_density(void)
Using this backend, calculate the reduced density (rho/rhoc)
Definition: AbstractState.h:422
virtual CoolPropDbl calc_hmolar(void)
Using this backend, calculate the molar enthalpy in J/mol.
Definition: AbstractState.h:147
double hmass(void)
Return the mass enthalpy in J/kg.
Definition: AbstractState.h:1092
virtual CoolPropDbl calc_d3alphar_dDelta_dTau2(void)
Using this backend, calculate the residual Helmholtz energy term (dimensionless) ...
Definition: AbstractState.h:295
CachedElement _viscosity
Transport properties.
Definition: AbstractState.h:113
virtual void update_with_guesses(CoolProp::input_pairs input_pair, double Value1, double Value2, const GuessesStructure &guesses)
Update the state using two state variables and providing guess values Some or all of the guesses will...
Definition: AbstractState.h:837
virtual CoolPropDbl calc_cpmolar_idealgas(void)
Using this backend, calculate the ideal gas molar constant-pressure specific heat in J/mol/K...
Definition: AbstractState.h:171
double gibbsmass(void)
Return the Gibbs energy in J/kg.
Definition: AbstractState.h:1150
Subcritical liquid.
Definition: DataStructures.h:178
virtual CoolPropDbl calc_d3alphar_dTau3(void)
Using this backend, calculate the residual Helmholtz energy term (dimensionless) ...
Definition: AbstractState.h:299
virtual void calc_specify_phase(phases phase)
Using this backend, specify the phase to be used for all further calculations.
Definition: AbstractState.h:570
virtual CoolPropDbl calc_acentric_factor(void)
Using this backend, calculate the acentric factor.
Definition: AbstractState.h:223
Supercritical gas (p < pc, T > Tc)
Definition: DataStructures.h:180
virtual std::vector< std::string > calc_fluid_names(void)
Using this backend, get a vector of fluid names.
Definition: AbstractState.h:578
phases
These are constants for the phases of the fluid.
Definition: DataStructures.h:176
double rhomass(void)
Return the mass density in kg/m^3.
Definition: AbstractState.h:1056
double smass_excess(void)
Return the molar entropy in J/kg/K.
Definition: AbstractState.h:1112
phases phase(void)
Get the phase of the state.
Definition: AbstractState.h:938
virtual CoolPropDbl calc_compressibility_factor(void)
Using this backend, calculate the compressibility factor Z .
Definition: AbstractState.h:447
CoolPropDbl d3alpha0_dTau3(void)
Return the term .
Definition: AbstractState.h:1485
virtual void calc_excess_properties(void)
Using this backend, calculate and cache the excess properties.
Definition: AbstractState.h:257
virtual std::string calc_description(void)
Using this backend, get the description of the fluid.
Definition: AbstractState.h:456
virtual CoolPropDbl calc_dalpha0_dTau(void)
Using this backend, calculate the ideal-gas Helmholtz energy term (dimensionless) ...
Definition: AbstractState.h:334
double cp0mass(void)
Return the mass constant pressure specific heat for ideal gas part only in J/kg/K.
Definition: AbstractState.h:1136
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!!!)
Definition: AbstractState.h:879
virtual CoolPropDbl calc_fraction_min(void)
Get the minimum fraction (mole, mass, volume) for incompressible fluid.
Definition: AbstractState.h:602
virtual CoolPropDbl calc_d4alphar_dDelta2_dTau2(void)
Using this backend, calculate the residual Helmholtz energy term (dimensionless) ...
Definition: AbstractState.h:312
double hmass_excess(void)
Return the excess mass enthalpy in J/kg.
Definition: AbstractState.h:1098
long _fluid_type
Some administrative variables.
Definition: AbstractState.h:82
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 ...
Definition: AbstractState.h:1397
void conformal_state(const std::string &reference_fluid, CoolPropDbl &T, CoolPropDbl &rhomolar)
Find the conformal state needed for ECS.
Definition: AbstractState.h:1439
virtual CoolPropDbl calc_health_hazard(void)
Using this backend, calculate the health hazard.
Definition: AbstractState.h:404
virtual CoolPropDbl calc_Cvirial(void)
Using this backend, calculate the third virial coefficient.
Definition: AbstractState.h:435
Unknown phase.
Definition: DataStructures.h:185
std::vector< CriticalState > all_critical_points(void)
Return the vector of critical points, including points that are unstable or correspond to negative pr...
Definition: AbstractState.h:960
An abstract AbstractState generator class.
Definition: AbstractState.h:1589
double cpmass(void)
Return the mass constant pressure specific heat in J/kg/K.
Definition: AbstractState.h:1130
CoolPropDbl d2alphar_dDelta_dTau(void)
Return the term .
Definition: AbstractState.h:1526
virtual CoolPropDbl calc_Tmin(void)
Using this backend, calculate the minimum temperature in K.
Definition: AbstractState.h:375
double umass_excess(void)
Return the excess internal energy in J/kg.
Definition: AbstractState.h:1124
virtual CoolPropDbl calc_T_reducing(void)
Using this backend, get the reducing point temperature in K.
Definition: AbstractState.h:474
The mother of all state classes.
Definition: AbstractState.h:78
virtual CoolPropDbl calc_reciprocal_reduced_temperature(void)
Using this backend, calculate the reciprocal reduced temperature (Tc/T)
Definition: AbstractState.h:426
virtual CoolPropDbl calc_d2alpha0_dDelta_dTau(void)
Using this backend, calculate the ideal-gas Helmholtz energy term (dimensionless) ...
Definition: AbstractState.h:338
CoolPropDbl d4alphar_dDelta3_dTau(void)
Return the term .
Definition: AbstractState.h:1561
virtual CoolPropDbl calc_GWP500(void)
Using this backend, calculate the 500-year global warming potential (GWP)
Definition: AbstractState.h:392
virtual std::string fluid_param_string(const std::string &)
Return a string from the backend for the mixture/fluid - backend dependent - could be CAS #...
Definition: AbstractState.h:849
CoolPropDbl d3alpha0_dDelta_dTau2(void)
Return the term .
Definition: AbstractState.h:1490
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!!!)
Definition: AbstractState.h:866
virtual const double get_fluid_constant(std::size_t i, parameters param) const
Get a constant for one of the fluids forming this mixture.
Definition: AbstractState.h:860
virtual CoolPropDbl calc_d2alphar_dDelta_dTau(void)
Using this backend, calculate the residual Helmholtz energy term (dimensionless) ...
Definition: AbstractState.h:279
double helmholtzmass(void)
Return the Helmholtz energy in J/kg.
Definition: AbstractState.h:1162
std::vector< CoolPropDbl > mole_fractions_liquid(void)
Get the mole fractions of the equilibrium liquid phase.
Definition: AbstractState.h:806
void set_T(CoolPropDbl T)
Set the internal variable T without a flash call (expert use only!)
Definition: AbstractState.h:732
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...
Definition: AbstractState.h:654
CoolPropDbl d3alphar_dTau3(void)
Return the term .
Definition: AbstractState.h:1551
CachedElement _rho_spline
Smoothing values.
Definition: AbstractState.h:129
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.
Definition: AbstractState.h:1371
virtual CoolPropDbl calc_pmax(void)
Using this backend, calculate the maximum pressure in Pa.
Definition: AbstractState.h:379
This structure holds values obtained while tracing the spinodal curve (most often in the process of f...
Definition: AbstractState.h:23
A data structure to hold the data for a phase envelope.
Definition: PhaseEnvelope.h:36
virtual CoolPropDbl calc_chemical_potential(std::size_t i)
Using this backend, calculate the chemical potential in J/mol.
Definition: AbstractState.h:247
virtual CoolPropDbl calc_Ttriple(void)
Using this backend, get the triple point temperature in K.
Definition: AbstractState.h:461
backend_families
The structure is taken directly from the AbstractState class.
Definition: DataStructures.h:438
At the critical point.
Definition: DataStructures.h:182
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...
Definition: AbstractState.h:1208
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!!!)
Definition: AbstractState.h:874
double saturated_liquid_keyed_output(parameters key)
Get an output from the saturated liquid state by key.
Definition: AbstractState.h:1039
CoolPropDbl first_saturation_deriv(parameters Of1, parameters Wrt1)
The first partial derivative along the saturation curve.
Definition: AbstractState.h:1283
std::string description()
Return the description - backend dependent.
Definition: AbstractState.h:1021
virtual CoolPropDbl calc_viscosity(void)
Using this backend, calculate the viscosity in Pa-s.
Definition: AbstractState.h:207
void conductivity_contributions(CoolPropDbl &dilute, CoolPropDbl &initial_density, CoolPropDbl &residual, CoolPropDbl &critical)
Return the thermal conductivity contributions, each in W/m/K.
Definition: AbstractState.h:1424
void specify_phase(phases phase)
Specify the phase for all further calculations with this state class.
Definition: AbstractState.h:942
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.
Definition: AbstractState.h:685
CachedElement _alpha0
Cached low-level elements for in-place calculation of other properties.
Definition: AbstractState.h:132
input_pairs
These are input pairs that can be used for the update function (in each pair, input keys are sorted a...
Definition: DataStructures.h:274
virtual CoolPropDbl calc_d3alpha0_dTau3(void)
Using this backend, calculate the ideal-gas Helmholtz energy term (dimensionless) ...
Definition: AbstractState.h:362
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:
Definition: AbstractState.h:899
void build_spinodal()
Construct the spinodal curve for the mixture (or pure fluid)
Definition: AbstractState.h:965
virtual void calc_phase_envelope(const std::string &type)
Using this backend, construct the phase envelope, the variable type describes the type of phase envel...
Definition: AbstractState.h:498
virtual bool available_in_high_level(void)
A function that says whether the backend instance can be instantiated in the high-level interface In ...
Definition: AbstractState.h:844
virtual void calc_unspecify_phase(void)
Using this backend, unspecify the phase.
Definition: AbstractState.h:574
Subcritical gas.
Definition: DataStructures.h:183
Definition: Exceptions.h:45
CachedElement _hmolar_excess
Excess properties.
Definition: AbstractState.h:121
CoolPropDbl d2alpha0_dTau2(void)
Return the term .
Definition: AbstractState.h:1480
virtual CoolPropDbl calc_Tmax(void)
Using this backend, calculate the maximum temperature in K.
Definition: AbstractState.h:371
virtual CoolPropDbl calc_conductivity(void)
Using this backend, calculate the thermal conductivity in W/m/K.
Definition: AbstractState.h:211
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.
Definition: AbstractState.h:1447
CoolPropDbl d3alphar_dDelta3(void)
Return the term .
Definition: AbstractState.h:1536
double helmholtzmass_excess(void)
Return the excess Helmholtz energy in J/kg.
Definition: AbstractState.h:1168
CoolPropDbl d3alphar_dDelta2_dTau(void)
Return the term .
Definition: AbstractState.h:1541
virtual CoolPropDbl calc_isobaric_expansion_coefficient(void)
Using this backend, calculate the isobaric expansion coefficient in 1/K.
Definition: AbstractState.h:199
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".
Definition: AbstractState.h:1203
virtual CoolPropDbl calc_molar_mass(void)
Using this backend, calculate the molar mass in kg/mol.
Definition: AbstractState.h:219
double Q(void)
Return the vapor quality (mol/mol); Q = 0 for saturated liquid.
Definition: AbstractState.h:1064
CachedElement _rhoLmolar
Two-Phase variables.
Definition: AbstractState.h:140
virtual CoolPropDbl calc_dBvirial_dT(void)
Using this backend, calculate the derivative dB/dT.
Definition: AbstractState.h:439
virtual CoolPropDbl calc_umolar(void)
Using this backend, calculate the molar internal energy in J/mol.
Definition: AbstractState.h:163
virtual CoolPropDbl calc_T(void)
Using this backend, get the temperature.
Definition: AbstractState.h:640
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.
Definition: AbstractState.h:649
CoolPropDbl d4alphar_dDelta_dTau3(void)
Return the term .
Definition: AbstractState.h:1571
virtual CoolPropDbl calc_d3alpha0_dDelta3(void)
Using this backend, calculate the ideal-gas Helmholtz energy term (dimensionless) ...
Definition: AbstractState.h:350
virtual CoolPropDbl calc_helmholtzmolar(void)
Using this backend, calculate the molar Helmholtz energy in J/mol.
Definition: AbstractState.h:187
CachedElement _molar_mass
Molar mass [mol/kg].
Definition: AbstractState.h:102
double saturated_vapor_keyed_output(parameters key)
Get an output from the saturated vapor state by key.
Definition: AbstractState.h:1043
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) ...
Definition: AbstractState.h:820
virtual CoolPropDbl calc_gas_constant(void)
Using this backend, calculate the universal gas constant in J/mol/K.
Definition: AbstractState.h:231
virtual void update_states(void)
Update the states after having changed the reference state for enthalpy and entropy.
Definition: AbstractState.h:549
virtual CoolPropDbl calc_d2alpha0_dDelta2(void)
Using this backend, calculate the ideal-gas Helmholtz energy term (dimensionless) ...
Definition: AbstractState.h:342
void criticality_contour_values(double &L1star, double &M1star)
Calculate the criticality contour values and .
Definition: AbstractState.h:975
virtual double get_fluid_parameter_double(const size_t i, const std::string ¶meter)
Double fluid parameter (currently the volume translation parameter for cubic)
Definition: AbstractState.h:907
virtual CoolPropDbl calc_isentropic_expansion_coefficient(void)
Using this backend, calculate the isentropic expansion coefficient .
Definition: AbstractState.h:203
CachedElement _gas_constant
Universal gas constant [J/mol/K].
Definition: AbstractState.h:105
virtual CoolPropDbl calc_pressure(void)
Using this backend, calculate the pressure in Pa.
Definition: AbstractState.h:227
double smass(void)
Return the molar entropy in J/kg/K.
Definition: AbstractState.h:1106
virtual CoolPropDbl calc_rhomolar(void)
Using this backend, get the molar density in mol/m^3.
Definition: AbstractState.h:644
CoolPropDbl d4alphar_dTau4(void)
Return the term .
Definition: AbstractState.h:1576
CoolPropDbl dalpha0_dDelta(void)
Return the term .
Definition: AbstractState.h:1460
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 densi...
Definition: AbstractState.h:782
double Prandtl(void)
Return the Prandtl number (dimensionless)
Definition: AbstractState.h:1430
virtual CoolPropDbl calc_alpha0(void)
Using this backend, calculate the ideal-gas Helmholtz energy term (dimensionless) ...
Definition: AbstractState.h:326
double volumemass_excess(void)
Return the excess volume in m^3/kg.
Definition: AbstractState.h:1174
virtual CoolPropDbl calc_rhomass_critical(void)
Using this backend, get the critical point mass density in kg/m^3 - Added for IF97Backend which is ma...
Definition: AbstractState.h:490
virtual CoolPropDbl calc_GWP20(void)
Using this backend, calculate the 20-year global warming potential (GWP)
Definition: AbstractState.h:384
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!!!)
Definition: AbstractState.h:883
virtual CoolPropDbl calc_d2alphar_dDelta2(void)
Using this backend, calculate the residual Helmholtz energy term (dimensionless) ...
Definition: AbstractState.h:275
std::vector< double > M1
The determinant of the scaled matrix for the second criticality condition.
Definition: AbstractState.h:26
virtual std::vector< CoolPropDbl > calc_fugacity_coefficients()
Using this backend, calculate the fugacity in Pa.
Definition: AbstractState.h:239
virtual CoolPropDbl calc_cpmolar(void)
Using this backend, calculate the molar constant-pressure specific heat in J/mol/K.
Definition: AbstractState.h:167
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!!!)
Definition: AbstractState.h:870
CoolPropDbl alpha0(void)
Return the term .
Definition: AbstractState.h:1455
virtual CoolPropDbl calc_gibbsmolar_residual(void)
Using this backend, calculate the residual molar Gibbs function in J/mol.
Definition: AbstractState.h:183
virtual std::string calc_name(void)
Using this backend, get the name of the fluid.
Definition: AbstractState.h:452
void unspecify_phase(void)
Unspecify the phase and go back to calculating it based on the inputs.
Definition: AbstractState.h:946
virtual CoolPropDbl calc_d4alphar_dDelta3_dTau(void)
Using this backend, calculate the residual Helmholtz energy term (dimensionless) ...
Definition: AbstractState.h:308
double p(void)
Return the pressure in Pa.
Definition: AbstractState.h:1060
virtual CoolPropDbl calc_d3alphar_dDelta3(void)
Using this backend, calculate the residual Helmholtz energy term (dimensionless) ...
Definition: AbstractState.h:287
virtual CoolPropDbl calc_d4alphar_dDelta_dTau3(void)
Using this backend, calculate the residual Helmholtz energy term (dimensionless) ...
Definition: AbstractState.h:316
double umass(void)
Return the mass internal energy in J/kg.
Definition: AbstractState.h:1118
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!!!)
Definition: AbstractState.h:887
virtual CoolPropDbl calc_speed_sound(void)
Using this backend, calculate the speed of sound in m/s.
Definition: AbstractState.h:191
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) ...
Definition: AbstractState.h:810
CoolPropDbl d3alpha0_dDelta2_dTau(void)
Return the term .
Definition: AbstractState.h:1495
Twophase.
Definition: DataStructures.h:184
phases _phase
The key for the phase from CoolProp::phases enum.
Definition: AbstractState.h:83
virtual CoolPropDbl calc_T_critical(void)
Using this backend, get the critical point temperature in K.
Definition: AbstractState.h:470
CachedElement _hmolar_residual
Residual properties.
Definition: AbstractState.h:118
virtual CoolPropDbl calc_p_reducing(void)
Using this backend, get the reducing point pressure in Pa.
Definition: AbstractState.h:482
virtual const CoolProp::SimpleState & get_reducing_state()
Get the state that is used in the equation of state or mixture model to reduce the state...
Definition: AbstractState.h:919
double _rhomolar
Bulk values.
Definition: AbstractState.h:108
virtual CoolPropDbl calc_alphar(void)
Using this backend, calculate the residual Helmholtz energy term (dimensionless) ...
Definition: AbstractState.h:263
CoolPropDbl second_partial_deriv(parameters Of1, parameters Wrt1, parameters Constant1, parameters Wrt2, parameters Constant2)
The second partial derivative in homogeneous phases.
Definition: AbstractState.h:1258
std::vector< double > y
molar composition of the vapor phase
Definition: AbstractState.h:43
const CoolProp::PhaseEnvelopeData & get_phase_envelope_data()
After having calculated the phase envelope, return the phase envelope data.
Definition: AbstractState.h:1388
virtual CoolPropDbl calc_d4alphar_dDelta4(void)
Using this backend, calculate the residual Helmholtz energy term (dimensionless) ...
Definition: AbstractState.h:304
This simple class holds the values for guesses for use in some solvers that have the ability to use g...
Definition: AbstractState.h:33
virtual CoolPropDbl calc_rhomolar_reducing(void)
Using this backend, get the reducing point molar density in mol/m^3.
Definition: AbstractState.h:494
const CoolProp::SimpleState & get_state(const std::string &state)
Get a desired state point - backend dependent.
Definition: AbstractState.h:924
virtual CoolPropDbl calc_dCvirial_dT(void)
Using this backend, calculate the derivative dC/dT.
Definition: AbstractState.h:443
CoolPropDbl alphar(void)
Return the term .
Definition: AbstractState.h:1506
CoolPropDbl d2alphar_dDelta2(void)
Return the term .
Definition: AbstractState.h:1521
virtual CoolPropDbl calc_p_triple(void)
Using this backend, get the triple point pressure in Pa.
Definition: AbstractState.h:465
double gibbsmass_excess(void)
Return the excess Gibbs energy in J/kg.
Definition: AbstractState.h:1156
virtual CoolPropDbl calc_GWP100(void)
Using this backend, calculate the 100-year global warming potential (GWP)
Definition: AbstractState.h:388
Definition: AbstractState.h:1602
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)
Definition: AbstractState.h:903
virtual CoolPropDbl calc_fugacity(std::size_t i)
Using this backend, calculate the fugacity in Pa.
Definition: AbstractState.h:243
CachedElement _rhoLanc
Ancillary values.
Definition: AbstractState.h:124
virtual CoolPropDbl calc_p_critical(void)
Using this backend, get the critical point pressure in Pa.
Definition: AbstractState.h:478
CoolPropDbl d2alphar_dTau2(void)
Return the term .
Definition: AbstractState.h:1531
virtual CoolPropDbl calc_PIP(void)
Using this backend, calculate the phase identification parameter (PIP)
Definition: AbstractState.h:251
CoolPropDbl d3alphar_dDelta_dTau2(void)
Return the term .
Definition: AbstractState.h:1546
virtual CoolPropDbl calc_cvmolar(void)
Using this backend, calculate the molar constant-volume specific heat in J/mol/K. ...
Definition: AbstractState.h:175
virtual CoolPropDbl calc_d3alphar_dDelta2_dTau(void)
Using this backend, calculate the residual Helmholtz energy term (dimensionless) ...
Definition: AbstractState.h:291
virtual const CoolProp::SimpleState & calc_state(const std::string &state)
Using this backend, calculate a phase given by the state string.
Definition: AbstractState.h:583
virtual CoolPropDbl calc_fugacity_coefficient(std::size_t i)
Using this backend, calculate the fugacity coefficient (dimensionless)
Definition: AbstractState.h:235
virtual const std::vector< CoolPropDbl > get_mass_fractions(void)
Get the mass fractions of the fluid.
Definition: AbstractState.h:828
virtual CoolPropDbl calc_d4alphar_dTau4(void)
Using this backend, calculate the residual Helmholtz energy term (dimensionless) ...
Definition: AbstractState.h:320
Supercritical (p > pc, T > Tc)
Definition: DataStructures.h:179
This file contains flash routines in which the state is unknown, and a solver of some kind must be us...
Definition: AbstractState.h:19
CoolPropDbl d2alpha0_dDelta_dTau(void)
Return the term .
Definition: AbstractState.h:1475
virtual CoolPropDbl calc_fraction_max(void)
Get the maximum fraction (mole, mass, volume) for incompressible fluid.
Definition: AbstractState.h:606
virtual CoolPropDbl calc_smolar_residual(void)
Using this backend, calculate the residual molar entropy in J/mol/K.
Definition: AbstractState.h:159
void viscosity_contributions(CoolPropDbl &dilute, CoolPropDbl &initial_density, CoolPropDbl &residual, CoolPropDbl &critical)
Return the viscosity contributions, each in Pa-s.
Definition: AbstractState.h:1418
CoolPropDbl dalphar_dTau(void)
Return the term .
Definition: AbstractState.h:1516
virtual void set_reference_stateS(const std::string &reference_state)
Set the reference state based on a string representation.
Definition: AbstractState.h:772
virtual CoolPropDbl calc_d3alpha0_dDelta2_dTau(void)
Using this backend, calculate the ideal-gas Helmholtz energy term (dimensionless) ...
Definition: AbstractState.h:354
CoolPropDbl d4alphar_dDelta4(void)
Return the term .
Definition: AbstractState.h:1556
parameters
Define some constants that will be used throughout These are constants for the input and output para...
Definition: DataStructures.h:64
virtual CoolPropDbl calc_gibbsmolar(void)
Using this backend, calculate the molar Gibbs function in J/mol.
Definition: AbstractState.h:179
virtual CoolPropDbl calc_d3alpha0_dDelta_dTau2(void)
Using this backend, calculate the ideal-gas Helmholtz energy term (dimensionless) ...
Definition: AbstractState.h:358
virtual phases calc_phase(void)
Using this backend, calculate the phase.
Definition: AbstractState.h:566
void register_backend(const backend_families &bf, shared_ptr< AbstractStateGenerator > gen)
Register a backend in the backend library (statically defined in AbstractState.cpp and not publicly a...
Definition: AbstractState.cpp:56
Definition: CachedElement.h:32
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.
Definition: AbstractState.h:1327
virtual CoolPropDbl calc_dalphar_dTau(void)
Using this backend, calculate the residual Helmholtz energy term (dimensionless) ...
Definition: AbstractState.h:271
virtual CoolPropDbl calc_smolar(void)
Using this backend, calculate the molar entropy in J/mol/K.
Definition: AbstractState.h:155
virtual CoolPropDbl calc_rhomolar_critical(void)
Using this backend, get the critical point molar density in mol/m^3.
Definition: AbstractState.h:486
virtual CoolPropDbl calc_isothermal_compressibility(void)
Using this backend, calculate the isothermal compressibility in 1/Pa.
Definition: AbstractState.h:195
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.
Definition: AbstractState.h:1347
CoolPropDbl second_saturation_deriv(parameters Of1, parameters Wrt1, parameters Wrt2)
The second partial derivative along the saturation curve.
Definition: AbstractState.h:1304
virtual CoolPropDbl calc_surface_tension(void)
Using this backend, calculate the surface tension in N/m.
Definition: AbstractState.h:215
virtual CoolPropDbl calc_dipole_moment(void)
Using this backend, calculate the dipole moment in C-m (1 D = 3.33564e-30 C-m)
Definition: AbstractState.h:412
void ideal_curve(const std::string &type, std::vector< double > &T, std::vector< double > &p)
Calculate an ideal curve for a pure fluid.
Definition: AbstractState.h:1219
double dipole_moment()
Return the dipole moment in C-m (1 D = 3.33564e-30 C-m)
Definition: AbstractState.h:1026
Definition: DataStructures.h:16
SpinodalData get_spinodal_data()
Get the data from the spinodal curve constructed in the call to build_spinodal()
Definition: AbstractState.h:970
std::vector< CoolPropDbl > mole_fractions_vapor(void)
Get the mole fractions of the equilibrium vapor phase.
Definition: AbstractState.h:816