CoolProp 7.2.0
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
121
124
128
131
133
136
145
148
151
152 // ----------------------------------------
153 // Property accessors to be optionally implemented by the backend
154 // for properties that are not always calculated
155 // ----------------------------------------
157 virtual CoolPropDbl calc_hmolar(void) {
158 throw NotImplementedError("calc_hmolar is not implemented for this backend");
159 };
162 throw NotImplementedError("calc_hmolar_residual is not implemented for this backend");
163 };
165 virtual CoolPropDbl calc_smolar(void) {
166 throw NotImplementedError("calc_smolar is not implemented for this backend");
167 };
170 throw NotImplementedError("calc_smolar_residual is not implemented for this backend");
171 };
173 virtual CoolPropDbl calc_neff(void) {
174 throw NotImplementedError("calc_neff is not implemented for this backend");
175 };
177 virtual CoolPropDbl calc_umolar(void) {
178 throw NotImplementedError("calc_umolar is not implemented for this backend");
179 };
181 virtual CoolPropDbl calc_cpmolar(void) {
182 throw NotImplementedError("calc_cpmolar is not implemented for this backend");
183 };
186 throw NotImplementedError("calc_cpmolar_idealgas is not implemented for this backend");
187 };
189 virtual CoolPropDbl calc_cvmolar(void) {
190 throw NotImplementedError("calc_cvmolar is not implemented for this backend");
191 };
194 throw NotImplementedError("calc_gibbsmolar is not implemented for this backend");
195 };
198 throw NotImplementedError("calc_gibbsmolar_residual is not implemented for this backend");
199 };
202 throw NotImplementedError("calc_helmholtzmolar is not implemented for this backend");
203 };
206 throw NotImplementedError("calc_speed_sound is not implemented for this backend");
207 };
210 throw NotImplementedError("calc_isothermal_compressibility is not implemented for this backend");
211 };
214 throw NotImplementedError("calc_isobaric_expansion_coefficient is not implemented for this backend");
215 };
218 throw NotImplementedError("calc_isentropic_expansion_coefficient is not implemented for this backend");
219 };
222 throw NotImplementedError("calc_viscosity is not implemented for this backend");
223 };
226 throw NotImplementedError("calc_conductivity is not implemented for this backend");
227 };
230 throw NotImplementedError("calc_surface_tension is not implemented for this backend");
231 };
234 throw NotImplementedError("calc_molar_mass is not implemented for this backend");
235 };
238 throw NotImplementedError("calc_acentric_factor is not implemented for this backend");
239 };
242 throw NotImplementedError("calc_pressure is not implemented for this backend");
243 };
246 throw NotImplementedError("calc_gas_constant is not implemented for this backend");
247 };
249 virtual CoolPropDbl calc_fugacity_coefficient(std::size_t i) {
250 throw NotImplementedError("calc_fugacity_coefficient is not implemented for this backend");
251 };
253 virtual std::vector<CoolPropDbl> calc_fugacity_coefficients() {
254 throw NotImplementedError("calc_fugacity_coefficients is not implemented for this backend");
255 };
257 virtual CoolPropDbl calc_fugacity(std::size_t i) {
258 throw NotImplementedError("calc_fugacity is not implemented for this backend");
259 };
261 virtual CoolPropDbl calc_chemical_potential(std::size_t i) {
262 throw NotImplementedError("calc_chemical_potential is not implemented for this backend");
263 };
265 virtual CoolPropDbl calc_PIP(void) {
266 throw NotImplementedError("calc_PIP is not implemented for this backend");
267 };
268
269 // Excess properties
271 virtual void calc_excess_properties(void) {
272 throw NotImplementedError("calc_excess_properties is not implemented for this backend");
273 };
274
275 // Derivatives of residual helmholtz energy
277 virtual CoolPropDbl calc_alphar(void) {
278 throw NotImplementedError("calc_alphar is not implemented for this backend");
279 };
282 throw NotImplementedError("calc_dalphar_dDelta is not implemented for this backend");
283 };
286 throw NotImplementedError("calc_dalphar_dTau is not implemented for this backend");
287 };
290 throw NotImplementedError("calc_d2alphar_dDelta2 is not implemented for this backend");
291 };
294 throw NotImplementedError("calc_d2alphar_dDelta_dTau is not implemented for this backend");
295 };
298 throw NotImplementedError("calc_d2alphar_dTau2 is not implemented for this backend");
299 };
302 throw NotImplementedError("calc_d3alphar_dDelta3 is not implemented for this backend");
303 };
306 throw NotImplementedError("calc_d3alphar_dDelta2_dTau is not implemented for this backend");
307 };
310 throw NotImplementedError("calc_d3alphar_dDelta_dTau2 is not implemented for this backend");
311 };
314 throw NotImplementedError("calc_d3alphar_dTau3 is not implemented for this backend");
315 };
316
319 throw NotImplementedError("calc_d4alphar_dDelta4 is not implemented for this backend");
320 };
323 throw NotImplementedError("calc_d4alphar_dDelta3_dTau is not implemented for this backend");
324 };
327 throw NotImplementedError("calc_d4alphar_dDelta2_dTau2 is not implemented for this backend");
328 };
331 throw NotImplementedError("calc_d4alphar_dDelta_dTau3 is not implemented for this backend");
332 };
335 throw NotImplementedError("calc_d4alphar_dTau4 is not implemented for this backend");
336 };
337
338 // Derivatives of ideal-gas helmholtz energy
340 virtual CoolPropDbl calc_alpha0(void) {
341 throw NotImplementedError("calc_alpha0 is not implemented for this backend");
342 };
345 throw NotImplementedError("calc_dalpha0_dDelta is not implemented for this backend");
346 };
349 throw NotImplementedError("calc_dalpha0_dTau is not implemented for this backend");
350 };
353 throw NotImplementedError("calc_d2alpha0_dDelta_dTau is not implemented for this backend");
354 };
357 throw NotImplementedError("calc_d2alpha0_dDelta2 is not implemented for this backend");
358 };
361 throw NotImplementedError("calc_d2alpha0_dTau2 is not implemented for this backend");
362 };
365 throw NotImplementedError("calc_d3alpha0_dDelta3 is not implemented for this backend");
366 };
369 throw NotImplementedError("calc_d3alpha0_dDelta2_dTau is not implemented for this backend");
370 };
373 throw NotImplementedError("calc_d3alpha0_dDelta_dTau2 is not implemented for this backend");
374 };
377 throw NotImplementedError("calc_d3alpha0_dTau3 is not implemented for this backend");
378 };
379
380 virtual void calc_reducing_state(void) {
381 throw NotImplementedError("calc_reducing_state is not implemented for this backend");
382 };
383
385 virtual CoolPropDbl calc_Tmax(void) {
386 throw NotImplementedError("calc_Tmax is not implemented for this backend");
387 };
389 virtual CoolPropDbl calc_Tmin(void) {
390 throw NotImplementedError("calc_Tmin is not implemented for this backend");
391 };
393 virtual CoolPropDbl calc_pmax(void) {
394 throw NotImplementedError("calc_pmax is not implemented for this backend");
395 };
396
398 virtual CoolPropDbl calc_GWP20(void) {
399 throw NotImplementedError("calc_GWP20 is not implemented for this backend");
400 };
402 virtual CoolPropDbl calc_GWP100(void) {
403 throw NotImplementedError("calc_GWP100 is not implemented for this backend");
404 };
406 virtual CoolPropDbl calc_GWP500(void) {
407 throw NotImplementedError("calc_GWP500 is not implemented for this backend");
408 };
410 virtual CoolPropDbl calc_ODP(void) {
411 throw NotImplementedError("calc_ODP is not implemented for this backend");
412 };
415 throw NotImplementedError("calc_flame_hazard is not implemented for this backend");
416 };
419 throw NotImplementedError("calc_health_hazard is not implemented for this backend");
420 };
423 throw NotImplementedError("calc_physical_hazard is not implemented for this backend");
424 };
427 throw NotImplementedError("calc_dipole_moment is not implemented for this backend");
428 };
429
433 virtual CoolPropDbl calc_second_partial_deriv(parameters Of1, parameters Wrt1, parameters Constant1, parameters Wrt2, parameters Constant2);
434
437 throw NotImplementedError("calc_reduced_density is not implemented for this backend");
438 };
441 throw NotImplementedError("calc_reciprocal_reduced_temperature is not implemented for this backend");
442 };
443
445 virtual CoolPropDbl calc_Bvirial(void) {
446 throw NotImplementedError("calc_Bvirial is not implemented for this backend");
447 };
449 virtual CoolPropDbl calc_Cvirial(void) {
450 throw NotImplementedError("calc_Cvirial is not implemented for this backend");
451 };
454 throw NotImplementedError("calc_dBvirial_dT is not implemented for this backend");
455 };
458 throw NotImplementedError("calc_dCvirial_dT is not implemented for this backend");
459 };
462 throw NotImplementedError("calc_compressibility_factor is not implemented for this backend");
463 };
464
466 virtual std::string calc_name(void) {
467 throw NotImplementedError("calc_name is not implemented for this backend");
468 };
470 virtual std::string calc_description(void) {
471 throw NotImplementedError("calc_description is not implemented for this backend");
472 };
473
475 virtual CoolPropDbl calc_Ttriple(void) {
476 throw NotImplementedError("calc_Ttriple is not implemented for this backend");
477 };
480 throw NotImplementedError("calc_p_triple is not implemented for this backend");
481 };
482
485 throw NotImplementedError("calc_T_critical is not implemented for this backend");
486 };
489 throw NotImplementedError("calc_T_reducing is not implemented for this backend");
490 };
493 throw NotImplementedError("calc_p_critical is not implemented for this backend");
494 };
497 throw NotImplementedError("calc_p_reducing is not implemented for this backend");
498 };
501 throw NotImplementedError("calc_rhomolar_critical is not implemented for this backend");
502 };
505 throw NotImplementedError("calc_rhomass_critical is not implemented for this backend");
506 };
509 throw NotImplementedError("calc_rhomolar_reducing is not implemented for this backend");
510 };
512 virtual void calc_phase_envelope(const std::string& type) {
513 throw NotImplementedError("calc_phase_envelope is not implemented for this backend");
514 };
516 virtual CoolPropDbl calc_rhomass(void) {
517 return rhomolar() * molar_mass();
518 }
519 virtual CoolPropDbl calc_hmass(void) {
520 return hmolar() / molar_mass();
521 }
523 return hmolar_excess() / molar_mass();
524 }
525 virtual CoolPropDbl calc_smass(void) {
526 return smolar() / molar_mass();
527 }
529 return smolar_excess() / molar_mass();
530 }
531 virtual CoolPropDbl calc_cpmass(void) {
532 return cpmolar() / molar_mass();
533 }
534 virtual CoolPropDbl calc_cp0mass(void) {
535 return cp0molar() / molar_mass();
536 }
537 virtual CoolPropDbl calc_cvmass(void) {
538 return cvmolar() / molar_mass();
539 }
540 virtual CoolPropDbl calc_umass(void) {
541 return umolar() / molar_mass();
542 }
544 return umolar_excess() / molar_mass();
545 }
547 return gibbsmolar() / molar_mass();
548 }
550 return gibbsmolar_excess() / molar_mass();
551 }
553 return helmholtzmolar() / molar_mass();
554 }
557 }
559 return volumemolar_excess() / molar_mass();
560 }
561
563 virtual void update_states(void) {
564 throw NotImplementedError("This backend does not implement update_states function");
565 };
566
567 virtual CoolPropDbl calc_melting_line(int param, int given, CoolPropDbl value) {
568 throw NotImplementedError("This backend does not implement calc_melting_line function");
569 };
570
575 virtual CoolPropDbl calc_saturation_ancillary(parameters param, int Q, parameters given, double value) {
576 throw NotImplementedError("This backend does not implement calc_saturation_ancillary");
577 };
578
580 virtual phases calc_phase(void) {
581 throw NotImplementedError("This backend does not implement calc_phase function");
582 };
585 throw NotImplementedError("This backend does not implement calc_specify_phase function");
586 };
588 virtual void calc_unspecify_phase(void) {
589 throw NotImplementedError("This backend does not implement calc_unspecify_phase function");
590 };
592 virtual std::vector<std::string> calc_fluid_names(void) {
593 throw NotImplementedError("This backend does not implement calc_fluid_names function");
594 };
597 virtual const CoolProp::SimpleState& calc_state(const std::string& state) {
598 throw NotImplementedError("calc_state is not implemented for this backend");
599 };
600
602 throw NotImplementedError("calc_phase_envelope_data is not implemented for this backend");
603 };
604
605 virtual std::vector<CoolPropDbl> calc_mole_fractions_liquid(void) {
606 throw NotImplementedError("calc_mole_fractions_liquid is not implemented for this backend");
607 };
608 virtual std::vector<CoolPropDbl> calc_mole_fractions_vapor(void) {
609 throw NotImplementedError("calc_mole_fractions_vapor is not implemented for this backend");
610 };
611 virtual const std::vector<CoolPropDbl> calc_mass_fractions(void) {
612 throw NotImplementedError("calc_mass_fractions is not implemented for this backend");
613 };
614
617 throw NotImplementedError("calc_fraction_min is not implemented for this backend");
618 };
621 throw NotImplementedError("calc_fraction_max is not implemented for this backend");
622 };
624 throw NotImplementedError("calc_T_freeze is not implemented for this backend");
625 };
626
628 throw NotImplementedError("calc_first_saturation_deriv is not implemented for this backend");
629 };
631 throw NotImplementedError("calc_second_saturation_deriv is not implemented for this backend");
632 };
634 throw NotImplementedError("calc_first_two_phase_deriv is not implemented for this backend");
635 };
637 throw NotImplementedError("calc_second_two_phase_deriv is not implemented for this backend");
638 };
640 throw NotImplementedError("calc_first_two_phase_deriv_splined is not implemented for this backend");
641 };
642
644 throw NotImplementedError("calc_saturated_liquid_keyed_output is not implemented for this backend");
645 };
647 throw NotImplementedError("calc_saturated_vapor_keyed_output is not implemented for this backend");
648 };
649 virtual void calc_ideal_curve(const std::string& type, std::vector<double>& T, std::vector<double>& p) {
650 throw NotImplementedError("calc_ideal_curve is not implemented for this backend");
651 };
652
654 virtual CoolPropDbl calc_T(void) {
655 return _T;
656 }
659 return _rhomolar;
660 }
661
663 virtual double calc_tangent_plane_distance(const double T, const double p, const std::vector<double>& w, const double rhomolar_guess) {
664 throw NotImplementedError("calc_tangent_plane_distance is not implemented for this backend");
665 };
666
668 virtual void calc_true_critical_point(double& T, double& rho) {
669 throw NotImplementedError("calc_true_critical_point is not implemented for this backend");
670 };
671
672 virtual void calc_conformal_state(const std::string& reference_fluid, CoolPropDbl& T, CoolPropDbl& rhomolar) {
673 throw NotImplementedError("calc_conformal_state is not implemented for this backend");
674 };
675
676 virtual void calc_viscosity_contributions(CoolPropDbl& dilute, CoolPropDbl& initial_density, CoolPropDbl& residual, CoolPropDbl& critical) {
677 throw NotImplementedError("calc_viscosity_contributions is not implemented for this backend");
678 };
679 virtual void calc_conductivity_contributions(CoolPropDbl& dilute, CoolPropDbl& initial_density, CoolPropDbl& residual, CoolPropDbl& critical) {
680 throw NotImplementedError("calc_conductivity_contributions is not implemented for this backend");
681 };
682 virtual std::vector<CriticalState> calc_all_critical_points(void) {
683 throw NotImplementedError("calc_all_critical_points is not implemented for this backend");
684 };
685 virtual void calc_build_spinodal() {
686 throw NotImplementedError("calc_build_spinodal is not implemented for this backend");
687 };
689 throw NotImplementedError("calc_get_spinodal_data is not implemented for this backend");
690 };
691 virtual void calc_criticality_contour_values(double& L1star, double& M1star) {
692 throw NotImplementedError("calc_criticality_contour_values is not implemented for this backend");
693 };
694
696 virtual void mass_to_molar_inputs(CoolProp::input_pairs& input_pair, CoolPropDbl& value1, CoolPropDbl& value2);
697
699 virtual void calc_change_EOS(const std::size_t i, const std::string& EOS_name) {
700 throw NotImplementedError("calc_change_EOS is not implemented for this backend");
701 };
702
703 public:
705 clear();
706 }
707 virtual ~AbstractState() {};
708
710
716 static AbstractState* factory(const std::string& backend, const std::string& fluid_names) {
717 return factory(backend, strsplit(fluid_names, '&'));
718 };
719
743 static AbstractState* factory(const std::string& backend, const std::vector<std::string>& fluid_names);
744
747 _T = T;
748 }
749
754 virtual std::string backend_name(void) = 0;
755
756 // The derived classes must implement this function to define whether they use mole fractions (true) or mass fractions (false)
757 virtual bool using_mole_fractions(void) = 0;
758 virtual bool using_mass_fractions(void) = 0;
759 virtual bool using_volu_fractions(void) = 0;
760
761 virtual void set_mole_fractions(const std::vector<CoolPropDbl>& mole_fractions) = 0;
762 virtual void set_mass_fractions(const std::vector<CoolPropDbl>& mass_fractions) = 0;
763 virtual void set_volu_fractions(const std::vector<CoolPropDbl>& mass_fractions) {
764 throw NotImplementedError("Volume composition has not been implemented.");
765 }
766
786 virtual void set_reference_stateS(const std::string& reference_state) {
788 "Setting reference state has not been implemented for this backend. Try using CoolProp::set_reference_stateD instead.");
789 }
790
796 virtual void set_reference_stateD(double T, double rhomolar, double hmolar0, double smolar0) {
798 "Setting reference state has not been implemented for this backend. Try using CoolProp::set_reference_stateD instead.");
799 }
800
801#ifndef COOLPROPDBL_MAPS_TO_DOUBLE
802 void set_mole_fractions(const std::vector<double>& mole_fractions) {
803 set_mole_fractions(std::vector<CoolPropDbl>(mole_fractions.begin(), mole_fractions.end()));
804 };
805 void set_mass_fractions(const std::vector<double>& mass_fractions) {
806 set_mass_fractions(std::vector<CoolPropDbl>(mass_fractions.begin(), mass_fractions.end()));
807 };
808 void set_volu_fractions(const std::vector<double>& volu_fractions) {
809 set_volu_fractions(std::vector<CoolPropDbl>(volu_fractions.begin(), volu_fractions.end()));
810 };
811#endif
812
813#ifdef EMSCRIPTEN
814 void set_mole_fractions_double(const std::vector<double>& mole_fractions) {
815 set_mole_fractions(std::vector<CoolPropDbl>(mole_fractions.begin(), mole_fractions.end()));
816 };
817#endif
818
820 std::vector<CoolPropDbl> mole_fractions_liquid(void) {
822 };
824 std::vector<double> mole_fractions_liquid_double(void) {
825 std::vector<CoolPropDbl> x = calc_mole_fractions_liquid();
826 return std::vector<double>(x.begin(), x.end());
827 };
828
830 std::vector<CoolPropDbl> mole_fractions_vapor(void) {
832 };
834 std::vector<double> mole_fractions_vapor_double(void) {
835 std::vector<CoolPropDbl> y = calc_mole_fractions_vapor();
836 return std::vector<double>(y.begin(), y.end());
837 };
838
840 virtual const std::vector<CoolPropDbl>& get_mole_fractions(void) = 0;
842 virtual const std::vector<CoolPropDbl> get_mass_fractions(void) {
843 return this->calc_mass_fractions();
844 };
845
847 virtual void update(CoolProp::input_pairs input_pair, double Value1, double Value2) = 0;
848
850 virtual void update_QT_pure_superanc(double Q, double T) {
851 throw NotImplementedError("update_QT_pure_superanc is not implemented for this backend");
852 };
853
856 virtual void update_with_guesses(CoolProp::input_pairs input_pair, double Value1, double Value2, const GuessesStructure& guesses) {
857 throw NotImplementedError("update_with_guesses is not implemented for this backend");
858 };
859
863 virtual bool available_in_high_level(void) {
864 return true;
865 }
866
868 virtual std::string fluid_param_string(const std::string&) {
869 throw NotImplementedError("fluid_param_string has not been implemented for this backend");
870 }
871
873 std::vector<std::string> fluid_names(void);
874
879 virtual const double get_fluid_constant(std::size_t i, parameters param) const {
880 throw NotImplementedError("get_fluid_constant is not implemented for this backend");
881 };
882
884 virtual void set_binary_interaction_double(const std::string& CAS1, const std::string& CAS2, const std::string& parameter, const double value) {
885 throw NotImplementedError("set_binary_interaction_double is not implemented for this backend");
886 };
888 virtual void set_binary_interaction_double(const std::size_t i, const std::size_t j, const std::string& parameter, const double value) {
889 throw NotImplementedError("set_binary_interaction_double is not implemented for this backend");
890 };
892 virtual void set_binary_interaction_string(const std::string& CAS1, const std::string& CAS2, const std::string& parameter,
893 const std::string& value) {
894 throw NotImplementedError("set_binary_interaction_string is not implemented for this backend");
895 };
897 virtual void set_binary_interaction_string(const std::size_t i, const std::size_t j, const std::string& parameter, const std::string& value) {
898 throw NotImplementedError("set_binary_interaction_string is not implemented for this backend");
899 };
901 virtual double get_binary_interaction_double(const std::string& CAS1, const std::string& CAS2, const std::string& parameter) {
902 throw NotImplementedError("get_binary_interaction_double is not implemented for this backend");
903 };
905 virtual double get_binary_interaction_double(const std::size_t i, const std::size_t j, const std::string& parameter) {
906 throw NotImplementedError("get_binary_interaction_double is not implemented for this backend");
907 };
909 virtual std::string get_binary_interaction_string(const std::string& CAS1, const std::string& CAS2, const std::string& parameter) {
910 throw NotImplementedError("get_binary_interaction_string is not implemented for this backend");
911 };
913 virtual void apply_simple_mixing_rule(std::size_t i, std::size_t j, const std::string& model) {
914 throw NotImplementedError("apply_simple_mixing_rule is not implemented for this backend");
915 };
917 virtual void set_cubic_alpha_C(const size_t i, const std::string& parameter, const double c1, const double c2, const double c3) {
918 throw ValueError("set_cubic_alpha_C only defined for cubic backends");
919 };
921 virtual void set_fluid_parameter_double(const size_t i, const std::string& parameter, const double value) {
922 throw ValueError("set_fluid_parameter_double only defined for cubic backends");
923 };
925 virtual double get_fluid_parameter_double(const size_t i, const std::string& parameter) {
926 throw ValueError("get_fluid_parameter_double only defined for cubic backends");
927 };
928
930 virtual bool clear();
932 virtual bool clear_comp_change();
933
938 return _reducing;
939 };
940
942 const CoolProp::SimpleState& get_state(const std::string& state) {
943 return calc_state(state);
944 };
945
947 double Tmin(void);
949 double Tmax(void);
951 double pmax(void);
953 double Ttriple(void);
954
956 phases phase(void) {
957 return calc_phase();
958 };
962 };
964 void unspecify_phase(void) {
966 };
967
969 double T_critical(void);
971 double p_critical(void);
973 double rhomolar_critical(void);
975 double rhomass_critical(void);
976
978 std::vector<CriticalState> all_critical_points(void) {
980 };
981
985 };
986
989 return calc_get_spinodal_data();
990 };
991
993 void criticality_contour_values(double& L1star, double& M1star) {
994 return calc_criticality_contour_values(L1star, M1star);
995 }
996
1020 double tangent_plane_distance(const double T, const double p, const std::vector<double>& w, const double rhomolar_guess = -1) {
1021 return calc_tangent_plane_distance(T, p, w, rhomolar_guess);
1022 };
1023
1025 double T_reducing(void);
1027 double rhomolar_reducing(void);
1029 double rhomass_reducing(void);
1030
1032 double p_triple(void);
1033
1035 std::string name() {
1036 return calc_name();
1037 };
1039 std::string description() {
1040 return calc_description();
1041 };
1042
1044 double dipole_moment() {
1045 return calc_dipole_moment();
1046 }
1047
1048 // ----------------------------------------
1049 // Bulk properties - temperature and density are directly calculated every time
1050 // All other parameters are calculated on an as-needed basis
1051 // ----------------------------------------
1053 double keyed_output(parameters key);
1055 double trivial_keyed_output(parameters key);
1059 };
1063 };
1064
1066 double T(void) {
1067 return calc_T();
1068 };
1070 double rhomolar(void) {
1071 return calc_rhomolar();
1072 };
1074 double rhomass(void) {
1075 return calc_rhomass();
1076 };
1078 double p(void) {
1079 return _p;
1080 };
1082 double Q(void) {
1083 return _Q;
1084 };
1086 double tau(void);
1088 double delta(void);
1090 double molar_mass(void);
1092 double acentric_factor(void);
1094 double gas_constant(void);
1096 double Bvirial(void);
1098 double dBvirial_dT(void);
1100 double Cvirial(void);
1102 double dCvirial_dT(void);
1104 double compressibility_factor(void);
1106 double hmolar(void);
1108 double hmolar_residual(void);
1110 double hmolar_idealgas(void);
1112 double hmass_idealgas(void);
1114 double hmass(void) {
1115 return calc_hmass();
1116 };
1118 double hmolar_excess(void);
1120 double hmass_excess(void) {
1121 return calc_hmass_excess();
1122 };
1124 double smolar(void);
1126 double smolar_residual(void);
1128 double smolar_idealgas(void);
1130 double smass_idealgas(void);
1132 double neff(void);
1134 double smass(void) {
1135 return calc_smass();
1136 };
1138 double smolar_excess(void);
1140 double smass_excess(void) {
1141 return calc_smass_excess();
1142 };
1144 double umolar(void);
1146 double umass(void) {
1147 return calc_umass();
1148 };
1150 double umolar_excess(void);
1152 double umass_excess(void) {
1153 return calc_umass_excess();
1154 };
1156 double umolar_idealgas(void);
1158 double umass_idealgas(void);
1160 double cpmolar(void);
1162 double cpmass(void) {
1163 return calc_cpmass();
1164 };
1166 double cp0molar(void);
1168 double cp0mass(void) {
1169 return calc_cp0mass();
1170 };
1172 double cvmolar(void);
1174 double cvmass(void) {
1175 return calc_cvmass();
1176 };
1178 double gibbsmolar(void);
1180 double gibbsmolar_residual(void);
1182 double gibbsmass(void) {
1183 return calc_gibbsmass();
1184 };
1186 double gibbsmolar_excess(void);
1188 double gibbsmass_excess(void) {
1189 return calc_gibbsmass_excess();
1190 };
1192 double helmholtzmolar(void);
1194 double helmholtzmass(void) {
1195 return calc_helmholtzmass();
1196 };
1198 double helmholtzmolar_excess(void);
1202 };
1204 double volumemolar_excess(void);
1206 double volumemass_excess(void) {
1207 return calc_volumemass_excess();
1208 };
1210 double speed_sound(void);
1212 double isothermal_compressibility(void);
1214 double isobaric_expansion_coefficient(void);
1218 double fugacity_coefficient(std::size_t i);
1220 std::vector<double> fugacity_coefficients();
1222 double fugacity(std::size_t i);
1224 double chemical_potential(std::size_t i);
1235 double PIP() {
1236 return calc_PIP();
1237 };
1238
1240 void true_critical_point(double& T, double& rho) {
1242 }
1243
1251 void ideal_curve(const std::string& type, std::vector<double>& T, std::vector<double>& p) {
1252 calc_ideal_curve(type, T, p);
1253 };
1254
1255 // ----------------------------------------
1256 // Partial derivatives
1257 // ----------------------------------------
1258
1264 return calc_first_partial_deriv(Of, Wrt, Constant);
1265 };
1266
1291 return calc_second_partial_deriv(Of1, Wrt1, Constant1, Wrt2, Constant2);
1292 };
1293
1316 return calc_first_saturation_deriv(Of1, Wrt1);
1317 };
1318
1337 return calc_second_saturation_deriv(Of1, Wrt1, Wrt2);
1338 };
1339
1360 return calc_first_two_phase_deriv(Of, Wrt, Constant);
1361 };
1362
1379 double second_two_phase_deriv(parameters Of, parameters Wrt1, parameters Constant1, parameters Wrt2, parameters Constant2) {
1380 return calc_second_two_phase_deriv(Of, Wrt1, Constant1, Wrt2, Constant2);
1381 };
1382
1403 double first_two_phase_deriv_splined(parameters Of, parameters Wrt, parameters Constant, double x_end) {
1404 return calc_first_two_phase_deriv_splined(Of, Wrt, Constant, x_end);
1405 };
1406
1407 // ----------------------------------------
1408 // Phase envelope for mixtures
1409 // ----------------------------------------
1410
1416 void build_phase_envelope(const std::string& type = "");
1421 return calc_phase_envelope_data();
1422 };
1423
1424 // ----------------------------------------
1425 // Ancillary equations
1426 // ----------------------------------------
1427
1429 virtual bool has_melting_line(void) {
1430 return false;
1431 };
1436 double melting_line(int param, int given, double value);
1442 double saturation_ancillary(parameters param, int Q, parameters given, double value);
1443
1444 // ----------------------------------------
1445 // Transport properties
1446 // ----------------------------------------
1448 double viscosity(void);
1450 void viscosity_contributions(CoolPropDbl& dilute, CoolPropDbl& initial_density, CoolPropDbl& residual, CoolPropDbl& critical) {
1451 calc_viscosity_contributions(dilute, initial_density, residual, critical);
1452 };
1454 double conductivity(void);
1456 void conductivity_contributions(CoolPropDbl& dilute, CoolPropDbl& initial_density, CoolPropDbl& residual, CoolPropDbl& critical) {
1457 calc_conductivity_contributions(dilute, initial_density, residual, critical);
1458 };
1460 double surface_tension(void);
1462 double Prandtl(void) {
1463 return cpmass() * viscosity() / conductivity();
1464 };
1471 void conformal_state(const std::string& reference_fluid, CoolPropDbl& T, CoolPropDbl& rhomolar) {
1472 return calc_conformal_state(reference_fluid, T, rhomolar);
1473 };
1474
1479 void change_EOS(const std::size_t i, const std::string& EOS_name) {
1480 calc_change_EOS(i, EOS_name);
1481 }
1482
1483 // ----------------------------------------
1484 // Helmholtz energy and derivatives
1485 // ----------------------------------------
1488 if (!_alpha0) _alpha0 = calc_alpha0();
1489 return _alpha0;
1490 };
1494 return _dalpha0_dDelta;
1495 };
1499 return _dalpha0_dTau;
1500 };
1504 return _d2alpha0_dDelta2;
1505 };
1509 return _d2alpha0_dDelta_dTau;
1510 };
1514 return _d2alpha0_dTau2;
1515 };
1519 return _d3alpha0_dTau3;
1520 };
1525 };
1530 };
1534 return _d3alpha0_dDelta3;
1535 };
1536
1539 if (!_alphar) _alphar = calc_alphar();
1540 return _alphar;
1541 };
1545 return _dalphar_dDelta;
1546 };
1550 return _dalphar_dTau;
1551 };
1555 return _d2alphar_dDelta2;
1556 };
1560 return _d2alphar_dDelta_dTau;
1561 };
1565 return _d2alphar_dTau2;
1566 };
1570 return _d3alphar_dDelta3;
1571 };
1576 };
1581 };
1585 return _d3alphar_dTau3;
1586 };
1590 return _d4alphar_dDelta4;
1591 };
1596 };
1601 };
1606 };
1610 return _d4alphar_dTau4;
1611 };
1612};
1613
1622{
1623 public:
1624 virtual AbstractState* get_AbstractState(const std::vector<std::string>& fluid_names) = 0;
1626};
1627
1631void register_backend(const backend_families& bf, shared_ptr<AbstractStateGenerator> gen);
1632
1633template <class T>
1635{
1636 public:
1638 register_backend(bf, shared_ptr<AbstractStateGenerator>(new T()));
1639 };
1640};
1641
1642} /* namespace CoolProp */
1643#endif /* ABSTRACTSTATE_H_ */