CoolProp 6.8.1dev
An open-source fluid property and humid air property database
AbstractState.h
Go to the documentation of this file.
1/*
2 * AbstractState.h
3 *
4 * Created on: 21 Dec 2013
5 * Author: jowr
6 */
7
8#ifndef ABSTRACTSTATE_H_
9#define ABSTRACTSTATE_H_
10
11#include "CachedElement.h"
12#include "Exceptions.h"
13#include "DataStructures.h"
14#include "PhaseEnvelope.h"
16
17#include <numeric>
18
19namespace CoolProp {
20
24{
25 public:
26 std::vector<double> tau,
29};
30
34{
35 public:
36 double T,
37 p,
43 std::vector<double> x,
44 y;
46 clear();
47 };
48 void clear() {
49 T = _HUGE;
50 p = _HUGE;
51 rhomolar = _HUGE;
52 hmolar = _HUGE;
53 smolar = _HUGE;
54 rhomolar_liq = _HUGE;
55 rhomolar_vap = _HUGE;
56 x.clear(), y.clear();
57 }
58};
59
61
79{
80 protected:
85
87 return (this->_phase == iphase_supercritical || this->_phase == iphase_supercritical_liquid || this->_phase == iphase_supercritical_gas);
88 }
89
90 bool isHomogeneousPhase(void) {
91 return (this->_phase == iphase_liquid || this->_phase == iphase_gas || isSupercriticalPhase() || this->_phase == iphase_critical_point);
92 }
93
94 bool isTwoPhase(void) {
95 return (this->_phase == iphase_twophase);
96 }
97
100
103
106
109
111 double _rhomolar, _T, _p, _Q;
112
114
117
119
122
125
128
130
133
139
141
144
145 // ----------------------------------------
146 // Property accessors to be optionally implemented by the backend
147 // for properties that are not always calculated
148 // ----------------------------------------
150 virtual CoolPropDbl calc_hmolar(void) {
151 throw NotImplementedError("calc_hmolar is not implemented for this backend");
152 };
155 throw NotImplementedError("calc_hmolar_residual is not implemented for this backend");
156 };
158 virtual CoolPropDbl calc_smolar(void) {
159 throw NotImplementedError("calc_smolar is not implemented for this backend");
160 };
163 throw NotImplementedError("calc_smolar_residual is not implemented for this backend");
164 };
166 virtual CoolPropDbl calc_neff(void) {
167 throw NotImplementedError("calc_neff is not implemented for this backend");
168 };
170 virtual CoolPropDbl calc_umolar(void) {
171 throw NotImplementedError("calc_umolar is not implemented for this backend");
172 };
174 virtual CoolPropDbl calc_cpmolar(void) {
175 throw NotImplementedError("calc_cpmolar is not implemented for this backend");
176 };
179 throw NotImplementedError("calc_cpmolar_idealgas is not implemented for this backend");
180 };
182 virtual CoolPropDbl calc_cvmolar(void) {
183 throw NotImplementedError("calc_cvmolar is not implemented for this backend");
184 };
187 throw NotImplementedError("calc_gibbsmolar is not implemented for this backend");
188 };
191 throw NotImplementedError("calc_gibbsmolar_residual is not implemented for this backend");
192 };
195 throw NotImplementedError("calc_helmholtzmolar is not implemented for this backend");
196 };
199 throw NotImplementedError("calc_speed_sound is not implemented for this backend");
200 };
203 throw NotImplementedError("calc_isothermal_compressibility is not implemented for this backend");
204 };
207 throw NotImplementedError("calc_isobaric_expansion_coefficient is not implemented for this backend");
208 };
211 throw NotImplementedError("calc_isentropic_expansion_coefficient is not implemented for this backend");
212 };
215 throw NotImplementedError("calc_viscosity is not implemented for this backend");
216 };
219 throw NotImplementedError("calc_conductivity is not implemented for this backend");
220 };
223 throw NotImplementedError("calc_surface_tension is not implemented for this backend");
224 };
227 throw NotImplementedError("calc_molar_mass is not implemented for this backend");
228 };
231 throw NotImplementedError("calc_acentric_factor is not implemented for this backend");
232 };
235 throw NotImplementedError("calc_pressure is not implemented for this backend");
236 };
239 throw NotImplementedError("calc_gas_constant is not implemented for this backend");
240 };
242 virtual CoolPropDbl calc_fugacity_coefficient(std::size_t i) {
243 throw NotImplementedError("calc_fugacity_coefficient is not implemented for this backend");
244 };
246 virtual std::vector<CoolPropDbl> calc_fugacity_coefficients() {
247 throw NotImplementedError("calc_fugacity_coefficients is not implemented for this backend");
248 };
250 virtual CoolPropDbl calc_fugacity(std::size_t i) {
251 throw NotImplementedError("calc_fugacity is not implemented for this backend");
252 };
254 virtual CoolPropDbl calc_chemical_potential(std::size_t i) {
255 throw NotImplementedError("calc_chemical_potential is not implemented for this backend");
256 };
258 virtual CoolPropDbl calc_PIP(void) {
259 throw NotImplementedError("calc_PIP is not implemented for this backend");
260 };
261
262 // Excess properties
264 virtual void calc_excess_properties(void) {
265 throw NotImplementedError("calc_excess_properties is not implemented for this backend");
266 };
267
268 // Derivatives of residual helmholtz energy
270 virtual CoolPropDbl calc_alphar(void) {
271 throw NotImplementedError("calc_alphar is not implemented for this backend");
272 };
275 throw NotImplementedError("calc_dalphar_dDelta is not implemented for this backend");
276 };
279 throw NotImplementedError("calc_dalphar_dTau is not implemented for this backend");
280 };
283 throw NotImplementedError("calc_d2alphar_dDelta2 is not implemented for this backend");
284 };
287 throw NotImplementedError("calc_d2alphar_dDelta_dTau is not implemented for this backend");
288 };
291 throw NotImplementedError("calc_d2alphar_dTau2 is not implemented for this backend");
292 };
295 throw NotImplementedError("calc_d3alphar_dDelta3 is not implemented for this backend");
296 };
299 throw NotImplementedError("calc_d3alphar_dDelta2_dTau is not implemented for this backend");
300 };
303 throw NotImplementedError("calc_d3alphar_dDelta_dTau2 is not implemented for this backend");
304 };
307 throw NotImplementedError("calc_d3alphar_dTau3 is not implemented for this backend");
308 };
309
312 throw NotImplementedError("calc_d4alphar_dDelta4 is not implemented for this backend");
313 };
316 throw NotImplementedError("calc_d4alphar_dDelta3_dTau is not implemented for this backend");
317 };
320 throw NotImplementedError("calc_d4alphar_dDelta2_dTau2 is not implemented for this backend");
321 };
324 throw NotImplementedError("calc_d4alphar_dDelta_dTau3 is not implemented for this backend");
325 };
328 throw NotImplementedError("calc_d4alphar_dTau4 is not implemented for this backend");
329 };
330
331 // Derivatives of ideal-gas helmholtz energy
333 virtual CoolPropDbl calc_alpha0(void) {
334 throw NotImplementedError("calc_alpha0 is not implemented for this backend");
335 };
338 throw NotImplementedError("calc_dalpha0_dDelta is not implemented for this backend");
339 };
342 throw NotImplementedError("calc_dalpha0_dTau is not implemented for this backend");
343 };
346 throw NotImplementedError("calc_d2alpha0_dDelta_dTau is not implemented for this backend");
347 };
350 throw NotImplementedError("calc_d2alpha0_dDelta2 is not implemented for this backend");
351 };
354 throw NotImplementedError("calc_d2alpha0_dTau2 is not implemented for this backend");
355 };
358 throw NotImplementedError("calc_d3alpha0_dDelta3 is not implemented for this backend");
359 };
362 throw NotImplementedError("calc_d3alpha0_dDelta2_dTau is not implemented for this backend");
363 };
366 throw NotImplementedError("calc_d3alpha0_dDelta_dTau2 is not implemented for this backend");
367 };
370 throw NotImplementedError("calc_d3alpha0_dTau3 is not implemented for this backend");
371 };
372
373 virtual void calc_reducing_state(void) {
374 throw NotImplementedError("calc_reducing_state is not implemented for this backend");
375 };
376
378 virtual CoolPropDbl calc_Tmax(void) {
379 throw NotImplementedError("calc_Tmax is not implemented for this backend");
380 };
382 virtual CoolPropDbl calc_Tmin(void) {
383 throw NotImplementedError("calc_Tmin is not implemented for this backend");
384 };
386 virtual CoolPropDbl calc_pmax(void) {
387 throw NotImplementedError("calc_pmax is not implemented for this backend");
388 };
389
391 virtual CoolPropDbl calc_GWP20(void) {
392 throw NotImplementedError("calc_GWP20 is not implemented for this backend");
393 };
395 virtual CoolPropDbl calc_GWP100(void) {
396 throw NotImplementedError("calc_GWP100 is not implemented for this backend");
397 };
399 virtual CoolPropDbl calc_GWP500(void) {
400 throw NotImplementedError("calc_GWP500 is not implemented for this backend");
401 };
403 virtual CoolPropDbl calc_ODP(void) {
404 throw NotImplementedError("calc_ODP is not implemented for this backend");
405 };
408 throw NotImplementedError("calc_flame_hazard is not implemented for this backend");
409 };
412 throw NotImplementedError("calc_health_hazard is not implemented for this backend");
413 };
416 throw NotImplementedError("calc_physical_hazard is not implemented for this backend");
417 };
420 throw NotImplementedError("calc_dipole_moment is not implemented for this backend");
421 };
422
426 virtual CoolPropDbl calc_second_partial_deriv(parameters Of1, parameters Wrt1, parameters Constant1, parameters Wrt2, parameters Constant2);
427
430 throw NotImplementedError("calc_reduced_density is not implemented for this backend");
431 };
434 throw NotImplementedError("calc_reciprocal_reduced_temperature is not implemented for this backend");
435 };
436
438 virtual CoolPropDbl calc_Bvirial(void) {
439 throw NotImplementedError("calc_Bvirial is not implemented for this backend");
440 };
442 virtual CoolPropDbl calc_Cvirial(void) {
443 throw NotImplementedError("calc_Cvirial is not implemented for this backend");
444 };
447 throw NotImplementedError("calc_dBvirial_dT is not implemented for this backend");
448 };
451 throw NotImplementedError("calc_dCvirial_dT is not implemented for this backend");
452 };
455 throw NotImplementedError("calc_compressibility_factor is not implemented for this backend");
456 };
457
459 virtual std::string calc_name(void) {
460 throw NotImplementedError("calc_name is not implemented for this backend");
461 };
463 virtual std::string calc_description(void) {
464 throw NotImplementedError("calc_description is not implemented for this backend");
465 };
466
468 virtual CoolPropDbl calc_Ttriple(void) {
469 throw NotImplementedError("calc_Ttriple is not implemented for this backend");
470 };
473 throw NotImplementedError("calc_p_triple is not implemented for this backend");
474 };
475
478 throw NotImplementedError("calc_T_critical is not implemented for this backend");
479 };
482 throw NotImplementedError("calc_T_reducing is not implemented for this backend");
483 };
486 throw NotImplementedError("calc_p_critical is not implemented for this backend");
487 };
490 throw NotImplementedError("calc_p_reducing is not implemented for this backend");
491 };
494 throw NotImplementedError("calc_rhomolar_critical is not implemented for this backend");
495 };
498 throw NotImplementedError("calc_rhomass_critical is not implemented for this backend");
499 };
502 throw NotImplementedError("calc_rhomolar_reducing is not implemented for this backend");
503 };
505 virtual void calc_phase_envelope(const std::string& type) {
506 throw NotImplementedError("calc_phase_envelope is not implemented for this backend");
507 };
509 virtual CoolPropDbl calc_rhomass(void) {
510 return rhomolar() * molar_mass();
511 }
512 virtual CoolPropDbl calc_hmass(void) {
513 return hmolar() / molar_mass();
514 }
516 return hmolar_excess() / molar_mass();
517 }
518 virtual CoolPropDbl calc_smass(void) {
519 return smolar() / molar_mass();
520 }
522 return smolar_excess() / molar_mass();
523 }
524 virtual CoolPropDbl calc_cpmass(void) {
525 return cpmolar() / molar_mass();
526 }
527 virtual CoolPropDbl calc_cp0mass(void) {
528 return cp0molar() / molar_mass();
529 }
530 virtual CoolPropDbl calc_cvmass(void) {
531 return cvmolar() / molar_mass();
532 }
533 virtual CoolPropDbl calc_umass(void) {
534 return umolar() / molar_mass();
535 }
537 return umolar_excess() / molar_mass();
538 }
540 return gibbsmolar() / molar_mass();
541 }
543 return gibbsmolar_excess() / molar_mass();
544 }
546 return helmholtzmolar() / molar_mass();
547 }
550 }
552 return volumemolar_excess() / molar_mass();
553 }
554
556 virtual void update_states(void) {
557 throw NotImplementedError("This backend does not implement update_states function");
558 };
559
560 virtual CoolPropDbl calc_melting_line(int param, int given, CoolPropDbl value) {
561 throw NotImplementedError("This backend does not implement calc_melting_line function");
562 };
563
568 virtual CoolPropDbl calc_saturation_ancillary(parameters param, int Q, parameters given, double value) {
569 throw NotImplementedError("This backend does not implement calc_saturation_ancillary");
570 };
571
573 virtual phases calc_phase(void) {
574 throw NotImplementedError("This backend does not implement calc_phase function");
575 };
578 throw NotImplementedError("This backend does not implement calc_specify_phase function");
579 };
581 virtual void calc_unspecify_phase(void) {
582 throw NotImplementedError("This backend does not implement calc_unspecify_phase function");
583 };
585 virtual std::vector<std::string> calc_fluid_names(void) {
586 throw NotImplementedError("This backend does not implement calc_fluid_names function");
587 };
590 virtual const CoolProp::SimpleState& calc_state(const std::string& state) {
591 throw NotImplementedError("calc_state is not implemented for this backend");
592 };
593
595 throw NotImplementedError("calc_phase_envelope_data is not implemented for this backend");
596 };
597
598 virtual std::vector<CoolPropDbl> calc_mole_fractions_liquid(void) {
599 throw NotImplementedError("calc_mole_fractions_liquid is not implemented for this backend");
600 };
601 virtual std::vector<CoolPropDbl> calc_mole_fractions_vapor(void) {
602 throw NotImplementedError("calc_mole_fractions_vapor is not implemented for this backend");
603 };
604 virtual const std::vector<CoolPropDbl> calc_mass_fractions(void) {
605 throw NotImplementedError("calc_mass_fractions is not implemented for this backend");
606 };
607
610 throw NotImplementedError("calc_fraction_min is not implemented for this backend");
611 };
614 throw NotImplementedError("calc_fraction_max is not implemented for this backend");
615 };
617 throw NotImplementedError("calc_T_freeze is not implemented for this backend");
618 };
619
621 throw NotImplementedError("calc_first_saturation_deriv is not implemented for this backend");
622 };
624 throw NotImplementedError("calc_second_saturation_deriv is not implemented for this backend");
625 };
627 throw NotImplementedError("calc_first_two_phase_deriv is not implemented for this backend");
628 };
630 throw NotImplementedError("calc_second_two_phase_deriv is not implemented for this backend");
631 };
633 throw NotImplementedError("calc_first_two_phase_deriv_splined is not implemented for this backend");
634 };
635
637 throw NotImplementedError("calc_saturated_liquid_keyed_output is not implemented for this backend");
638 };
640 throw NotImplementedError("calc_saturated_vapor_keyed_output is not implemented for this backend");
641 };
642 virtual void calc_ideal_curve(const std::string& type, std::vector<double>& T, std::vector<double>& p) {
643 throw NotImplementedError("calc_ideal_curve is not implemented for this backend");
644 };
645
647 virtual CoolPropDbl calc_T(void) {
648 return _T;
649 }
652 return _rhomolar;
653 }
654
656 virtual double calc_tangent_plane_distance(const double T, const double p, const std::vector<double>& w, const double rhomolar_guess) {
657 throw NotImplementedError("calc_tangent_plane_distance is not implemented for this backend");
658 };
659
661 virtual void calc_true_critical_point(double& T, double& rho) {
662 throw NotImplementedError("calc_true_critical_point is not implemented for this backend");
663 };
664
665 virtual void calc_conformal_state(const std::string& reference_fluid, CoolPropDbl& T, CoolPropDbl& rhomolar) {
666 throw NotImplementedError("calc_conformal_state is not implemented for this backend");
667 };
668
669 virtual void calc_viscosity_contributions(CoolPropDbl& dilute, CoolPropDbl& initial_density, CoolPropDbl& residual, CoolPropDbl& critical) {
670 throw NotImplementedError("calc_viscosity_contributions is not implemented for this backend");
671 };
672 virtual void calc_conductivity_contributions(CoolPropDbl& dilute, CoolPropDbl& initial_density, CoolPropDbl& residual, CoolPropDbl& critical) {
673 throw NotImplementedError("calc_conductivity_contributions is not implemented for this backend");
674 };
675 virtual std::vector<CriticalState> calc_all_critical_points(void) {
676 throw NotImplementedError("calc_all_critical_points is not implemented for this backend");
677 };
678 virtual void calc_build_spinodal() {
679 throw NotImplementedError("calc_build_spinodal is not implemented for this backend");
680 };
682 throw NotImplementedError("calc_get_spinodal_data is not implemented for this backend");
683 };
684 virtual void calc_criticality_contour_values(double& L1star, double& M1star) {
685 throw NotImplementedError("calc_criticality_contour_values is not implemented for this backend");
686 };
687
689 virtual void mass_to_molar_inputs(CoolProp::input_pairs& input_pair, CoolPropDbl& value1, CoolPropDbl& value2);
690
692 virtual void calc_change_EOS(const std::size_t i, const std::string& EOS_name) {
693 throw NotImplementedError("calc_change_EOS is not implemented for this backend");
694 };
695
696 public:
698 clear();
699 }
700 virtual ~AbstractState(){};
701
703
709 static AbstractState* factory(const std::string& backend, const std::string& fluid_names) {
710 return factory(backend, strsplit(fluid_names, '&'));
711 };
712
736 static AbstractState* factory(const std::string& backend, const std::vector<std::string>& fluid_names);
737
740 _T = T;
741 }
742
747 virtual std::string backend_name(void) = 0;
748
749 // The derived classes must implement this function to define whether they use mole fractions (true) or mass fractions (false)
750 virtual bool using_mole_fractions(void) = 0;
751 virtual bool using_mass_fractions(void) = 0;
752 virtual bool using_volu_fractions(void) = 0;
753
754 virtual void set_mole_fractions(const std::vector<CoolPropDbl>& mole_fractions) = 0;
755 virtual void set_mass_fractions(const std::vector<CoolPropDbl>& mass_fractions) = 0;
756 virtual void set_volu_fractions(const std::vector<CoolPropDbl>& mass_fractions) {
757 throw NotImplementedError("Volume composition has not been implemented.");
758 }
759
779 virtual void set_reference_stateS(const std::string& reference_state) {
781 "Setting reference state has not been implemented for this backend. Try using CoolProp::set_reference_stateD instead.");
782 }
783
789 virtual void set_reference_stateD(double T, double rhomolar, double hmolar0, double smolar0) {
791 "Setting reference state has not been implemented for this backend. Try using CoolProp::set_reference_stateD instead.");
792 }
793
794#ifndef COOLPROPDBL_MAPS_TO_DOUBLE
795 void set_mole_fractions(const std::vector<double>& mole_fractions) {
796 set_mole_fractions(std::vector<CoolPropDbl>(mole_fractions.begin(), mole_fractions.end()));
797 };
798 void set_mass_fractions(const std::vector<double>& mass_fractions) {
799 set_mass_fractions(std::vector<CoolPropDbl>(mass_fractions.begin(), mass_fractions.end()));
800 };
801 void set_volu_fractions(const std::vector<double>& volu_fractions) {
802 set_volu_fractions(std::vector<CoolPropDbl>(volu_fractions.begin(), volu_fractions.end()));
803 };
804#endif
805
806#ifdef EMSCRIPTEN
807 void set_mole_fractions_double(const std::vector<double>& mole_fractions) {
808 set_mole_fractions(std::vector<CoolPropDbl>(mole_fractions.begin(), mole_fractions.end()));
809 };
810#endif
811
813 std::vector<CoolPropDbl> mole_fractions_liquid(void) {
815 };
817 std::vector<double> mole_fractions_liquid_double(void) {
818 std::vector<CoolPropDbl> x = calc_mole_fractions_liquid();
819 return std::vector<double>(x.begin(), x.end());
820 };
821
823 std::vector<CoolPropDbl> mole_fractions_vapor(void) {
825 };
827 std::vector<double> mole_fractions_vapor_double(void) {
828 std::vector<CoolPropDbl> y = calc_mole_fractions_vapor();
829 return std::vector<double>(y.begin(), y.end());
830 };
831
833 virtual const std::vector<CoolPropDbl>& get_mole_fractions(void) = 0;
835 virtual const std::vector<CoolPropDbl> get_mass_fractions(void) {
836 return this->calc_mass_fractions();
837 };
838
840 virtual void update(CoolProp::input_pairs input_pair, double Value1, double Value2) = 0;
841
843 virtual void update_QT_pure_superanc(double Q, double T){
844 throw NotImplementedError("update_QT_pure_superanc is not implemented for this backend");
845 };
846
849 virtual void update_with_guesses(CoolProp::input_pairs input_pair, double Value1, double Value2, const GuessesStructure& guesses) {
850 throw NotImplementedError("update_with_guesses is not implemented for this backend");
851 };
852
856 virtual bool available_in_high_level(void) {
857 return true;
858 }
859
861 virtual std::string fluid_param_string(const std::string&) {
862 throw NotImplementedError("fluid_param_string has not been implemented for this backend");
863 }
864
866 std::vector<std::string> fluid_names(void);
867
872 virtual const double get_fluid_constant(std::size_t i, parameters param) const {
873 throw NotImplementedError("get_fluid_constant is not implemented for this backend");
874 };
875
877 virtual void set_binary_interaction_double(const std::string& CAS1, const std::string& CAS2, const std::string& parameter, const double value) {
878 throw NotImplementedError("set_binary_interaction_double is not implemented for this backend");
879 };
881 virtual void set_binary_interaction_double(const std::size_t i, const std::size_t j, const std::string& parameter, const double value) {
882 throw NotImplementedError("set_binary_interaction_double is not implemented for this backend");
883 };
885 virtual void set_binary_interaction_string(const std::string& CAS1, const std::string& CAS2, const std::string& parameter,
886 const std::string& value) {
887 throw NotImplementedError("set_binary_interaction_string is not implemented for this backend");
888 };
890 virtual void set_binary_interaction_string(const std::size_t i, const std::size_t j, const std::string& parameter, const std::string& value) {
891 throw NotImplementedError("set_binary_interaction_string is not implemented for this backend");
892 };
894 virtual double get_binary_interaction_double(const std::string& CAS1, const std::string& CAS2, const std::string& parameter) {
895 throw NotImplementedError("get_binary_interaction_double is not implemented for this backend");
896 };
898 virtual double get_binary_interaction_double(const std::size_t i, const std::size_t j, const std::string& parameter) {
899 throw NotImplementedError("get_binary_interaction_double is not implemented for this backend");
900 };
902 virtual std::string get_binary_interaction_string(const std::string& CAS1, const std::string& CAS2, const std::string& parameter) {
903 throw NotImplementedError("get_binary_interaction_string is not implemented for this backend");
904 };
906 virtual void apply_simple_mixing_rule(std::size_t i, std::size_t j, const std::string& model) {
907 throw NotImplementedError("apply_simple_mixing_rule is not implemented for this backend");
908 };
910 virtual void set_cubic_alpha_C(const size_t i, const std::string& parameter, const double c1, const double c2, const double c3) {
911 throw ValueError("set_cubic_alpha_C only defined for cubic backends");
912 };
914 virtual void set_fluid_parameter_double(const size_t i, const std::string& parameter, const double value) {
915 throw ValueError("set_fluid_parameter_double only defined for cubic backends");
916 };
918 virtual double get_fluid_parameter_double(const size_t i, const std::string& parameter) {
919 throw ValueError("get_fluid_parameter_double only defined for cubic backends");
920 };
921
923 virtual bool clear();
925 virtual bool clear_comp_change();
926
931 return _reducing;
932 };
933
935 const CoolProp::SimpleState& get_state(const std::string& state) {
936 return calc_state(state);
937 };
938
940 double Tmin(void);
942 double Tmax(void);
944 double pmax(void);
946 double Ttriple(void);
947
949 phases phase(void) {
950 return calc_phase();
951 };
955 };
957 void unspecify_phase(void) {
959 };
960
962 double T_critical(void);
964 double p_critical(void);
966 double rhomolar_critical(void);
968 double rhomass_critical(void);
969
971 std::vector<CriticalState> all_critical_points(void) {
973 };
974
978 };
979
982 return calc_get_spinodal_data();
983 };
984
986 void criticality_contour_values(double& L1star, double& M1star) {
987 return calc_criticality_contour_values(L1star, M1star);
988 }
989
1013 double tangent_plane_distance(const double T, const double p, const std::vector<double>& w, const double rhomolar_guess = -1) {
1014 return calc_tangent_plane_distance(T, p, w, rhomolar_guess);
1015 };
1016
1018 double T_reducing(void);
1020 double rhomolar_reducing(void);
1022 double rhomass_reducing(void);
1023
1025 double p_triple(void);
1026
1028 std::string name() {
1029 return calc_name();
1030 };
1032 std::string description() {
1033 return calc_description();
1034 };
1035
1037 double dipole_moment() {
1038 return calc_dipole_moment();
1039 }
1040
1041 // ----------------------------------------
1042 // Bulk properties - temperature and density are directly calculated every time
1043 // All other parameters are calculated on an as-needed basis
1044 // ----------------------------------------
1046 double keyed_output(parameters key);
1048 double trivial_keyed_output(parameters key);
1052 };
1056 };
1057
1059 double T(void) {
1060 return calc_T();
1061 };
1063 double rhomolar(void) {
1064 return calc_rhomolar();
1065 };
1067 double rhomass(void) {
1068 return calc_rhomass();
1069 };
1071 double p(void) {
1072 return _p;
1073 };
1075 double Q(void) {
1076 return _Q;
1077 };
1079 double tau(void);
1081 double delta(void);
1083 double molar_mass(void);
1085 double acentric_factor(void);
1087 double gas_constant(void);
1089 double Bvirial(void);
1091 double dBvirial_dT(void);
1093 double Cvirial(void);
1095 double dCvirial_dT(void);
1097 double compressibility_factor(void);
1099 double hmolar(void);
1101 double hmolar_residual(void);
1103 double hmass(void) {
1104 return calc_hmass();
1105 };
1107 double hmolar_excess(void);
1109 double hmass_excess(void) {
1110 return calc_hmass_excess();
1111 };
1113 double smolar(void);
1115 double smolar_residual(void);
1117 double neff(void);
1119 double smass(void) {
1120 return calc_smass();
1121 };
1123 double smolar_excess(void);
1125 double smass_excess(void) {
1126 return calc_smass_excess();
1127 };
1129 double umolar(void);
1131 double umass(void) {
1132 return calc_umass();
1133 };
1135 double umolar_excess(void);
1137 double umass_excess(void) {
1138 return calc_umass_excess();
1139 };
1141 double cpmolar(void);
1143 double cpmass(void) {
1144 return calc_cpmass();
1145 };
1147 double cp0molar(void);
1149 double cp0mass(void) {
1150 return calc_cp0mass();
1151 };
1153 double cvmolar(void);
1155 double cvmass(void) {
1156 return calc_cvmass();
1157 };
1159 double gibbsmolar(void);
1161 double gibbsmolar_residual(void);
1163 double gibbsmass(void) {
1164 return calc_gibbsmass();
1165 };
1167 double gibbsmolar_excess(void);
1169 double gibbsmass_excess(void) {
1170 return calc_gibbsmass_excess();
1171 };
1173 double helmholtzmolar(void);
1175 double helmholtzmass(void) {
1176 return calc_helmholtzmass();
1177 };
1179 double helmholtzmolar_excess(void);
1183 };
1185 double volumemolar_excess(void);
1187 double volumemass_excess(void) {
1188 return calc_volumemass_excess();
1189 };
1191 double speed_sound(void);
1193 double isothermal_compressibility(void);
1195 double isobaric_expansion_coefficient(void);
1199 double fugacity_coefficient(std::size_t i);
1201 std::vector<double> fugacity_coefficients();
1203 double fugacity(std::size_t i);
1205 double chemical_potential(std::size_t i);
1216 double PIP() {
1217 return calc_PIP();
1218 };
1219
1221 void true_critical_point(double& T, double& rho) {
1223 }
1224
1232 void ideal_curve(const std::string& type, std::vector<double>& T, std::vector<double>& p) {
1233 calc_ideal_curve(type, T, p);
1234 };
1235
1236 // ----------------------------------------
1237 // Partial derivatives
1238 // ----------------------------------------
1239
1245 return calc_first_partial_deriv(Of, Wrt, Constant);
1246 };
1247
1272 return calc_second_partial_deriv(Of1, Wrt1, Constant1, Wrt2, Constant2);
1273 };
1274
1297 return calc_first_saturation_deriv(Of1, Wrt1);
1298 };
1299
1318 return calc_second_saturation_deriv(Of1, Wrt1, Wrt2);
1319 };
1320
1341 return calc_first_two_phase_deriv(Of, Wrt, Constant);
1342 };
1343
1360 double second_two_phase_deriv(parameters Of, parameters Wrt1, parameters Constant1, parameters Wrt2, parameters Constant2) {
1361 return calc_second_two_phase_deriv(Of, Wrt1, Constant1, Wrt2, Constant2);
1362 };
1363
1384 double first_two_phase_deriv_splined(parameters Of, parameters Wrt, parameters Constant, double x_end) {
1385 return calc_first_two_phase_deriv_splined(Of, Wrt, Constant, x_end);
1386 };
1387
1388 // ----------------------------------------
1389 // Phase envelope for mixtures
1390 // ----------------------------------------
1391
1397 void build_phase_envelope(const std::string& type = "");
1402 return calc_phase_envelope_data();
1403 };
1404
1405 // ----------------------------------------
1406 // Ancillary equations
1407 // ----------------------------------------
1408
1410 virtual bool has_melting_line(void) {
1411 return false;
1412 };
1417 double melting_line(int param, int given, double value);
1423 double saturation_ancillary(parameters param, int Q, parameters given, double value);
1424
1425 // ----------------------------------------
1426 // Transport properties
1427 // ----------------------------------------
1429 double viscosity(void);
1431 void viscosity_contributions(CoolPropDbl& dilute, CoolPropDbl& initial_density, CoolPropDbl& residual, CoolPropDbl& critical) {
1432 calc_viscosity_contributions(dilute, initial_density, residual, critical);
1433 };
1435 double conductivity(void);
1437 void conductivity_contributions(CoolPropDbl& dilute, CoolPropDbl& initial_density, CoolPropDbl& residual, CoolPropDbl& critical) {
1438 calc_conductivity_contributions(dilute, initial_density, residual, critical);
1439 };
1441 double surface_tension(void);
1443 double Prandtl(void) {
1444 return cpmass() * viscosity() / conductivity();
1445 };
1452 void conformal_state(const std::string& reference_fluid, CoolPropDbl& T, CoolPropDbl& rhomolar) {
1453 return calc_conformal_state(reference_fluid, T, rhomolar);
1454 };
1455
1460 void change_EOS(const std::size_t i, const std::string& EOS_name) {
1461 calc_change_EOS(i, EOS_name);
1462 }
1463
1464 // ----------------------------------------
1465 // Helmholtz energy and derivatives
1466 // ----------------------------------------
1469 if (!_alpha0) _alpha0 = calc_alpha0();
1470 return _alpha0;
1471 };
1475 return _dalpha0_dDelta;
1476 };
1480 return _dalpha0_dTau;
1481 };
1485 return _d2alpha0_dDelta2;
1486 };
1490 return _d2alpha0_dDelta_dTau;
1491 };
1495 return _d2alpha0_dTau2;
1496 };
1500 return _d3alpha0_dTau3;
1501 };
1506 };
1511 };
1515 return _d3alpha0_dDelta3;
1516 };
1517
1520 if (!_alphar) _alphar = calc_alphar();
1521 return _alphar;
1522 };
1526 return _dalphar_dDelta;
1527 };
1531 return _dalphar_dTau;
1532 };
1536 return _d2alphar_dDelta2;
1537 };
1541 return _d2alphar_dDelta_dTau;
1542 };
1546 return _d2alphar_dTau2;
1547 };
1551 return _d3alphar_dDelta3;
1552 };
1557 };
1562 };
1566 return _d3alphar_dTau3;
1567 };
1571 return _d4alphar_dDelta4;
1572 };
1577 };
1582 };
1587 };
1591 return _d4alphar_dTau4;
1592 };
1593};
1594
1603{
1604 public:
1605 virtual AbstractState* get_AbstractState(const std::vector<std::string>& fluid_names) = 0;
1607};
1608
1612void register_backend(const backend_families& bf, shared_ptr<AbstractStateGenerator> gen);
1613
1614template <class T>
1616{
1617 public:
1619 register_backend(bf, shared_ptr<AbstractStateGenerator>(new T()));
1620 };
1621};
1622
1623} /* namespace CoolProp */
1624#endif /* ABSTRACTSTATE_H_ */