teqp 0.19.1
Loading...
Searching...
No Matches
pcsaft.hpp
Go to the documentation of this file.
1/***
2
3 \brief This file contains the contributions that can be composed together to form SAFT models
4
5*/
6
7#pragma once
8
9#include "nlohmann/json.hpp"
10#include "teqp/types.hpp"
11#include "teqp/exceptions.hpp"
12#include "teqp/constants.hpp"
13#include "teqp/json_tools.hpp"
15#include <optional>
16
17namespace teqp {
18namespace PCSAFT {
19
21struct SAFTCoeffs {
22 std::string name;
23 double m = -1,
26 std::string BibTeXKey;
27 double mustar2 = 0,
28 nmu = 0,
29 Qstar2 = 0,
30 nQ = 0;
31};
32
35 std::map<std::string, SAFTCoeffs> coeffs;
36public:
38 insert_normal_fluid("Methane", 1.0000, 3.7039, 150.03, "Gross-IECR-2001");
39 insert_normal_fluid("Ethane", 1.6069, 3.5206, 191.42, "Gross-IECR-2001");
40 insert_normal_fluid("Propane", 2.0020, 3.6184, 208.11, "Gross-IECR-2001");
41 }
42 void insert_normal_fluid(const std::string& name, double m, const double sigma_Angstrom, const double epsilon_over_k, const std::string& BibTeXKey) {
43 SAFTCoeffs coeff;
44 coeff.name = name;
45 coeff.m = m;
46 coeff.sigma_Angstrom = sigma_Angstrom;
47 coeff.epsilon_over_k = epsilon_over_k;
48 coeff.BibTeXKey = BibTeXKey;
49 coeffs.insert(std::pair<std::string, SAFTCoeffs>(name, coeff));
50 }
51 const auto& get_normal_fluid(const std::string& name) {
52 auto it = coeffs.find(name);
53 if (it != coeffs.end()) {
54 return it->second;
55 }
56 else {
57 throw std::invalid_argument("Bad name:" + name);
58 }
59 }
60 auto get_coeffs(const std::vector<std::string>& names){
61 std::vector<SAFTCoeffs> c;
62 for (auto n : names){
63 c.push_back(get_normal_fluid(n));
64 }
65 return c;
66 }
67};
68
72template <typename Eta, typename Mbar>
73auto C1(const Eta& eta, Mbar mbar) {
74 return forceeval(1.0 / (1.0
75 + mbar * (8.0 * eta - 2.0 * eta * eta) / pow(1.0 - eta, 4)
76 + (1.0 - mbar) * (20.0 * eta - 27.0 * eta * eta + 12.0 * pow(eta, 3) - 2.0 * pow(eta, 4)) / pow((1.0 - eta) * (2.0 - eta), 2)));
77}
79template <typename Eta, typename Mbar>
80auto C2(const Eta& eta, Mbar mbar) {
81 return forceeval(-pow(C1(eta, mbar), 2) * (
82 mbar * (-4.0 * eta * eta + 20.0 * eta + 8.0) / pow(1.0 - eta, 5)
83 + (1.0 - mbar) * (2.0 * eta * eta * eta + 12.0 * eta * eta - 48.0 * eta + 40.0) / pow((1.0 - eta) * (2.0 - eta), 3)
84 ));
85}
87template<typename TYPE>
88auto get_a(TYPE mbar) {
89 static Eigen::ArrayXd a_0 = (Eigen::ArrayXd(7) << 0.9105631445, 0.6361281449, 2.6861347891, -26.547362491, 97.759208784, -159.59154087, 91.297774084).finished();
90 static Eigen::ArrayXd a_1 = (Eigen::ArrayXd(7) << -0.3084016918, 0.1860531159, -2.5030047259, 21.419793629, -65.255885330, 83.318680481, -33.746922930).finished();
91 static Eigen::ArrayXd a_2 = (Eigen::ArrayXd(7) << -0.0906148351, 0.4527842806, 0.5962700728, -1.7241829131, -4.1302112531, 13.776631870, -8.6728470368).finished();
92 return forceeval(a_0.cast<TYPE>().array() + ((mbar - 1.0) / mbar) * a_1.cast<TYPE>().array() + ((mbar - 1.0) / mbar * (mbar - 2.0) / mbar) * a_2.cast<TYPE>().array()).eval();
93}
95template<typename TYPE>
96auto get_b(TYPE mbar) {
97 // See https://stackoverflow.com/a/35170514/1360263
98 static Eigen::ArrayXd b_0 = (Eigen::ArrayXd(7) << 0.7240946941, 2.2382791861, -4.0025849485, -21.003576815, 26.855641363, 206.55133841, -355.60235612).finished();
99 static Eigen::ArrayXd b_1 = (Eigen::ArrayXd(7) << -0.5755498075, 0.6995095521, 3.8925673390, -17.215471648, 192.67226447, -161.82646165, -165.20769346).finished();
100 static Eigen::ArrayXd b_2 = (Eigen::ArrayXd(7) << 0.0976883116, -0.2557574982, -9.1558561530, 20.642075974, -38.804430052, 93.626774077, -29.666905585).finished();
101 return forceeval(b_0.cast<TYPE>().array() + (mbar - 1.0) / mbar * b_1.cast<TYPE>().array() + (mbar - 1.0) / mbar * (mbar - 2.0) / mbar * b_2.cast<TYPE>().array()).eval();
102}
104template<typename VecType>
105auto get_alphar_hs(const VecType& zeta) {
106 // The limit of alphar_hs in the case of density going to zero is still zero,
107 // but the way it goes to zero is subtle
108 if (getbaseval(zeta[3]) == 0){
109 return forceeval(4.0*zeta[3]);
110 }
111 auto Upsilon = 1.0 - zeta[3];
112 return forceeval(1.0 / zeta[0] * (3.0 * zeta[1] * zeta[2] / Upsilon
113 + zeta[2] * zeta[2] * zeta[2] / zeta[3] / Upsilon / Upsilon
114 + (zeta[2] * zeta[2] * zeta[2] / (zeta[3] * zeta[3]) - zeta[0]) * log(1.0 - zeta[3])
115 ));
116}
117
119template<typename zVecType, typename dVecType>
120auto gij_HS(const zVecType& zeta, const dVecType& d,
121 std::size_t i, std::size_t j) {
122 auto Upsilon = 1.0 - zeta[3];
123 return forceeval(1.0 / (Upsilon)+d[i] * d[j] / (d[i] + d[j]) * 3.0 * zeta[2] / pow(Upsilon, 2)
124 + pow(d[i] * d[j] / (d[i] + d[j]), 2) * 2.0 * pow(zeta[2], 2) / pow(Upsilon, 3));
125}
127template <typename Eta, typename MbarType>
128auto get_I1(const Eta& eta, MbarType mbar) {
129 auto avec = get_a(mbar);
130 Eta summer_I1 = 0.0, summer_etadI1deta = 0.0;
131 for (std::size_t i = 0; i < 7; ++i) {
132 auto increment = avec(i) * powi(eta, static_cast<int>(i));
133 summer_I1 = summer_I1 + increment;
134 summer_etadI1deta = summer_etadI1deta + increment * (i + 1.0);
135 }
136 return std::make_tuple(forceeval(summer_I1), forceeval(summer_etadI1deta));
137}
139template <typename Eta, typename MbarType>
140auto get_I2(const Eta& eta, MbarType mbar) {
141 auto bvec = get_b(mbar);
142 Eta summer_I2 = 0.0 * eta, summer_etadI2deta = 0.0 * eta;
143 for (std::size_t i = 0; i < 7; ++i) {
144 auto increment = bvec(i) * powi(eta, static_cast<int>(i));
145 summer_I2 = summer_I2 + increment;
146 summer_etadI2deta = summer_etadI2deta + increment * (i + 1.0);
147 }
148 return std::make_tuple(forceeval(summer_I2), forceeval(summer_etadI2deta));
149}
150
154template<typename VecType1, typename NType>
155auto powvec(const VecType1& v1, NType n) {
156 auto o = v1;
157 for (auto i = 0; i < v1.size(); ++i) {
158 o[i] = pow(v1[i], n);
159 }
160 return o;
161}
162
166template<typename VecType1, typename VecType2, typename VecType3>
167auto sumproduct(const VecType1& v1, const VecType2& v2, const VecType3& v3) {
168 using ResultType = typename std::common_type_t<decltype(v1[0]), decltype(v2[0]), decltype(v3[0])>;
169 return forceeval((v1.template cast<ResultType>().array() * v2.template cast<ResultType>().array() * v3.template cast<ResultType>().array()).sum());
170}
171
173template<typename NumType, typename ProductType>
174class SAFTCalc {
175public:
176 // Just temperature dependent things
177 Eigen::ArrayX<NumType> d;
178
179 // These things also have composition dependence
182};
183
184/***
185 * \brief This class provides the evaluation of the hard chain contribution from classic PC-SAFT
186 */
188
189protected:
190 const Eigen::ArrayX<double> m,
194 const Eigen::ArrayXXd kmat;
195
196public:
197 PCSAFTHardChainContribution(const Eigen::ArrayX<double> &m, const Eigen::ArrayX<double> &mminus1, const Eigen::ArrayX<double> &sigma_Angstrom, const Eigen::ArrayX<double> &epsilon_over_k, const Eigen::ArrayXXd &kmat)
199
201
202 template<typename TTYPE, typename RhoType, typename VecType>
203 auto eval(const TTYPE& T, const RhoType& rhomolar, const VecType& mole_fractions) const {
204
205 Eigen::Index N = m.size();
206
207 if (mole_fractions.size() != N) {
208 throw std::invalid_argument("Length of mole_fractions (" + std::to_string(mole_fractions.size()) + ") is not the length of components (" + std::to_string(N) + ")");
209 }
210
211 using TRHOType = std::common_type_t<std::decay_t<TTYPE>, std::decay_t<RhoType>, std::decay_t<decltype(mole_fractions[0])>, std::decay_t<decltype(m[0])>>;
212
214 c.m2_epsilon_sigma3_bar = static_cast<TRHOType>(0.0);
215 c.m2_epsilon2_sigma3_bar = static_cast<TRHOType>(0.0);
216 c.d.resize(N);
217 for (auto i = 0L; i < N; ++i) {
218 c.d[i] = sigma_Angstrom[i]*(1.0 - 0.12 * exp(-3.0*epsilon_over_k[i]/T)); // [A]
219 for (auto j = 0; j < N; ++j) {
220 // Eq. A.5
221 auto sigma_ij = 0.5 * sigma_Angstrom[i] + 0.5 * sigma_Angstrom[j];
222 auto eij_over_k = sqrt(epsilon_over_k[i] * epsilon_over_k[j]) * (1.0 - kmat(i,j));
223 c.m2_epsilon_sigma3_bar = c.m2_epsilon_sigma3_bar + mole_fractions[i] * mole_fractions[j] * m[i] * m[j] * eij_over_k / T * pow(sigma_ij, 3);
224 c.m2_epsilon2_sigma3_bar = c.m2_epsilon2_sigma3_bar + mole_fractions[i] * mole_fractions[j] * m[i] * m[j] * pow(eij_over_k / T, 2) * pow(sigma_ij, 3);
225 }
226 }
227 auto mbar = (mole_fractions.template cast<TRHOType>().array()*m.template cast<TRHOType>().array()).sum();
228
230 RhoType rho_A3 = rhomolar * N_A * 1e-30; //[molecules (not moles)/A^3]
231
232 constexpr double MY_PI = EIGEN_PI;
233 double pi6 = (MY_PI / 6.0);
234
236 using ta = std::common_type_t<decltype(pi6), decltype(m[0]), decltype(c.d[0]), decltype(rho_A3)>;
237 std::vector<ta> zeta(4);
238 for (std::size_t n = 0; n < 4; ++n) {
239 // Eqn A.8
240 auto dn = pow(c.d, static_cast<int>(n));
241 TRHOType xmdn = forceeval((mole_fractions.template cast<TRHOType>().array()*m.template cast<TRHOType>().array()*dn.template cast<TRHOType>().array()).sum());
242 zeta[n] = forceeval(pi6*rho_A3*xmdn);
243 }
244
246 auto eta = zeta[3];
247
248 auto [I1, etadI1deta] = get_I1(eta, mbar);
249 auto [I2, etadI2deta] = get_I2(eta, mbar);
250
251 // Hard chain contribution from G&S
252 using tt = std::common_type_t<decltype(zeta[0]), decltype(c.d[0])>;
253 Eigen::ArrayX<tt> lngii_hs(mole_fractions.size());
254 for (auto i = 0; i < lngii_hs.size(); ++i) {
255 lngii_hs[i] = log(gij_HS(zeta, c.d, i, i));
256 }
257 auto alphar_hc = forceeval(mbar * get_alphar_hs(zeta) - sumproduct(mole_fractions, mminus1, lngii_hs)); // Eq. A.4
258
259 // Dispersive contribution
260 auto alphar_disp = forceeval(-2 * MY_PI * rho_A3 * I1 * c.m2_epsilon_sigma3_bar - MY_PI * rho_A3 * mbar * C1(eta, mbar) * I2 * c.m2_epsilon2_sigma3_bar);
261
262 using eta_t = decltype(eta);
263 using hc_t = decltype(alphar_hc);
264 using disp_t = decltype(alphar_disp);
265 struct PCSAFTHardChainContributionTerms{
266 eta_t eta;
267 hc_t alphar_hc;
268 disp_t alphar_disp;
269 };
270 return PCSAFTHardChainContributionTerms{forceeval(eta), forceeval(alphar_hc), forceeval(alphar_disp)};
271 }
272};
273
281public:
284protected:
285 Eigen::ArrayX<double> m,
289 std::vector<std::string> names, bibtex;
290 Eigen::ArrayXXd kmat;
291
293 std::optional<PCSAFTDipolarContribution> dipolar; // Can be present or not
294 std::optional<PCSAFTQuadrupolarContribution> quadrupolar; // Can be present or not
295
296 void check_kmat(Eigen::Index N) {
297 if (kmat.cols() != kmat.rows()) {
298 throw teqp::InvalidArgument("kmat rows and columns are not identical");
299 }
300 if (kmat.cols() == 0) {
301 kmat.resize(N, N); kmat.setZero();
302 }
303 else if (kmat.cols() != N) {
304 throw teqp::InvalidArgument("kmat needs to be a square matrix the same size as the number of components");
305 }
306 };
307 auto get_coeffs_from_names(const std::vector<std::string> &the_names){
308 PCSAFTLibrary library;
309 return library.get_coeffs(the_names);
310 }
311 auto build_hardchain(const std::vector<SAFTCoeffs> &coeffs){
312 check_kmat(coeffs.size());
313
314 m.resize(coeffs.size());
315 mminus1.resize(coeffs.size());
316 sigma_Angstrom.resize(coeffs.size());
317 epsilon_over_k.resize(coeffs.size());
318 names.resize(coeffs.size());
319 bibtex.resize(coeffs.size());
320 auto i = 0;
321 for (const auto &coeff : coeffs) {
322 m[i] = coeff.m;
323 mminus1[i] = m[i] - 1;
324 sigma_Angstrom[i] = coeff.sigma_Angstrom;
325 epsilon_over_k[i] = coeff.epsilon_over_k;
326 names[i] = coeff.name;
327 bibtex[i] = coeff.BibTeXKey;
328 i++;
329 }
331 }
332 auto extract_names(const std::vector<SAFTCoeffs> &coeffs){
333 std::vector<std::string> names_;
334 for (const auto& c: coeffs){
335 names_.push_back(c.name);
336 }
337 return names_;
338 }
339 auto build_dipolar(const std::vector<SAFTCoeffs> &coeffs) -> std::optional<PCSAFTDipolarContribution>{
340 Eigen::ArrayXd mustar2(coeffs.size()), nmu(coeffs.size());
341 auto i = 0;
342 for (const auto &coeff : coeffs) {
343 mustar2[i] = coeff.mustar2;
344 nmu[i] = coeff.nmu;
345 i++;
346 }
347 if ((mustar2*nmu).cwiseAbs().sum() == 0){
348 return std::nullopt; // No dipolar contribution is present
349 }
350 // The dispersive and hard chain initialization has already happened at this point
352 }
353 auto build_quadrupolar(const std::vector<SAFTCoeffs> &coeffs) -> std::optional<PCSAFTQuadrupolarContribution>{
354 // The dispersive and hard chain initialization has already happened at this point
355 Eigen::ArrayXd Qstar2(coeffs.size()), nQ(coeffs.size());
356 auto i = 0;
357 for (const auto &coeff : coeffs) {
358 Qstar2[i] = coeff.Qstar2;
359 nQ[i] = coeff.nQ;
360 i++;
361 }
362 if ((Qstar2*nQ).cwiseAbs().sum() == 0){
363 return std::nullopt; // No quadrupolar contribution is present
364 }
366 }
367public:
368 PCSAFTMixture(const std::vector<std::string> &names, const Eigen::ArrayXXd& kmat = {}) : PCSAFTMixture(get_coeffs_from_names(names), kmat){};
369 PCSAFTMixture(const std::vector<SAFTCoeffs> &coeffs, const Eigen::ArrayXXd &kmat = {}) : names(extract_names(coeffs)), kmat(kmat), hardchain(build_hardchain(coeffs)), dipolar(build_dipolar(coeffs)), quadrupolar(build_quadrupolar(coeffs)) {};
370
371// PCSAFTMixture( const PCSAFTMixture& ) = delete; // non construction-copyable
372 PCSAFTMixture& operator=( const PCSAFTMixture& ) = delete; // non copyable
373
374 auto get_m() const { return m; }
375 auto get_sigma_Angstrom() const { return sigma_Angstrom; }
376 auto get_epsilon_over_k_K() const { return epsilon_over_k; }
377 auto get_kmat() const { return kmat; }
378 auto get_names() const { return names;}
379 auto get_BibTeXKeys() const { return bibtex;}
380
381 auto print_info() {
382 std::string s = std::string("i m sigma / A e/kB / K \n ++++++++++++++") + "\n";
383 for (auto i = 0; i < m.size(); ++i) {
384 s += std::to_string(i) + " " + std::to_string(m[i]) + " " + std::to_string(sigma_Angstrom[i]) + " " + std::to_string(epsilon_over_k[i]) + "\n";
385 }
386 return s;
387 }
388
389 template<typename VecType>
390 double max_rhoN(const double T, const VecType& mole_fractions) const {
391 auto N = mole_fractions.size();
392 Eigen::ArrayX<double> d(N);
393 for (auto i = 0; i < N; ++i) {
394 d[i] = sigma_Angstrom[i] * (1.0 - 0.12 * exp(-3.0 * epsilon_over_k[i] / T));
395 }
396 return 6 * 0.74 / EIGEN_PI / (mole_fractions*m*powvec(d, 3)).sum()*1e30; // particles/m^3
397 }
398
399 template<class VecType>
400 auto R(const VecType& molefrac) const {
401 return get_R_gas<decltype(molefrac[0])>();
402 }
403
404 template<typename TTYPE, typename RhoType, typename VecType>
405 auto alphar(const TTYPE& T, const RhoType& rhomolar, const VecType& mole_fractions) const {
406 // First values for the chain with dispersion (always included)
407 auto vals = hardchain.eval(T, rhomolar, mole_fractions);
408 auto alphar = forceeval(vals.alphar_hc + vals.alphar_disp);
409
410 auto rho_A3 = forceeval(rhomolar*N_A*1e-30);
411 // If dipole is present, add its contribution
412 if (dipolar){
413 auto valsdip = dipolar.value().eval(T, rho_A3, vals.eta, mole_fractions);
414 alphar += valsdip.alpha;
415 }
416 // If quadrupole is present, add its contribution
417 if (quadrupolar){
418 auto valsquad = quadrupolar.value().eval(T, rho_A3, vals.eta, mole_fractions);
419 alphar += valsquad.alpha;
420 }
421 return forceeval(alphar);
422 }
423};
424
426inline auto PCSAFTfactory(const nlohmann::json& spec) {
427 std::optional<Eigen::ArrayXXd> kmat;
428 if (spec.contains("kmat") && spec.at("kmat").is_array() && spec.at("kmat").size() > 0){
429 kmat = build_square_matrix(spec["kmat"]);
430 }
431
432 if (spec.contains("names")){
433 std::vector<std::string> names = spec["names"];
434 if (kmat && static_cast<std::size_t>(kmat.value().rows()) != names.size()){
435 throw teqp::InvalidArgument("Provided length of names of " + std::to_string(names.size()) + " does not match the dimension of the kmat of " + std::to_string(kmat.value().rows()));
436 }
437 return PCSAFTMixture(names, kmat.value_or(Eigen::ArrayXXd{}));
438 }
439 else if (spec.contains("coeffs")){
440 std::vector<SAFTCoeffs> coeffs;
441 for (auto j : spec["coeffs"]) {
442 SAFTCoeffs c;
443 c.name = j.at("name");
444 c.m = j.at("m");
445 c.sigma_Angstrom = j.at("sigma_Angstrom");
446 c.epsilon_over_k = j.at("epsilon_over_k");
447 c.BibTeXKey = j.at("BibTeXKey");
448 if (j.contains("(mu^*)^2") && j.contains("nmu")){
449 c.mustar2 = j.at("(mu^*)^2");
450 c.nmu = j.at("nmu");
451 }
452 if (j.contains("(Q^*)^2") && j.contains("nQ")){
453 c.Qstar2 = j.at("(Q^*)^2");
454 c.nQ = j.at("nQ");
455 }
456 coeffs.push_back(c);
457 }
458 if (kmat && static_cast<std::size_t>(kmat.value().rows()) != coeffs.size()){
459 throw teqp::InvalidArgument("Provided length of coeffs of " + std::to_string(coeffs.size()) + " does not match the dimension of the kmat of " + std::to_string(kmat.value().rows()));
460 }
461 return PCSAFTMixture(coeffs, kmat.value_or(Eigen::ArrayXXd{}));
462 }
463 else{
464 throw std::invalid_argument("you must provide names or coeffs, but not both");
465 }
466}
467
468} /* namespace PCSAFT */
469}; // namespace teqp
PCSAFTHardChainContribution & operator=(const PCSAFTHardChainContribution &)=delete
const Eigen::ArrayX< double > sigma_Angstrom
Definition pcsaft.hpp:192
auto eval(const TTYPE &T, const RhoType &rhomolar, const VecType &mole_fractions) const
Definition pcsaft.hpp:203
PCSAFTHardChainContribution(const Eigen::ArrayX< double > &m, const Eigen::ArrayX< double > &mminus1, const Eigen::ArrayX< double > &sigma_Angstrom, const Eigen::ArrayX< double > &epsilon_over_k, const Eigen::ArrayXXd &kmat)
Definition pcsaft.hpp:197
const Eigen::ArrayX< double > epsilon_over_k
depth of pair potential divided by Boltzman constant
Definition pcsaft.hpp:193
const Eigen::ArrayXXd kmat
binary interaction parameter matrix
Definition pcsaft.hpp:194
const Eigen::ArrayX< double > mminus1
m-1
Definition pcsaft.hpp:191
const Eigen::ArrayX< double > m
number of segments
Definition pcsaft.hpp:190
Manager class for PCSAFT coefficients.
Definition pcsaft.hpp:34
void insert_normal_fluid(const std::string &name, double m, const double sigma_Angstrom, const double epsilon_over_k, const std::string &BibTeXKey)
Definition pcsaft.hpp:42
const auto & get_normal_fluid(const std::string &name)
Definition pcsaft.hpp:51
auto get_coeffs(const std::vector< std::string > &names)
Definition pcsaft.hpp:60
PCSAFTMixture(const std::vector< std::string > &names, const Eigen::ArrayXXd &kmat={})
Definition pcsaft.hpp:368
std::vector< std::string > names
Definition pcsaft.hpp:289
auto R(const VecType &molefrac) const
Definition pcsaft.hpp:400
auto build_hardchain(const std::vector< SAFTCoeffs > &coeffs)
Definition pcsaft.hpp:311
Eigen::ArrayXXd kmat
binary interaction parameter matrix
Definition pcsaft.hpp:290
auto build_dipolar(const std::vector< SAFTCoeffs > &coeffs) -> std::optional< PCSAFTDipolarContribution >
Definition pcsaft.hpp:339
auto alphar(const TTYPE &T, const RhoType &rhomolar, const VecType &mole_fractions) const
Definition pcsaft.hpp:405
Eigen::ArrayX< double > epsilon_over_k
depth of pair potential divided by Boltzman constant
Definition pcsaft.hpp:288
auto get_epsilon_over_k_K() const
Definition pcsaft.hpp:376
Eigen::ArrayX< double > sigma_Angstrom
Definition pcsaft.hpp:287
SAFTpolar::QuadrupolarContributionGross PCSAFTQuadrupolarContribution
Definition pcsaft.hpp:283
Eigen::ArrayX< double > m
number of segments
Definition pcsaft.hpp:285
SAFTpolar::DipolarContributionGrossVrabec PCSAFTDipolarContribution
Definition pcsaft.hpp:282
std::optional< PCSAFTDipolarContribution > dipolar
Definition pcsaft.hpp:293
auto build_quadrupolar(const std::vector< SAFTCoeffs > &coeffs) -> std::optional< PCSAFTQuadrupolarContribution >
Definition pcsaft.hpp:353
auto get_BibTeXKeys() const
Definition pcsaft.hpp:379
auto get_coeffs_from_names(const std::vector< std::string > &the_names)
Definition pcsaft.hpp:307
std::vector< std::string > bibtex
Definition pcsaft.hpp:289
PCSAFTMixture(const std::vector< SAFTCoeffs > &coeffs, const Eigen::ArrayXXd &kmat={})
Definition pcsaft.hpp:369
std::optional< PCSAFTQuadrupolarContribution > quadrupolar
Definition pcsaft.hpp:294
void check_kmat(Eigen::Index N)
Definition pcsaft.hpp:296
Eigen::ArrayX< double > mminus1
m-1
Definition pcsaft.hpp:286
PCSAFTHardChainContribution hardchain
Definition pcsaft.hpp:292
auto get_sigma_Angstrom() const
Definition pcsaft.hpp:375
auto extract_names(const std::vector< SAFTCoeffs > &coeffs)
Definition pcsaft.hpp:332
double max_rhoN(const double T, const VecType &mole_fractions) const
Definition pcsaft.hpp:390
PCSAFTMixture & operator=(const PCSAFTMixture &)=delete
Parameters for model evaluation.
Definition pcsaft.hpp:174
ProductType m2_epsilon2_sigma3_bar
Eq. A. 13.
Definition pcsaft.hpp:181
Eigen::ArrayX< NumType > d
Definition pcsaft.hpp:177
ProductType m2_epsilon_sigma3_bar
Eq. A. 12.
Definition pcsaft.hpp:180
auto powvec(const VecType1 &v1, NType n)
Definition pcsaft.hpp:155
auto get_I2(const Eta &eta, MbarType mbar)
Eqn. A.17, Eqn. A.30.
Definition pcsaft.hpp:140
auto get_I1(const Eta &eta, MbarType mbar)
Eqn. A.16, Eqn. A.29.
Definition pcsaft.hpp:128
auto gij_HS(const zVecType &zeta, const dVecType &d, std::size_t i, std::size_t j)
Term from Eqn. A.7.
Definition pcsaft.hpp:120
auto get_a(TYPE mbar)
Eqn. A.18.
Definition pcsaft.hpp:88
auto C1(const Eta &eta, Mbar mbar)
Definition pcsaft.hpp:73
auto C2(const Eta &eta, Mbar mbar)
Eqn. A.31.
Definition pcsaft.hpp:80
auto PCSAFTfactory(const nlohmann::json &spec)
A JSON-based factory function for the PC-SAFT model.
Definition pcsaft.hpp:426
auto sumproduct(const VecType1 &v1, const VecType2 &v2, const VecType3 &v3)
Definition pcsaft.hpp:167
auto get_b(TYPE mbar)
Eqn. A.19.
Definition pcsaft.hpp:96
auto get_alphar_hs(const VecType &zeta)
Residual contribution to alphar from hard-sphere (Eqn. A.6)
Definition pcsaft.hpp:105
auto build_square_matrix
auto getbaseval(const T &expr)
Definition types.hpp:84
T powi(const T &x, int n)
From Ulrich Deiters.
Definition types.hpp:133
auto pow(const double &x, const double &e)
Definition types.hpp:189
auto forceeval(T &&expr)
Definition types.hpp:46
Coefficients for one fluid.
Definition pcsaft.hpp:21
double epsilon_over_k
[K] depth of pair potential divided by Boltzman constant
Definition pcsaft.hpp:25
double Qstar2
nondimensional, the reduced quadrupole squared
Definition pcsaft.hpp:29
std::string BibTeXKey
The BibTeXKey for the reference for these coefficients.
Definition pcsaft.hpp:26
std::string name
Name of fluid.
Definition pcsaft.hpp:22
double mustar2
nondimensional, the reduced dipole moment squared
Definition pcsaft.hpp:27
double m
number of segments
Definition pcsaft.hpp:23
double nQ
number of quadrupolar segments
Definition pcsaft.hpp:30
double nmu
number of dipolar segments
Definition pcsaft.hpp:28
double sigma_Angstrom
[A] segment diameter
Definition pcsaft.hpp:24