6#include "nlohmann/json.hpp"
17template<
typename X>
auto POW2(
X x) {
return x * x; };
18template<
typename X>
auto POW3(
X x) {
return x *
POW2(x); };
26template<
typename BType,
typename TType,
typename RhoType,
typename VecType>
27inline auto get_DeltaAB_pure(radial_dist dist,
double epsABi,
double betaABi, BType b_cubic, TType RT, RhoType rhomolar,
const VecType& ) {
29 using eta_type = std::common_type_t<
decltype(rhomolar),
decltype(b_cubic)>;
35 case radial_dist::CS: {
37 eta = (rhomolar / 4.0) * b_cubic;
38 g_vm_ref = (2.0 - eta) / (2.0 *
POW3(1.0 - eta));
41 case radial_dist::KG: {
43 eta = (rhomolar / 4.0) * b_cubic;
44 g_vm_ref = 1.0 / (1.0 - 1.9 * eta);
52 throw std::invalid_argument(
"Bad radial_dist");
57 auto DeltaAiBj =
forceeval(g_vm_ref*(exp(epsABi/RT) - 1.0)*b_cubic* betaABi);
67template<
typename BType,
typename TType,
typename RhoType,
typename VecType>
68inline auto XA_calc_pure(
int N_sites, association_classes scheme, radial_dist dist,
double epsABi,
double betaABi,
const BType b_cubic,
const TType RT,
const RhoType rhomolar,
const VecType& molefrac) {
72 using result_type = std::common_type_t<
decltype(RT),
decltype(rhomolar),
decltype(molefrac[0])>;
73 Eigen::Array<result_type, Eigen::Dynamic, Eigen::Dynamic> XA;
74 XA.resize(N_sites, 1);
78 auto DeltaAiBj =
get_DeltaAB_pure(dist, epsABi, betaABi, b_cubic, RT, rhomolar, molefrac);
80 if (scheme == association_classes::a1A) {
82 XA(0, 0) =
forceeval((-1.0 + sqrt(1.0 + 4.0 * rhomolar * DeltaAiBj)) / (2.0 * rhomolar * DeltaAiBj));
84 else if (scheme == association_classes::a2B) {
86 XA(0, 0) =
forceeval((-1.0 + sqrt(1.0 + 4.0 * rhomolar * DeltaAiBj)) / (2.0 * rhomolar * DeltaAiBj));
89 else if (scheme == association_classes::a3B) {
91 XA(0, 0) =
forceeval((-(1.0 - rhomolar * DeltaAiBj) + sqrt(
POW2(1.0 + rhomolar * DeltaAiBj) + 4.0 * rhomolar * DeltaAiBj)) / (4.0 * rhomolar * DeltaAiBj));
93 XA(2, 0) = 2.0*XA(0, 0) - 1.0;
95 else if (scheme == association_classes::a4C) {
97 XA(0, 0) =
forceeval((-1.0 + sqrt(1.0 + 8.0 * rhomolar * DeltaAiBj)) / (4.0 * rhomolar * DeltaAiBj));
102 else if (scheme == association_classes::not_associating) {
109 throw std::invalid_argument(
"Bad scheme");
119 throw std::invalid_argument(
"bad cubic flag:" + s);
126 const std::valarray<double> a0, bi, c1, Tc;
127 double delta_1, delta_2;
129 const std::optional<std::vector<std::vector<double>>> kmat;
131 CPACubic(
cubic_flag flag,
const std::valarray<double> &a0,
const std::valarray<double> &bi,
const std::valarray<double> &c1,
const std::valarray<double> &Tc,
const double R_gas,
const std::optional<std::vector<std::vector<double>>> & kmat = std::nullopt) : a0(a0), bi(bi), c1(c1), Tc(Tc), R_gas(R_gas), kmat(kmat) {
134 { delta_1 = 1 + sqrt(2); delta_2 = 1 - sqrt(2);
break; }
136 { delta_1 = 1; delta_2 = 0;
break; }
138 throw std::invalid_argument(
"Bad cubic flag");
142 std::size_t
size()
const {
return a0.size(); }
144 template<
typename VecType>
145 auto R(
const VecType& )
const {
return R_gas; }
147 template<
typename TType>
149 return forceeval(a0[i] *
POW2(1.0 + c1[i]*(1.0 - sqrt(T / Tc[i]))));
152 template<
typename TType,
typename VecType>
153 auto get_ab(
const TType T,
const VecType& molefrac)
const {
154 using return_type = std::common_type_t<
decltype(T),
decltype(molefrac[0])>;
155 return_type asummer = 0.0, bsummer = 0.0;
156 for (
auto i = 0U; i < molefrac.size(); ++i) {
157 bsummer += molefrac[i] * bi[i];
159 for (
auto j = 0U; j < molefrac.size(); ++j) {
161 double kij = (kmat) ? kmat.value()[i][j] : 0.0;
162 auto a_ij = (1.0 - kij) * sqrt(ai * aj);
163 asummer += molefrac[i] * molefrac[j] * a_ij;
166 return std::make_tuple(asummer, bsummer);
169 template<
typename TType,
typename RhoType,
typename VecType>
170 auto alphar(
const TType T,
const RhoType rhomolar,
const VecType& molefrac)
const {
171 auto [a_cubic, b_cubic] =
get_ab(T, molefrac);
173 -log(1.0 - b_cubic * rhomolar)
174 -a_cubic/R_gas/T*log((delta_1*b_cubic*rhomolar + 1.0) / (delta_2*b_cubic*rhomolar + 1.0)) / b_cubic / (delta_1 - delta_2)
183 const std::vector<association_classes> classes;
184 const radial_dist dist;
185 const std::valarray<double> epsABi, betaABi, bi;
186 const std::vector<int> N_sites;
189 auto get_N_sites(
const std::vector<association_classes> &the_classes) {
190 std::vector<int> N_sites_out;
191 auto get_N = [](
auto cl) {
193 case association_classes::a1A:
return 1;
194 case association_classes::a2B:
return 2;
195 case association_classes::a3B:
return 3;
196 case association_classes::a4C:
return 4;
197 default:
throw std::invalid_argument(
"Bad association class");
200 for (
auto cl : the_classes) {
201 N_sites_out.push_back(get_N(cl));
207 CPAAssociation(
const std::vector<association_classes>& classes,
const radial_dist dist,
const std::valarray<double> &epsABi,
const std::valarray<double> &betaABi,
const std::valarray<double> &bi,
double R_gas)
208 : classes(classes), dist(dist), epsABi(epsABi), betaABi(betaABi), bi(bi), N_sites(get_N_sites(classes)), R_gas(R_gas) {};
210 nlohmann::json
get_assoc_calcs(
double T,
double rhomolar,
const Eigen::ArrayXd& mole_fractions)
const{
212 auto fromArrayX = [](
const Eigen::ArrayXd &x){std::valarray<double>n(x.size());
for (
auto i = 0U; i < n.size(); ++i){ n[i] = x[i];}
return n;};
215 auto b_cubic = (Eigen::Map<const Eigen::ArrayXd>(&bi[0], bi.size())*mole_fractions).sum();
219 auto XA =
XA_calc_pure(N_sites[0], classes[0], dist, epsABi[0], betaABi[0], b_cubic, RT, rhomolar, mole_fractions);
223 {
"X_A", fromArrayX(XA)},
224 {
"note",
"X_A is the fraction of non-bonded sites for each site type"}
228 template<
typename TType,
typename RhoType,
typename VecType>
229 auto alphar(
const TType& T,
const RhoType& rhomolar,
const VecType& molefrac)
const {
231 auto b_cubic = (Eigen::Map<const Eigen::ArrayXd>(&bi[0], bi.size())*molefrac).sum();
235 auto XA =
XA_calc_pure(N_sites[0], classes[0], dist, epsABi[0], betaABi[0], b_cubic, RT, rhomolar, molefrac);
237 using return_type = std::common_type_t<
decltype(T),
decltype(rhomolar),
decltype(molefrac[0])>;
238 return_type alpha_r_asso = 0.0;
240 for (
auto xi : molefrac){
241 auto XAi = XA.col(i);
242 alpha_r_asso +=
forceeval(xi * (log(XAi) - XAi / 2).sum());
243 alpha_r_asso += xi*
static_cast<double>(N_sites[i])/2;
250template <
typename Cubic,
typename Assoc>
256 template<
class VecType>
257 auto R(
const VecType& molefrac)
const {
258 return cubic.R(molefrac);
266 template<
typename TType,
typename RhoType,
typename VecType>
267 auto alphar(
const TType& T,
const RhoType& rhomolar,
const VecType& molefrac)
const {
268 if (
static_cast<std::size_t
>(molefrac.size()) !=
cubic.size()){
273 auto alpha_r_cubic =
cubic.alphar(T, rhomolar, molefrac);
276 auto alpha_r_assoc =
assoc.alphar(T, rhomolar, molefrac);
278 return forceeval(alpha_r_cubic + alpha_r_assoc);
283 using vartype = std::variant<CPAAssociation, association::Association>;
288 template<
typename TType,
typename RhoType,
typename MoleFracsType>
289 auto alphar(
const TType& T,
const RhoType& rhomolar,
const MoleFracsType& molefracs)
const{
290 return std::visit([&](
auto& h){
return h.alphar(T, rhomolar, molefracs); },
holder);
293 auto get_assoc_calcs(
double T,
double rhomolar,
const Eigen::ArrayXd& mole_fractions)
const {
294 return std::visit([&](
auto& h){
return h.get_assoc_calcs(T, rhomolar, mole_fractions); },
holder);
302 auto build_cubic = [](
const auto& j) {
303 auto N = j[
"pures"].size();
304 std::valarray<double> a0i(N), bi(N), c1(N), Tc(N);
305 std::vector<std::vector<double>> kmat;
306 if (j.contains(
"kmat")){
308 std::string kmaterr =
"The kmat is the wrong size. It should be square with dimension " + std::to_string(N);
309 if (kmat.size() != N){
313 for (
auto& krow: kmat){
314 if(krow.size() != N){
321 kmat.resize(N);
for (
auto i = 0U; i < N; ++i){ kmat[i].resize(N);
for (
auto j = 0U; j < N; ++j){kmat[i][j] = 0.0;} }
325 for (
auto p : j[
"pures"]) {
326 a0i[i] = p[
"a0i / Pa m^6/mol^2"];
327 bi[i] = p[
"bi / m^3/mol"];
336 auto N = j[
"pures"].size();
337 if (N == 1 && j.at(
"pures").at(0).contains(
"class") ){
340 std::vector<association_classes> classes;
341 radial_dist dist = get_radial_dist(j.at(
"radial_dist"));
342 std::valarray<double> epsABi(N), betaABi(N), bi(N);
344 for (
auto p : j.at(
"pures")) {
345 epsABi[i] = p.at(
"epsABi / J/mol");
346 betaABi[i] = p.at(
"betaABi");
347 bi[i] = p.at(
"bi / m^3/mol");
348 classes.push_back(get_association_classes(p.at(
"class")));
355 Eigen::ArrayXd b_m3mol(N), beta(N), epsilon_Jmol(N);
357 opt.
radial_dist = get_radial_dist(j.at(
"radial_dist"));
358 if (j.contains(
"options")){
359 opt = j.at(
"options");
362 std::vector<std::vector<std::string>> molecule_sites;
364 std::set<std::string> unique_site_types;
365 for (
auto p : j[
"pures"]) {
366 epsilon_Jmol[i] = p.at(
"epsABi / J/mol");
367 beta[i] = p.at(
"betaABi");
368 b_m3mol[i] = p.at(
"bi / m^3/mol");
369 molecule_sites.push_back(p.at(
"sites"));
370 for (
auto & s : molecule_sites.back()){
371 unique_site_types.insert(s);
375 if (j.contains(
"options") && j.at(
"options").contains(
"interaction_partners")){
378 if (unique_site_types.count(k) == 0){
381 for (
auto& partner : partners){
382 if (unique_site_types.count(partner) == 0){
390 for (
auto& site1 : unique_site_types){
391 std::vector<std::string> partners;
392 for (
auto& site2: unique_site_types){
394 partners.push_back(site2);
404 return CPAEOS(build_cubic(j), build_assoc_pure( j));
CPAAssociation(const std::vector< association_classes > &classes, const radial_dist dist, const std::valarray< double > &epsABi, const std::valarray< double > &betaABi, const std::valarray< double > &bi, double R_gas)
auto alphar(const TType &T, const RhoType &rhomolar, const VecType &molefrac) const
nlohmann::json get_assoc_calcs(double T, double rhomolar, const Eigen::ArrayXd &mole_fractions) const
CPACubic(cubic_flag flag, const std::valarray< double > &a0, const std::valarray< double > &bi, const std::valarray< double > &c1, const std::valarray< double > &Tc, const double R_gas, const std::optional< std::vector< std::vector< double > > > &kmat=std::nullopt)
auto alphar(const TType T, const RhoType rhomolar, const VecType &molefrac) const
auto get_ai(TType T, int i) const
auto R(const VecType &) const
auto get_ab(const TType T, const VecType &molefrac) const
auto alphar(const TType &T, const RhoType &rhomolar, const VecType &molefrac) const
auto R(const VecType &molefrac) const
CPAEOS(Cubic &&cubic, Assoc &&assoc)
auto get_DeltaAB_pure(radial_dist dist, double epsABi, double betaABi, BType b_cubic, TType RT, RhoType rhomolar, const VecType &)
Function that calculates the association binding strength between site A of molecule i and site B on ...
auto XA_calc_pure(int N_sites, association_classes scheme, radial_dist dist, double epsABi, double betaABi, const BType b_cubic, const TType RT, const RhoType rhomolar, const VecType &molefrac)
auto CPAfactory(const nlohmann::json &j)
auto get_cubic_flag(const std::string &s)
auto get_association_classes(const std::string &s)
auto get_radial_dist(const std::string &s)
std::variant< CPAAssociation, association::Association > vartype
AssociationVariantWrapper(const vartype &holder)
auto get_assoc_calcs(double T, double rhomolar, const Eigen::ArrayXd &mole_fractions) const
auto alphar(const TType &T, const RhoType &rhomolar, const MoleFracsType &molefracs) const
std::map< std::string, std::vector< std::string > > interaction_partners
association::radial_dist radial_dist