teqp 0.19.1
Loading...
Searching...
No Matches
teqpcpp.hpp
Go to the documentation of this file.
1#pragma once
2#include <memory>
3#include <typeindex>
4#include <optional>
5
6#include <Eigen/Dense>
7#include "nlohmann/json.hpp"
8
9// The only headers that can be included here are
10// ones that define and use POD (plain ole' data) types
14
15using EArray2 = Eigen::Array<double, 2, 1>;
16using EArrayd = Eigen::ArrayX<double>;
17using EArray33d = Eigen::Array<double, 3, 3>;
18using REArrayd = Eigen::Ref<const EArrayd>;
19using EMatrixd = Eigen::Array<double, Eigen::Dynamic, Eigen::Dynamic>;
20using REMatrixd = Eigen::Ref<const Eigen::Array<double, Eigen::Dynamic, Eigen::Dynamic>>;
21
22#define ARXY_args \
23 X(0,0) \
24 X(0,1) \
25 X(0,2) \
26 X(0,3) \
27 X(0,4) \
28 X(1,0) \
29 X(1,1) \
30 X(1,2) \
31 X(1,3) \
32 X(1,4) \
33 X(2,0) \
34 X(2,1) \
35 X(2,2) \
36 X(2,3) \
37 X(2,4)
38
39#define AR0N_args \
40 X(0) \
41 X(1) \
42 X(2) \
43 X(3) \
44 X(4) \
45 X(5) \
46 X(6)
47
48// Functions that return a double, take T and rhovec as arguments
49#define ISOCHORIC_double_args \
50 X(get_pr) \
51 X(get_splus) \
52 X(get_dpdT_constrhovec)
53
54#define ISOCHORIC_array_args \
55 X(build_Psir_gradient_autodiff) \
56 X(get_chempotVLE_autodiff) \
57 X(get_dchempotdT_autodiff) \
58 X(get_fugacity_coefficients) \
59 X(get_partial_molar_volumes) \
60 X(build_d2PsirdTdrhoi_autodiff) \
61 X(get_dpdrhovec_constT)
62
63#define ISOCHORIC_matrix_args \
64 X(build_Psir_Hessian_autodiff) \
65 X(build_Psi_Hessian_autodiff)
66
67#define ISOCHORIC_multimatrix_args \
68 X(build_Psir_fgradHessian_autodiff)
69
70namespace teqp {
71 namespace cppinterface {
72
88 public:
89
90 virtual ~AbstractModel() = default;
91
92 virtual const std::type_index& get_type_index() const = 0;
93
94 virtual double get_R(const EArrayd&) const = 0;
95 double R(const EArrayd& x) const { return get_R(x); };
96
97 virtual double get_Arxy(const int, const int, const double, const double, const EArrayd&) const = 0;
98
99 // Here X-Macros are used to create functions like get_Ar00, get_Ar01, ....
100 #define X(i,j) virtual double get_Ar ## i ## j(const double T, const double rho, const REArrayd& molefrac) const = 0;
102 #undef X
103 // And like get_Ar01n, get_Ar02n, ....
104 #define X(i) virtual EArrayd get_Ar0 ## i ## n(const double T, const double rho, const REArrayd& molefrac) const = 0;
106 #undef X
107
108
109 // Virial derivatives
110 virtual double get_B2vir(const double T, const EArrayd& z) const = 0;
111 virtual std::map<int, double> get_Bnvir(const int Nderiv, const double T, const EArrayd& z) const = 0;
112 virtual double get_B12vir(const double T, const EArrayd& z) const = 0;
113 virtual double get_dmBnvirdTm(const int Nderiv, const int NTderiv, const double T, const EArrayd& z) const = 0;
114
115 // Composition derivatives
116 virtual double get_ATrhoXi(const double T, const int NT, const double rhomolar, int ND, const EArrayd& molefrac, const int i, const int NXi) const = 0;
117 virtual double get_ATrhoXiXj(const double T, const int NT, const double rhomolar, int ND, const EArrayd& molefrac, const int i, const int NXi, const int j, const int NXj) const = 0;
118 virtual double get_ATrhoXiXjXk(const double T, const int NT, const double rhomolar, int ND, const EArrayd& molefrac, const int i, const int NXi, const int j, const int NXj, const int k, const int NXk) const = 0;
119
120 virtual double get_AtaudeltaXi(const double tau, const int Ntau, const double delta, int Ndelta, const EArrayd& molefrac, const int i, const int NXi) const = 0;
121 virtual double get_AtaudeltaXiXj(const double tau, const int Ntau, const double delta, int Ndelta, const EArrayd& molefrac, const int i, const int NXi, const int j, const int NXj) const = 0;
122 virtual double get_AtaudeltaXiXjXk(const double tau, const int Ntau, const double delta, int Ndelta, const EArrayd& molefrac, const int i, const int NXi, const int j, const int NXj, const int k, const int NXk) const = 0;
123
124 // Derivatives from isochoric thermodynamics (all have the same signature whithin each block)
125 #define X(f) virtual double f(const double T, const EArrayd& rhovec) const = 0;
127 #undef X
128 #define X(f) virtual EArrayd f(const double T, const EArrayd& rhovec) const = 0;
130 #undef X
131 #define X(f) virtual EMatrixd f(const double T, const EArrayd& rhovec) const = 0;
133 #undef X
134 #define X(f) virtual std::tuple<double, Eigen::ArrayXd, Eigen::MatrixXd> f(const double T, const EArrayd& rhovec) const = 0;
136 #undef X
137 virtual Eigen::ArrayXd get_Psir_sigma_derivs(const double T, const EArrayd& rhovec, const EArrayd& v) const = 0;
138
139 double get_neff(const double, const double, const EArrayd&) const;
140
141 virtual EArray33d get_deriv_mat2(const double T, double rho, const EArrayd& z ) const = 0;
142
143 std::tuple<double, double> solve_pure_critical(const double T, const double rho, const std::optional<nlohmann::json>& = std::nullopt) const ;
144 EArray2 extrapolate_from_critical(const double Tc, const double rhoc, const double Tgiven, const std::optional<Eigen::ArrayXd>& molefracs = std::nullopt) const;
145 std::tuple<EArrayd, EMatrixd> get_pure_critical_conditions_Jacobian(const double T, const double rho, const std::optional<std::size_t>& alternative_pure_index, const std::optional<std::size_t>& alternative_length) const;
146
147 EArray2 pure_VLE_T(const double T, const double rhoL, const double rhoV, int maxiter, const std::optional<Eigen::ArrayXd>& molefracs = std::nullopt) const;
148 double dpsatdT_pure(const double T, const double rhoL, const double rhoV) const;
149
150 virtual std::tuple<EArrayd, EArrayd> get_drhovecdp_Tsat(const double T, const REArrayd& rhovecL, const REArrayd& rhovecV) const;
151 virtual std::tuple<EArrayd, EArrayd> get_drhovecdT_psat(const double T, const REArrayd& rhovecL, const REArrayd& rhovecV) const;
152 virtual double get_dpsat_dTsat_isopleth(const double T, const REArrayd& rhovecL, const REArrayd& rhovecV) const;
153 virtual nlohmann::json trace_VLE_isotherm_binary(const double T0, const EArrayd& rhovec0, const EArrayd& rhovecV0, const std::optional<TVLEOptions> & = std::nullopt) const;
154 virtual nlohmann::json trace_VLE_isobar_binary(const double p, const double T0, const EArrayd& rhovecL0, const EArrayd& rhovecV0, const std::optional<PVLEOptions> & = std::nullopt) const;
155 virtual std::tuple<VLE_return_code,EArrayd,EArrayd> mix_VLE_Tx(const double T, const REArrayd& rhovecL0, const REArrayd& rhovecV0, const REArrayd& xspec, const double atol, const double reltol, const double axtol, const double relxtol, const int maxiter) const;
156 virtual MixVLEReturn mix_VLE_Tp(const double T, const double pgiven, const REArrayd& rhovecL0, const REArrayd& rhovecV0, const std::optional<MixVLETpFlags> &flags = std::nullopt) const;
157 virtual std::tuple<VLE_return_code,double,EArrayd,EArrayd> mixture_VLE_px(const double p_spec, const REArrayd& xmolar_spec, const double T0, const REArrayd& rhovecL0, const REArrayd& rhovecV0, const std::optional<MixVLEpxFlags>& flags = std::nullopt) const;
158
159 std::tuple<VLLE::VLLE_return_code,EArrayd,EArrayd,EArrayd> mix_VLLE_T(const double T, const REArrayd& rhovecVinit, const REArrayd& rhovecL1init, const REArrayd& rhovecL2init, const double atol, const double reltol, const double axtol, const double relxtol, const int maxiter) const;
160 std::vector<nlohmann::json> find_VLLE_T_binary(const std::vector<nlohmann::json>& traces, const std::optional<VLLE::VLLEFinderOptions> options = std::nullopt) const;
161 std::vector<nlohmann::json> find_VLLE_p_binary(const std::vector<nlohmann::json>& traces, const std::optional<VLLE::VLLEFinderOptions> options = std::nullopt) const;
162 nlohmann::json trace_VLLE_binary(const double T, const REArrayd& rhovecV, const REArrayd& rhovecL1, const REArrayd& rhovecL2, const std::optional<VLLE::VLLETracerOptions> options) const;
163
164 virtual nlohmann::json trace_critical_arclength_binary(const double T0, const EArrayd& rhovec0, const std::optional<std::string>& = std::nullopt, const std::optional<TCABOptions> & = std::nullopt) const;
165 virtual EArrayd get_drhovec_dT_crit(const double T, const REArrayd& rhovec) const;
166 virtual double get_dp_dT_crit(const double T, const REArrayd& rhovec) const;
167 virtual EArray2 get_criticality_conditions(const double T, const REArrayd& rhovec) const;
168 virtual EigenData eigen_problem(const double T, const REArrayd& rhovec, const std::optional<REArrayd>& = std::nullopt) const;
169 virtual double get_minimum_eigenvalue_Psi_Hessian(const double T, const REArrayd& rhovec) const;
170
171 };
172
173 // Generic JSON-based interface where the model description is encoded as JSON
174 std::unique_ptr<AbstractModel> make_model(const nlohmann::json &, bool validate = true);
175
176 // Expose specialized factory functions for different models
177 // Mostly these are just adapter functions that prepare some
178 // JSON and pass it to the make_model function
179 // ....
180 std::unique_ptr<AbstractModel> make_multifluid_model(
181 const std::vector<std::string>& components,
182 const std::string& coolprop_root,
183 const std::string& BIPcollectionpath = {},
184 const nlohmann::json& flags = {},
185 const std::string& departurepath = {}
186 );
187
188 std::unique_ptr<AbstractModel> build_model_ptr(const nlohmann::json& json, bool validate = true);
189
191 nlohmann::json get_model_schema(const std::string& kind);
192
193
194 using ModelPointerFactoryFunction = std::function<std::unique_ptr<teqp::cppinterface::AbstractModel>(const nlohmann::json &j)>;
195
201 void add_model_pointer_factory_function(const std::string& key, ModelPointerFactoryFunction& func);
202
203 }
204}
std::vector< nlohmann::json > find_VLLE_T_binary(const std::vector< nlohmann::json > &traces, const std::optional< VLLE::VLLEFinderOptions > options=std::nullopt) const
virtual double get_AtaudeltaXi(const double tau, const int Ntau, const double delta, int Ndelta, const EArrayd &molefrac, const int i, const int NXi) const =0
std::tuple< VLLE::VLLE_return_code, EArrayd, EArrayd, EArrayd > mix_VLLE_T(const double T, const REArrayd &rhovecVinit, const REArrayd &rhovecL1init, const REArrayd &rhovecL2init, const double atol, const double reltol, const double axtol, const double relxtol, const int maxiter) const
virtual std::map< int, double > get_Bnvir(const int Nderiv, const double T, const EArrayd &z) const =0
virtual Eigen::ArrayXd get_Psir_sigma_derivs(const double T, const EArrayd &rhovec, const EArrayd &v) const =0
virtual double get_ATrhoXiXjXk(const double T, const int NT, const double rhomolar, int ND, const EArrayd &molefrac, const int i, const int NXi, const int j, const int NXj, const int k, const int NXk) const =0
double dpsatdT_pure(const double T, const double rhoL, const double rhoV) const
virtual std::tuple< VLE_return_code, double, EArrayd, EArrayd > mixture_VLE_px(const double p_spec, const REArrayd &xmolar_spec, const double T0, const REArrayd &rhovecL0, const REArrayd &rhovecV0, const std::optional< MixVLEpxFlags > &flags=std::nullopt) const
virtual double get_R(const EArrayd &) const =0
virtual double get_ATrhoXiXj(const double T, const int NT, const double rhomolar, int ND, const EArrayd &molefrac, const int i, const int NXi, const int j, const int NXj) const =0
virtual nlohmann::json trace_VLE_isobar_binary(const double p, const double T0, const EArrayd &rhovecL0, const EArrayd &rhovecV0, const std::optional< PVLEOptions > &=std::nullopt) const
virtual double get_dpsat_dTsat_isopleth(const double T, const REArrayd &rhovecL, const REArrayd &rhovecV) const
virtual double get_B2vir(const double T, const EArrayd &z) const =0
virtual double get_AtaudeltaXiXjXk(const double tau, const int Ntau, const double delta, int Ndelta, const EArrayd &molefrac, const int i, const int NXi, const int j, const int NXj, const int k, const int NXk) const =0
double get_neff(const double, const double, const EArrayd &) const
nlohmann::json trace_VLLE_binary(const double T, const REArrayd &rhovecV, const REArrayd &rhovecL1, const REArrayd &rhovecL2, const std::optional< VLLE::VLLETracerOptions > options) const
virtual EArray33d get_deriv_mat2(const double T, double rho, const EArrayd &z) const =0
virtual EArrayd get_drhovec_dT_crit(const double T, const REArrayd &rhovec) const
virtual double get_minimum_eigenvalue_Psi_Hessian(const double T, const REArrayd &rhovec) const
virtual EArray2 get_criticality_conditions(const double T, const REArrayd &rhovec) const
EArray2 extrapolate_from_critical(const double Tc, const double rhoc, const double Tgiven, const std::optional< Eigen::ArrayXd > &molefracs=std::nullopt) const
virtual MixVLEReturn mix_VLE_Tp(const double T, const double pgiven, const REArrayd &rhovecL0, const REArrayd &rhovecV0, const std::optional< MixVLETpFlags > &flags=std::nullopt) const
virtual nlohmann::json trace_VLE_isotherm_binary(const double T0, const EArrayd &rhovec0, const EArrayd &rhovecV0, const std::optional< TVLEOptions > &=std::nullopt) const
virtual double get_ATrhoXi(const double T, const int NT, const double rhomolar, int ND, const EArrayd &molefrac, const int i, const int NXi) const =0
virtual std::tuple< VLE_return_code, EArrayd, EArrayd > mix_VLE_Tx(const double T, const REArrayd &rhovecL0, const REArrayd &rhovecV0, const REArrayd &xspec, const double atol, const double reltol, const double axtol, const double relxtol, const int maxiter) const
double R(const EArrayd &x) const
Definition teqpcpp.hpp:95
virtual EigenData eigen_problem(const double T, const REArrayd &rhovec, const std::optional< REArrayd > &=std::nullopt) const
EArray2 pure_VLE_T(const double T, const double rhoL, const double rhoV, int maxiter, const std::optional< Eigen::ArrayXd > &molefracs=std::nullopt) const
virtual double get_AtaudeltaXiXj(const double tau, const int Ntau, const double delta, int Ndelta, const EArrayd &molefrac, const int i, const int NXi, const int j, const int NXj) const =0
std::tuple< double, double > solve_pure_critical(const double T, const double rho, const std::optional< nlohmann::json > &=std::nullopt) const
virtual double get_dmBnvirdTm(const int Nderiv, const int NTderiv, const double T, const EArrayd &z) const =0
virtual const std::type_index & get_type_index() const =0
virtual double get_dp_dT_crit(const double T, const REArrayd &rhovec) const
virtual std::tuple< EArrayd, EArrayd > get_drhovecdp_Tsat(const double T, const REArrayd &rhovecL, const REArrayd &rhovecV) const
std::tuple< EArrayd, EMatrixd > get_pure_critical_conditions_Jacobian(const double T, const double rho, const std::optional< std::size_t > &alternative_pure_index, const std::optional< std::size_t > &alternative_length) const
virtual double get_Arxy(const int, const int, const double, const double, const EArrayd &) const =0
std::vector< nlohmann::json > find_VLLE_p_binary(const std::vector< nlohmann::json > &traces, const std::optional< VLLE::VLLEFinderOptions > options=std::nullopt) const
virtual nlohmann::json trace_critical_arclength_binary(const double T0, const EArrayd &rhovec0, const std::optional< std::string > &=std::nullopt, const std::optional< TCABOptions > &=std::nullopt) const
virtual double get_B12vir(const double T, const EArrayd &z) const =0
virtual std::tuple< EArrayd, EArrayd > get_drhovecdT_psat(const double T, const REArrayd &rhovecL, const REArrayd &rhovecV) const
std::unique_ptr< AbstractModel > make_model(const nlohmann::json &, bool validate=true)
std::function< std::unique_ptr< teqp::cppinterface::AbstractModel >(const nlohmann::json &j)> ModelPointerFactoryFunction
Definition teqpcpp.hpp:194
std::unique_ptr< AbstractModel > make_multifluid_model(const std::vector< std::string > &components, const std::string &coolprop_root, const std::string &BIPcollectionpath={}, const nlohmann::json &flags={}, const std::string &departurepath={})
void add_model_pointer_factory_function(const std::string &key, ModelPointerFactoryFunction &func)
This function allows you to inject your own model factory function into the set of factory functions ...
std::unique_ptr< AbstractModel > build_model_ptr(const nlohmann::json &json, bool validate=true)
nlohmann::json get_model_schema(const std::string &kind)
Return the schema for the given model kind.
#define AR0N_args
Definition teqpcpp.hpp:39
#define ISOCHORIC_matrix_args
Definition teqpcpp.hpp:63
#define ISOCHORIC_double_args
Definition teqpcpp.hpp:49
Eigen::Ref< const Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic > > REMatrixd
Definition teqpcpp.hpp:20
Eigen::Ref< const EArrayd > REArrayd
Definition teqpcpp.hpp:18
Eigen::Array< double, 3, 3 > EArray33d
Definition teqpcpp.hpp:17
Eigen::Array< double, 2, 1 > EArray2
Definition teqpcpp.hpp:15
Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic > EMatrixd
Definition teqpcpp.hpp:19
#define ISOCHORIC_multimatrix_args
Definition teqpcpp.hpp:67
#define ISOCHORIC_array_args
Definition teqpcpp.hpp:54
Eigen::ArrayX< double > EArrayd
Definition teqpcpp.hpp:16
#define ARXY_args
Definition teqpcpp.hpp:22