Cantera  4.0.0a1
Loading...
Searching...
No Matches
RedlichKwongMFTP.cpp
Go to the documentation of this file.
1//! @file RedlichKwongMFTP.cpp
2
3// This file is part of Cantera. See License.txt in the top-level directory or
4// at https://cantera.org/license.txt for license and copyright information.
5
11
12#include <boost/algorithm/string.hpp>
13#include <algorithm>
14
15namespace Cantera
16{
17
18const double RedlichKwongMFTP::omega_a = 4.27480233540E-01;
19const double RedlichKwongMFTP::omega_b = 8.66403499650E-02;
20const double RedlichKwongMFTP::omega_vc = 3.33333333333333E-01;
21
22RedlichKwongMFTP::RedlichKwongMFTP(const string& infile, const string& id_)
23{
24 initThermoFile(infile, id_);
25}
26
27void RedlichKwongMFTP::setSpeciesCoeffs(const string& species,
28 double a0, double a1, double b)
29{
30 size_t k = speciesIndex(species, true);
31
32 if (a1 != 0.0) {
33 m_formTempParam = 1; // expression is temperature-dependent
34 }
35
36 size_t counter = k + m_kk * k;
37 a_coeff_vec(0, counter) = a0;
38 a_coeff_vec(1, counter) = a1;
39
40 // standard mixing rule for cross-species interaction term
41 for (size_t j = 0; j < m_kk; j++) {
42 if (k == j) {
43 continue;
44 }
45
46 // a_coeff_vec(0) is initialized to NaN to mark uninitialized species
47 if (isnan(a_coeff_vec(0, j + m_kk * j))) {
48 // The diagonal element of the jth species has not yet been defined.
49 continue;
50 } else if (isnan(a_coeff_vec(0, j + m_kk * k))) {
51 // Only use the mixing rules if the off-diagonal element has not already been defined by a
52 // user-specified crossFluidParameters entry:
53 double a0kj = sqrt(a_coeff_vec(0, j + m_kk * j) * a0);
54 double a1kj = sqrt(a_coeff_vec(1, j + m_kk * j) * a1);
55 a_coeff_vec(0, j + m_kk * k) = a0kj;
56 a_coeff_vec(1, j + m_kk * k) = a1kj;
57 a_coeff_vec(0, k + m_kk * j) = a0kj;
58 a_coeff_vec(1, k + m_kk * j) = a1kj;
59 }
60 }
61 a_coeff_vec.getRow(0, a_vec_Curr_);
62 b_vec_Curr_[k] = b;
63}
64
65void RedlichKwongMFTP::setBinaryCoeffs(const string& species_i, const string& species_j,
66 double a0, double a1)
67{
68 size_t ki = speciesIndex(species_i, true);
69 size_t kj = speciesIndex(species_j, true);
70
71 if (a1 != 0.0) {
72 m_formTempParam = 1; // expression is temperature-dependent
73 }
74
75 m_binaryParameters[species_i][species_j] = {a0, a1};
76 m_binaryParameters[species_j][species_i] = {a0, a1};
77 size_t counter1 = ki + m_kk * kj;
78 size_t counter2 = kj + m_kk * ki;
79 a_coeff_vec(0, counter1) = a_coeff_vec(0, counter2) = a0;
80 a_coeff_vec(1, counter1) = a_coeff_vec(1, counter2) = a1;
81 a_vec_Curr_[counter1] = a_vec_Curr_[counter2] = a0;
82}
83
84// ------------Molar Thermodynamic Properties -------------------------
85
87{
89 double cpref = GasConstant * mean_X(m_cp0_R);
90 if (m_b_current == 0.0) {
91 return cpref;
92 }
93 double TKelvin = temperature();
94 double sqt = sqrt(TKelvin);
95 double mv = molarVolume();
96 double vpb = mv + m_b_current;
98 double dadt = da_dt();
99 double fac = TKelvin * dadt - 3.0 * m_a_current / 2.0;
100 double dHdT_V = (cpref + mv * dpdT_ - GasConstant - 1.0 / (2.0 * m_b_current * TKelvin * sqt) * log(vpb/mv) * fac
101 +1.0/(m_b_current * sqt) * log(vpb/mv) * (-0.5 * dadt));
102 return dHdT_V - (mv + TKelvin * dpdT_ / dpdV_) * dpdT_;
103}
104
106{
108 double cvref = GasConstant * (mean_X(m_cp0_R) - 1.0);
109 if (m_b_current == 0.0) {
110 return cvref;
111 }
112 double TKelvin = temperature();
113 double sqt = sqrt(TKelvin);
114 double mv = molarVolume();
115 double vpb = mv + m_b_current;
116 double dadt = da_dt();
117 double fac = TKelvin * dadt - 3.0 * m_a_current / 2.0;
118 return (cvref - 1.0/(2.0 * m_b_current * TKelvin * sqt) * log(vpb/mv)*fac
119 +1.0/(m_b_current * sqt) * log(vpb/mv)*(-0.5*dadt));
120}
121
123{
125
126 // Get a copy of the private variables stored in the State object
127 double T = temperature();
128 double molarV = meanMolecularWeight() / density();
129 double pp = GasConstant * T/(molarV - m_b_current) - m_a_current/(sqrt(T) * molarV * (molarV + m_b_current));
130 return pp;
131}
132
134{
136 return 1.0 / m_workS[k];
137}
138
140{
141 checkArraySize("RedlichKwongMFTP::getActivityCoefficients", ac.size(), m_kk);
142 for (size_t k = 0; k < m_kk; k++) {
143 ac[k] = 1.0;
144 }
145 if (m_b_current == 0.0) {
146 return;
147 }
148 double mv = molarVolume();
149 double sqt = sqrt(temperature());
150 double vpb = mv + m_b_current;
151 double vmb = mv - m_b_current;
152
153 for (size_t k = 0; k < m_kk; k++) {
154 m_pp[k] = 0.0;
155 for (size_t i = 0; i < m_kk; i++) {
156 size_t counter = k + m_kk*i;
157 m_pp[k] += moleFractions_[i] * a_vec_Curr_[counter];
158 }
159 }
160 double pres = pressure();
161
162 for (size_t k = 0; k < m_kk; k++) {
163 ac[k] = (- RT() * log(pres * mv / RT())
164 + RT() * log(mv / vmb)
165 + RT() * b_vec_Curr_[k] / vmb
166 - 2.0 * m_pp[k] / (m_b_current * sqt) * log(vpb/mv)
167 + m_a_current * b_vec_Curr_[k] / (m_b_current * m_b_current * sqt) * log(vpb/mv)
168 - m_a_current / (m_b_current * sqt) * (b_vec_Curr_[k]/vpb)
169 );
170 }
171 for (size_t k = 0; k < m_kk; k++) {
172 ac[k] = exp(ac[k]/RT());
173 }
174}
175
176// ---- Partial Molar Properties of the Solution -----------------
177
178void RedlichKwongMFTP::getChemPotentials(span<double> mu) const
179{
181 for (size_t k = 0; k < m_kk; k++) {
182 double xx = std::max(SmallNumber, moleFraction(k));
183 mu[k] += RT() * log(xx);
184 }
185 if (m_b_current == 0.0) {
186 return;
187 }
188
189 double mv = molarVolume();
190 double sqt = sqrt(temperature());
191 double vpb = mv + m_b_current;
192 double vmb = mv - m_b_current;
193
194 for (size_t k = 0; k < m_kk; k++) {
195 m_pp[k] = 0.0;
196 for (size_t i = 0; i < m_kk; i++) {
197 size_t counter = k + m_kk*i;
198 m_pp[k] += moleFractions_[i] * a_vec_Curr_[counter];
199 }
200 }
201 double pres = pressure();
202
203 for (size_t k = 0; k < m_kk; k++) {
204 mu[k] += (- RT() * log(pres * mv / RT())
205 + RT() * log(mv / vmb)
206 + RT() * b_vec_Curr_[k] / vmb
207 - 2.0 * m_pp[k] / (m_b_current * sqt) * log(vpb/mv)
208 + m_a_current * b_vec_Curr_[k] / (m_b_current * m_b_current * sqt) * log(vpb/mv)
209 - m_a_current / (m_b_current * sqt) * (b_vec_Curr_[k]/vpb)
210 );
211 }
212}
213
215{
216 // First we get the reference state contributions
217 getEnthalpy_RT_ref(hbar);
218 scale(hbar.begin(), hbar.end(), hbar.begin(), RT());
219 if (m_b_current == 0.0) {
220 return;
221 }
222
223 // We calculate dpdni_
224 double TKelvin = temperature();
225 double mv = molarVolume();
226 double sqt = sqrt(TKelvin);
227 double vpb = mv + m_b_current;
228 double vmb = mv - m_b_current;
229 for (size_t k = 0; k < m_kk; k++) {
230 m_pp[k] = 0.0;
231 for (size_t i = 0; i < m_kk; i++) {
232 size_t counter = k + m_kk*i;
233 m_pp[k] += moleFractions_[i] * a_vec_Curr_[counter];
234 }
235 }
236 for (size_t k = 0; k < m_kk; k++) {
237 dpdni_[k] = RT()/vmb + RT() * b_vec_Curr_[k] / (vmb * vmb) - 2.0 * m_pp[k] / (sqt * mv * vpb)
238 + m_a_current * b_vec_Curr_[k]/(sqt * mv * vpb * vpb);
239 }
240 double dadt = da_dt();
241 double fac = TKelvin * dadt - 3.0 * m_a_current / 2.0;
242
243 for (size_t k = 0; k < m_kk; k++) {
244 m_workS[k] = 0.0;
245 for (size_t i = 0; i < m_kk; i++) {
246 size_t counter = k + m_kk*i;
247 m_workS[k] += 2.0 * moleFractions_[i] * TKelvin * a_coeff_vec(1,counter) - 3.0 * moleFractions_[i] * a_vec_Curr_[counter];
248 }
249 }
250
252 double fac2 = mv + TKelvin * dpdT_ / dpdV_;
253 for (size_t k = 0; k < m_kk; k++) {
254 double hE_v = (mv * dpdni_[k] - RT() - b_vec_Curr_[k]/ (m_b_current * m_b_current * sqt) * log(vpb/mv)*fac
255 + 1.0 / (m_b_current * sqt) * log(vpb/mv) * m_workS[k]
256 + b_vec_Curr_[k] / vpb / (m_b_current * sqt) * fac);
257 hbar[k] = hbar[k] + hE_v;
258 hbar[k] -= fac2 * dpdni_[k];
259 }
260}
261
263{
264 getEntropy_R_ref(sbar);
265 scale(sbar.begin(), sbar.end(), sbar.begin(), GasConstant);
266 double TKelvin = temperature();
267 double sqt = sqrt(TKelvin);
268 double mv = molarVolume();
269 double pres = pressure();
270 double logPres = log(pres / refPressure());
271
272 for (size_t k = 0; k < m_kk; k++) {
273 double xx = std::max(SmallNumber, moleFraction(k));
274 sbar[k] -= GasConstant * (log(xx) + logPres);
275 }
276 if (m_b_current == 0.0) {
277 return;
278 }
279 for (size_t k = 0; k < m_kk; k++) {
280 m_pp[k] = 0.0;
281 for (size_t i = 0; i < m_kk; i++) {
282 size_t counter = k + m_kk*i;
283 m_pp[k] += moleFractions_[i] * a_vec_Curr_[counter];
284 }
285 }
286 for (size_t k = 0; k < m_kk; k++) {
287 m_workS[k] = 0.0;
288 for (size_t i = 0; i < m_kk; i++) {
289 size_t counter = k + m_kk*i;
290 m_workS[k] += moleFractions_[i] * a_coeff_vec(1,counter);
291 }
292 }
293
294 double dadt = da_dt();
295 double fac = dadt - m_a_current / (2.0 * TKelvin);
296 double vmb = mv - m_b_current;
297 double vpb = mv + m_b_current;
298 for (size_t k = 0; k < m_kk; k++) {
299 sbar[k] += (GasConstant * log(pres * mv / RT())
301 - GasConstant * log(mv/vmb)
302 - GasConstant * b_vec_Curr_[k]/vmb
303 - m_pp[k]/(m_b_current * TKelvin * sqt) * log(vpb/mv)
304 + 2.0 * m_workS[k]/(m_b_current * sqt) * log(vpb/mv)
305 - b_vec_Curr_[k] / (m_b_current * m_b_current * sqt) * log(vpb/mv) * fac
306 + 1.0 / (m_b_current * sqt) * b_vec_Curr_[k] / vpb * fac);
307 }
308
310 getPartialMolarVolumes(m_partialMolarVolumes);
311 for (size_t k = 0; k < m_kk; k++) {
312 sbar[k] += m_partialMolarVolumes[k] * dpdT_;
313 }
314}
315
317{
318 // u_k = h_k - P * v_k
319 getPartialMolarVolumes(m_partialMolarVolumes);
321 double p = pressure();
322 for (size_t k = 0; k < nSpecies(); k++) {
323 ubar[k] -= p * m_partialMolarVolumes[k];
324 }
325}
326
328{
329 checkArraySize("RedlichKwongMFTP::getPartialMolarIntEnergies_TV", utilde.size(), m_kk);
331
332 double T = temperature();
333 double sqt = sqrt(T);
334 double mv = molarVolume();
335 double b = m_b_current;
336 double a = m_a_current;
337 double dadt = da_dt();
338
339 for (size_t k = 0; k < m_kk; k++) {
340 utilde[k] = RT() * (m_h0_RT[k] - 1.0);
341 }
342
343 if (fabs(b) < SmallNumber) {
344 return;
345 }
346
347 // A_k = sum_i x_i a_ki and A'_k = dA_k/dT = sum_i x_i a_ki,1
348 for (size_t k = 0; k < m_kk; k++) {
349 m_pp[k] = 0.0;
350 m_dAkdT[k] = 0.0;
351 for (size_t i = 0; i < m_kk; i++) {
352 size_t counter = k + m_kk * i;
353 m_pp[k] += moleFractions_[i] * a_vec_Curr_[counter];
354 m_dAkdT[k] += moleFractions_[i] * a_coeff_vec(1, counter);
355 }
356 }
357
358 double vpb = mv + b;
359 double logv = log(vpb / mv);
360 double F = T * dadt - 1.5 * a;
361 double C = -logv / b + 1.0 / vpb;
362 double pref = 1.0 / (b * sqt);
363
364 for (size_t k = 0; k < m_kk; k++) {
365 double Sk = 2.0 * T * m_dAkdT[k] - 3.0 * m_pp[k];
366 double ures = pref * (logv * Sk + b_vec_Curr_[k] * F * C);
367 utilde[k] += ures;
368 }
369}
370
371void RedlichKwongMFTP::getPartialMolarCp(span<double> cpbar) const
372{
373 checkArraySize("RedlichKwongMFTP::getPartialMolarCp", cpbar.size(), m_kk);
375
376 double T = temperature();
377 double sqt = sqrt(T);
378 double v = molarVolume();
379 double b = m_b_current;
380 double a = m_a_current;
381 double dadt = da_dt();
382 double d2adt2 = 0.0; // current RK model uses a(T)=a0+a1*T
383
384 for (size_t k = 0; k < m_kk; k++) {
385 cpbar[k] = GasConstant * m_cp0_R[k];
386 }
387
388 if (fabs(b) < SmallNumber) {
389 return;
390 }
391
392 // Mixture pressure derivatives
394 double pT = dpdT_;
395 double pV = dpdV_;
396 double dvdT_P = -pT / pV;
397
398 double vpb = v + b;
399 double vmb = v - b;
400
401 // P_TT, P_TV and P_VV for v''(T)|P from EOS:
402 // v'' = -(P_TT + 2 P_TV v' + P_VV v'^2) / P_V
403 double pTT = (-d2adt2 + dadt / T - 3.0 * a / (4.0 * T * T)) / (sqt * v * vpb);
404 double pTV = -GasConstant / (vmb * vmb)
405 + (dadt - a / (2.0 * T)) * (2.0 * v + b) / (sqt * v * v * vpb * vpb);
406 double pVV = 2.0 * GasConstant * T / (vmb * vmb * vmb)
407 - 2.0 * a * (3.0 * v * v + 3.0 * b * v + b * b)
408 / (sqt * v * v * v * vpb * vpb * vpb);
409 double d2vdT2_P = -(pTT + 2.0 * pTV * dvdT_P + pVV * dvdT_P * dvdT_P) / pV;
410
411 // Species-dependent A_k and dA_k/dT
412 for (size_t k = 0; k < m_kk; k++) {
413 m_pp[k] = 0.0;
414 m_dAkdT[k] = 0.0;
415 for (size_t i = 0; i < m_kk; i++) {
416 size_t counter = k + m_kk * i;
417 m_pp[k] += moleFractions_[i] * a_vec_Curr_[counter];
418 m_dAkdT[k] += moleFractions_[i] * a_coeff_vec(1, counter);
419 }
420 }
421
422 double F = T * dadt - 1.5 * a;
423 double dF_dT = -0.5 * dadt + T * d2adt2;
424 double logvpv = log(vpb / v);
425 double termLPrime = -b * dvdT_P / (v * vpb);
426
427 for (size_t k = 0; k < m_kk; k++) {
428 double bk = b_vec_Curr_[k];
429 double Ak = m_pp[k];
430 double dAk = m_dAkdT[k];
431 double Sk = 2.0 * T * dAk - 3.0 * Ak;
432 double dSk_dT = -dAk; // for linear a(T)
433
434 // dp/dn_k at constant T,V,n_j
435 double dpdni = RT() / vmb + RT() * bk / (vmb * vmb)
436 - 2.0 * Ak / (sqt * v * vpb)
437 + a * bk / (sqt * v * vpb * vpb);
438
439 // d/dT|P [dp/dn_k]
440 double d_dpdni_dT = GasConstant / vmb
441 - GasConstant * T * dvdT_P / (vmb * vmb)
442 + GasConstant * bk / (vmb * vmb)
443 - 2.0 * GasConstant * T * bk * dvdT_P / (vmb * vmb * vmb)
444 - 2.0 * dAk / (sqt * v * vpb)
445 + Ak / (T * sqt * v * vpb)
446 + 2.0 * Ak * (2.0 * v + b) * dvdT_P / (sqt * v * v * vpb * vpb)
447 + bk * dadt / (sqt * v * vpb * vpb)
448 - bk * a / (2.0 * T * sqt * v * vpb * vpb)
449 - bk * a * (3.0 * v + b) * dvdT_P / (sqt * v * v * vpb * vpb * vpb);
450
451 // d/dT|P of departure terms used in getPartialMolarEnthalpies
452 double dterm1 = -bk / (b * b * sqt)
453 * (termLPrime * F + logvpv * dF_dT - logvpv * F / (2.0 * T));
454 double dterm2 = 1.0 / (b * sqt)
455 * (termLPrime * Sk + logvpv * dSk_dT - logvpv * Sk / (2.0 * T));
456 double dterm3 = bk / (b * sqt)
457 * (dF_dT / vpb - F * dvdT_P / (vpb * vpb) - F / (2.0 * T * vpb));
458
459 // cp_k = d hbar_k / dT at constant P and composition
460 cpbar[k] += -GasConstant
461 + (dvdT_P + T * d2vdT2_P) * dpdni
462 + T * dvdT_P * d_dpdni_dT
463 + dterm1 + dterm2 + dterm3;
464 }
465}
466
467void RedlichKwongMFTP::getPartialMolarCv_TV(span<double> cvtilde) const
468{
469 checkArraySize("RedlichKwongMFTP::getPartialMolarCv_TV", cvtilde.size(), m_kk);
471
472 double T = temperature();
473 double sqt = sqrt(T);
474 double mv = molarVolume();
475 double b = m_b_current;
476 double a = m_a_current;
477 double dadt = da_dt();
478
479 for (size_t k = 0; k < m_kk; k++) {
480 cvtilde[k] = GasConstant * (m_cp0_R[k] - 1.0);
481 }
482
483 if (fabs(b) < SmallNumber) {
484 return;
485 }
486
487 // A_k = sum_i x_i a_ki and A'_k = dA_k/dT = sum_i x_i a_ki,1
488 for (size_t k = 0; k < m_kk; k++) {
489 m_pp[k] = 0.0;
490 m_dAkdT[k] = 0.0;
491 for (size_t i = 0; i < m_kk; i++) {
492 size_t counter = k + m_kk * i;
493 m_pp[k] += moleFractions_[i] * a_vec_Curr_[counter];
494 m_dAkdT[k] += moleFractions_[i] * a_coeff_vec(1, counter);
495 }
496 }
497
498 // Residual contribution for
499 // \tilde{u}_k = (\partial U / \partial n_k)_{T,V,n_j}
500 // with linear a(T): a_ij = a_ij,0 + a_ij,1 T
501 double vpb = mv + b;
502 double logv = log(vpb / mv);
503 double F = T * dadt - 1.5 * a;
504 double C = -logv / b + 1.0 / vpb;
505 double pref = 1.0 / (b * sqt);
506
507 for (size_t k = 0; k < m_kk; k++) {
508 double Sk = 2.0 * T * m_dAkdT[k] - 3.0 * m_pp[k];
509 double ures = pref * (logv * Sk + b_vec_Curr_[k] * F * C);
510 double dresdT = pref * (-logv * m_dAkdT[k] - 0.5 * b_vec_Curr_[k] * dadt * C)
511 - 0.5 * ures / T;
512 cvtilde[k] += dresdT;
513 }
514}
515
516void RedlichKwongMFTP::getPartialMolarVolumes(span<double> vbar) const
517{
518 checkArraySize("RedlichKwongMFTP::getPartialMolarVolumes", vbar.size(), m_kk);
519 for (size_t k = 0; k < m_kk; k++) {
520 m_pp[k] = 0.0;
521 for (size_t i = 0; i < m_kk; i++) {
522 size_t counter = k + m_kk*i;
523 m_pp[k] += moleFractions_[i] * a_vec_Curr_[counter];
524 }
525 }
526 for (size_t k = 0; k < m_kk; k++) {
527 m_workS[k] = 0.0;
528 for (size_t i = 0; i < m_kk; i++) {
529 size_t counter = k + m_kk*i;
530 m_workS[k] += moleFractions_[i] * a_coeff_vec(1,counter);
531 }
532 }
533
534 double sqt = sqrt(temperature());
535 double mv = molarVolume();
536 double vmb = mv - m_b_current;
537 double vpb = mv + m_b_current;
538 for (size_t k = 0; k < m_kk; k++) {
539 double num = (RT() + RT() * m_b_current/ vmb + RT() * b_vec_Curr_[k] / vmb
540 + RT() * m_b_current * b_vec_Curr_[k] /(vmb * vmb)
541 - 2.0 * m_pp[k] / (sqt * vpb)
542 + m_a_current * b_vec_Curr_[k] / (sqt * vpb * vpb)
543 );
544 double denom = (pressure() + RT() * m_b_current/(vmb * vmb) - m_a_current / (sqt * vpb * vpb)
545 );
546 vbar[k] = num / denom;
547 }
548}
549
550bool RedlichKwongMFTP::addSpecies(shared_ptr<Species> spec)
551{
552 bool added = MixtureFugacityTP::addSpecies(spec);
553 if (added) {
554 a_vec_Curr_.resize(m_kk * m_kk, 0.0);
555
556 // Initialize a_vec and b_vec to NaN, to screen for species with
557 // pureFluidParameters which are undefined in the input file:
558 b_vec_Curr_.push_back(NAN);
559 a_coeff_vec.resize(2, m_kk * m_kk, NAN);
560
561 m_pp.push_back(0.0);
562 m_dAkdT.push_back(0.0);
563 m_coeffSource.push_back(CoeffSource::EoS);
564 m_partialMolarVolumes.push_back(0.0);
565 dpdni_.push_back(0.0);
566 }
567 return added;
568}
569
571{
572 // Contents of 'critical-properties.yaml', loaded later if needed
573 AnyMap critPropsDb;
574 std::unordered_map<string, AnyMap*> dbSpecies;
575
576 for (auto& [name, species] : m_species) {
577 auto& data = species->input;
578 size_t k = speciesIndex(name, true);
579 if (!isnan(a_coeff_vec(0, k + m_kk * k))) {
580 continue;
581 }
582 bool foundCoeffs = false;
583
584 if (data.hasKey("equation-of-state") &&
585 data["equation-of-state"].hasMapWhere("model", "Redlich-Kwong"))
586 {
587 // Read a and b coefficients from species 'input' information (that is, as
588 // specified in a YAML input file) specific to the Redlich-Kwong model
589 auto eos = data["equation-of-state"].getMapWhere(
590 "model", "Redlich-Kwong");
591
592 if (eos.hasKey("a") && eos.hasKey("b")) {
593 double a0 = 0, a1 = 0;
594 if (eos["a"].isScalar()) {
595 a0 = eos.convert("a", "Pa*m^6/kmol^2*K^0.5");
596 } else {
597 auto avec = eos["a"].asVector<AnyValue>(2);
598 a0 = eos.units().convert(avec[0], "Pa*m^6/kmol^2*K^0.5");
599 a1 = eos.units().convert(avec[1], "Pa*m^6/kmol^2/K^0.5");
600 }
601 double b = eos.convert("b", "m^3/kmol");
602 foundCoeffs = true;
603 setSpeciesCoeffs(name, a0, a1, b);
604 m_coeffSource[k] = CoeffSource::EoS;
605 }
606
607 if (eos.hasKey("binary-a")) {
608 AnyMap& binary_a = eos["binary-a"].as<AnyMap>();
609 const UnitSystem& units = binary_a.units();
610 for (auto& [name2, coeff] : binary_a) {
611 double a0 = 0, a1 = 0;
612 if (coeff.isScalar()) {
613 a0 = units.convert(coeff, "Pa*m^6/kmol^2*K^0.5");
614 } else {
615 auto avec = coeff.asVector<AnyValue>(2);
616 a0 = units.convert(avec[0], "Pa*m^6/kmol^2*K^0.5");
617 a1 = units.convert(avec[1], "Pa*m^6/kmol^2/K^0.5");
618 }
619 setBinaryCoeffs(name, name2, a0, a1);
620 }
621 }
622
623 if (foundCoeffs) {
624 continue;
625 }
626 }
627
628 // Coefficients have not been populated from model-specific input
629 double Tc = NAN, Pc = NAN;
630 if (data.hasKey("critical-parameters")) {
631 // Use critical state information stored in the species entry to
632 // calculate a and b
633 auto& critProps = data["critical-parameters"].as<AnyMap>();
634 Tc = critProps.convert("critical-temperature", "K");
635 Pc = critProps.convert("critical-pressure", "Pa");
636 m_coeffSource[k] = CoeffSource::CritProps;
637 } else {
638 // Search 'crit-properties.yaml' to find Tc and Pc. Load data if needed
639 if (critPropsDb.empty()) {
640 critPropsDb = AnyMap::fromYamlFile("critical-properties.yaml");
641 dbSpecies = critPropsDb["species"].asMap("name");
642 }
643
644 // All names in critical-properties.yaml are upper case
645 auto ucName = boost::algorithm::to_upper_copy(name);
646 if (dbSpecies.count(ucName)) {
647 auto& spec = *dbSpecies.at(ucName);
648 auto& critProps = spec["critical-parameters"].as<AnyMap>();
649 Tc = critProps.convert("critical-temperature", "K");
650 Pc = critProps.convert("critical-pressure", "Pa");
651 m_coeffSource[k] = CoeffSource::Database;
652 }
653 }
654
655 // Check if critical properties were found in either location
656 if (!isnan(Tc)) {
657 // Assuming no temperature dependence (that is, a1 = 0)
658 double a = omega_a * pow(GasConstant, 2) * pow(Tc, 2.5) / Pc;
659 double b = omega_b * GasConstant * Tc / Pc;
660 setSpeciesCoeffs(name, a, 0.0, b);
661 } else {
662 throw InputFileError("RedlichKwongMFTP::initThermo", data,
663 "No critical property or Redlich-Kwong parameters found "
664 "for species {}.", name);
665 }
666 }
667}
668
670 AnyMap& speciesNode) const
671{
673 size_t k = speciesIndex(name, true);
674 if (m_coeffSource[k] == CoeffSource::EoS) {
675 auto& eosNode = speciesNode["equation-of-state"].getMapWhere(
676 "model", "Redlich-Kwong", true);
677
678 size_t counter = k + m_kk * k;
679 if (a_coeff_vec(1, counter) != 0.0) {
680 vector<AnyValue> coeffs(2);
681 coeffs[0].setQuantity(a_coeff_vec(0, counter), "Pa*m^6/kmol^2*K^0.5");
682 coeffs[1].setQuantity(a_coeff_vec(1, counter), "Pa*m^6/kmol^2/K^0.5");
683 eosNode["a"] = std::move(coeffs);
684 } else {
685 eosNode["a"].setQuantity(a_coeff_vec(0, counter),
686 "Pa*m^6/kmol^2*K^0.5");
687 }
688 eosNode["b"].setQuantity(b_vec_Curr_[k], "m^3/kmol");
689 } else if (m_coeffSource[k] == CoeffSource::CritProps) {
690 auto& critProps = speciesNode["critical-parameters"];
691 double a = a_coeff_vec(0, k + m_kk * k);
692 double b = b_vec_Curr_[k];
693 double Tc = pow(a * omega_b / (b * omega_a * GasConstant), 2.0/3.0);
694 double Pc = omega_b * GasConstant * Tc / b;
695 critProps["critical-temperature"].setQuantity(Tc, "K");
696 critProps["critical-pressure"].setQuantity(Pc, "Pa");
697 }
698
699 if (m_binaryParameters.count(name)) {
700 auto& eosNode = speciesNode["equation-of-state"].getMapWhere(
701 "model", "Redlich-Kwong", true);
702 AnyMap bin_a;
703 for (const auto& [name2, coeffs] : m_binaryParameters.at(name)) {
704 if (coeffs.second == 0) {
705 bin_a[name2].setQuantity(coeffs.first, "Pa*m^6/kmol^2*K^0.5");
706 } else {
707 vector<AnyValue> C(2);
708 C[0].setQuantity(coeffs.first, "Pa*m^6/kmol^2*K^0.5");
709 C[1].setQuantity(coeffs.second, "Pa*m^6/kmol^2/K^0.5");
710 bin_a[name2] = std::move(C);
711 }
712 }
713 eosNode["binary-a"] = std::move(bin_a);
714 }
715}
716
718{
719 if (m_b_current == 0.0) {
720 return 0.0;
721 }
722 // note this agrees with tpx
723 double rho = density();
724 double mmw = meanMolecularWeight();
725 double molarV = mmw / rho;
726 double hh = m_b_current / molarV;
727 double zz = z();
728 double dadt = da_dt();
729 double T = temperature();
730 double sqT = sqrt(T);
731 double fac = dadt - m_a_current / (2.0 * T);
732 double sresid_mol_R = log(zz*(1.0 - hh)) + log(1.0 + hh) * fac / (sqT * GasConstant * m_b_current);
733 return GasConstant * sresid_mol_R;
734}
735
737{
738 if (m_b_current == 0.0) {
739 return 0.0;
740 }
741 // note this agrees with tpx
742 double rho = density();
743 double mmw = meanMolecularWeight();
744 double molarV = mmw / rho;
745 double hh = m_b_current / molarV;
746 double zz = z();
747 double dadt = da_dt();
748 double T = temperature();
749 double sqT = sqrt(T);
750 double fac = T * dadt - 3.0 *m_a_current / (2.0);
751 return GasConstant * T * (zz - 1.0) + fac * log(1.0 + hh) / (sqT * m_b_current);
752}
753
754double RedlichKwongMFTP::liquidVolEst(double TKelvin, double& presGuess) const
755{
756 double v = m_b_current * 1.1;
757 double atmp;
758 double btmp;
759 calculateAB(TKelvin, atmp, btmp);
760 double pres = std::max(psatEst(TKelvin), presGuess);
761 double Vroot[3];
762 bool foundLiq = false;
763 int m = 0;
764 while (m < 100 && !foundLiq) {
765 int nsol = solveCubic(TKelvin, pres, atmp, btmp, Vroot);
766 if (nsol == 1 || nsol == 2) {
767 double pc = critPressure();
768 if (pres > pc) {
769 foundLiq = true;
770 }
771 pres *= 1.04;
772 } else {
773 foundLiq = true;
774 }
775 }
776
777 if (foundLiq) {
778 v = Vroot[0];
779 presGuess = pres;
780 } else {
781 v = -1.0;
782 }
783 return v;
784}
785
786double RedlichKwongMFTP::densityCalc(double TKelvin, double presPa, int phaseRequested,
787 double rhoguess)
788{
789 // It's necessary to set the temperature so that m_a_current is set correctly.
790 setTemperature(TKelvin);
791 double tcrit = critTemperature();
792 double mmw = meanMolecularWeight();
793 if (rhoguess == -1.0) {
794 if (phaseRequested >= FLUID_LIQUID_0) {
795 double lqvol = liquidVolEst(TKelvin, presPa);
796 rhoguess = mmw / lqvol;
797 }
798 } else {
799 // Assume the Gas phase initial guess, if nothing is specified to the routine
800 rhoguess = presPa * mmw / (GasConstant * TKelvin);
801 }
802
803
804 double volguess = mmw / rhoguess;
805 NSolns_ = solveCubic(TKelvin, presPa, m_a_current, m_b_current, Vroot_);
806
807 // Ensure root ordering used for branch selection only contains physical roots.
808 double vmin = std::max(0.0, m_b_current * (1.0 + 1e-12));
809 vector<double> physicalRoots;
810 for (double root : Vroot_) {
811 if (std::isfinite(root) && root > vmin) {
812 physicalRoots.push_back(root);
813 }
814 }
815 std::sort(physicalRoots.begin(), physicalRoots.end());
816 if (physicalRoots.empty()) {
817 return -1.0;
818 } else if (physicalRoots.size() == 1) {
819 // Preserve branch-request semantics for single-root states:
820 // return -2 when the requested phase is inconsistent with the
821 // single branch identified by solveCubic() sign convention.
822 if ((phaseRequested == FLUID_GAS && NSolns_ < 0)
823 || (phaseRequested >= FLUID_LIQUID_0 && NSolns_ > 0))
824 {
825 return -2.0;
826 }
827 // Otherwise, accept the physically admissible root directly.
828 return mmw / physicalRoots[0];
829 } else if (physicalRoots.size() == 2) {
830 Vroot_[0] = physicalRoots[0];
831 Vroot_[1] = 0.5 * (physicalRoots[0] + physicalRoots[1]);
832 Vroot_[2] = physicalRoots[1];
833 } else {
834 Vroot_[0] = physicalRoots[0];
835 Vroot_[1] = physicalRoots[1];
836 Vroot_[2] = physicalRoots[2];
837 }
838
839 double molarVolLast = Vroot_[0];
840 if (NSolns_ >= 2) {
841 if (phaseRequested >= FLUID_LIQUID_0) {
842 molarVolLast = Vroot_[0];
843 } else if (phaseRequested == FLUID_GAS || phaseRequested == FLUID_SUPERCRIT) {
844 molarVolLast = Vroot_[2];
845 } else {
846 if (volguess > Vroot_[1]) {
847 molarVolLast = Vroot_[2];
848 } else {
849 molarVolLast = Vroot_[0];
850 }
851 }
852 } else if (NSolns_ == 1) {
853 if (phaseRequested == FLUID_GAS || phaseRequested == FLUID_SUPERCRIT || phaseRequested == FLUID_UNDEFINED) {
854 molarVolLast = Vroot_[0];
855 } else {
856 return -2.0;
857 }
858 } else if (NSolns_ == -1) {
859 if (phaseRequested >= FLUID_LIQUID_0 || phaseRequested == FLUID_UNDEFINED || phaseRequested == FLUID_SUPERCRIT) {
860 molarVolLast = Vroot_[0];
861 } else if (TKelvin > tcrit) {
862 molarVolLast = Vroot_[0];
863 } else {
864 return -2.0;
865 }
866 } else {
867 molarVolLast = Vroot_[0];
868 return -1.0;
869 }
870 return mmw / molarVolLast;
871}
872
873double RedlichKwongMFTP::dpdVCalc(double TKelvin, double molarVol, double& presCalc) const
874{
875 double sqt = sqrt(TKelvin);
876 presCalc = GasConstant * TKelvin / (molarVol - m_b_current)
877 - m_a_current / (sqt * molarVol * (molarVol + m_b_current));
878
879 double vpb = molarVol + m_b_current;
880 double vmb = molarVol - m_b_current;
881 double dpdv = (- GasConstant * TKelvin / (vmb * vmb)
882 + m_a_current * (2 * molarVol + m_b_current) / (sqt * molarVol * molarVol * vpb * vpb));
883 return dpdv;
884}
885
887{
889 return -1 / (molarVolume() * dpdV_);
890}
891
893{
895 return -dpdT_ / (molarVolume() * dpdV_);
896}
897
899{
901 return temperature() * dpdT_ - pressure();
902}
903
905{
907 return molarVolume() * sqrt(-cp_mole() / cv_mole() * dpdV_ / meanMolecularWeight());
908}
909
911{
912 double TKelvin = temperature();
913 double mv = molarVolume();
914 double pres;
915
916 dpdV_ = dpdVCalc(TKelvin, mv, pres);
917 double sqt = sqrt(TKelvin);
918 double vpb = mv + m_b_current;
919 double vmb = mv - m_b_current;
920 double dadt = da_dt();
921 double fac = dadt - m_a_current/(2.0 * TKelvin);
922 dpdT_ = (GasConstant / vmb - fac / (sqt * mv * vpb));
923}
924
926{
927 double temp = temperature();
928 if (m_formTempParam == 1) {
929 for (size_t i = 0; i < m_kk; i++) {
930 for (size_t j = 0; j < m_kk; j++) {
931 size_t counter = i * m_kk + j;
932 a_vec_Curr_[counter] = a_coeff_vec(0,counter) + a_coeff_vec(1,counter) * temp;
933 }
934 }
935 }
936
937 m_b_current = 0.0;
938 m_a_current = 0.0;
939 for (size_t i = 0; i < m_kk; i++) {
940 m_b_current += moleFractions_[i] * b_vec_Curr_[i];
941 for (size_t j = 0; j < m_kk; j++) {
942 m_a_current += a_vec_Curr_[i * m_kk + j] * moleFractions_[i] * moleFractions_[j];
943 }
944 }
945 if (isnan(m_b_current)) {
946 // One or more species do not have specified coefficients.
947 fmt::memory_buffer b;
948 for (size_t k = 0; k < m_kk; k++) {
949 if (isnan(b_vec_Curr_[k])) {
950 if (b.size() > 0) {
951 fmt_append(b, ", {}", speciesName(k));
952 } else {
953 fmt_append(b, "{}", speciesName(k));
954 }
955 }
956 }
957 throw CanteraError("RedlichKwongMFTP::updateMixingExpressions",
958 "Missing Redlich-Kwong coefficients for species: {}", to_string(b));
959 }
960}
961
962void RedlichKwongMFTP::calculateAB(double temp, double& aCalc, double& bCalc) const
963{
964 bCalc = 0.0;
965 aCalc = 0.0;
966 if (m_formTempParam == 1) {
967 for (size_t i = 0; i < m_kk; i++) {
968 bCalc += moleFractions_[i] * b_vec_Curr_[i];
969 for (size_t j = 0; j < m_kk; j++) {
970 size_t counter = i * m_kk + j;
971 double a_vec_Curr = a_coeff_vec(0,counter) + a_coeff_vec(1,counter) * temp;
972 aCalc += a_vec_Curr * moleFractions_[i] * moleFractions_[j];
973 }
974 }
975 } else {
976 for (size_t i = 0; i < m_kk; i++) {
977 bCalc += moleFractions_[i] * b_vec_Curr_[i];
978 for (size_t j = 0; j < m_kk; j++) {
979 size_t counter = i * m_kk + j;
980 double a_vec_Curr = a_coeff_vec(0,counter);
981 aCalc += a_vec_Curr * moleFractions_[i] * moleFractions_[j];
982 }
983 }
984 }
985}
986
987double RedlichKwongMFTP::da_dt() const
988{
989 double dadT = 0.0;
990 if (m_formTempParam == 1) {
991 for (size_t i = 0; i < m_kk; i++) {
992 for (size_t j = 0; j < m_kk; j++) {
993 size_t counter = i * m_kk + j;
994 dadT+= a_coeff_vec(1,counter) * moleFractions_[i] * moleFractions_[j];
995 }
996 }
997 }
998 return dadT;
999}
1000
1001void RedlichKwongMFTP::calcCriticalConditions(double& pc, double& tc, double& vc) const
1002{
1003 double a0 = 0.0;
1004 double aT = 0.0;
1005 for (size_t i = 0; i < m_kk; i++) {
1006 for (size_t j = 0; j <m_kk; j++) {
1007 size_t counter = i + m_kk * j;
1008 a0 += moleFractions_[i] * moleFractions_[j] * a_coeff_vec(0, counter);
1009 aT += moleFractions_[i] * moleFractions_[j] * a_coeff_vec(1, counter);
1010 }
1011 }
1012 double a = m_a_current;
1013 double b = m_b_current;
1014 if (m_formTempParam != 0) {
1015 a = a0;
1016 }
1017 if (b <= 0.0) {
1018 tc = 1000000.;
1019 pc = 1.0E13;
1020 vc = omega_vc * GasConstant * tc / pc;
1021 return;
1022 }
1023 if (a <= 0.0) {
1024 tc = 0.0;
1025 pc = 0.0;
1026 vc = 2.0 * b;
1027 return;
1028 }
1029 double tmp = a * omega_b / (b * omega_a * GasConstant);
1030 double pp = 2./3.;
1031 double sqrttc, f, dfdt, deltatc;
1032
1033 if (m_formTempParam == 0) {
1034 tc = pow(tmp, pp);
1035 } else {
1036 tc = pow(tmp, pp);
1037 for (int j = 0; j < 10; j++) {
1038 sqrttc = sqrt(tc);
1039 f = omega_a * b * GasConstant * tc * sqrttc / omega_b - aT * tc - a0;
1040 dfdt = 1.5 * omega_a * b * GasConstant * sqrttc / omega_b - aT;
1041 deltatc = - f / dfdt;
1042 tc += deltatc;
1043 }
1044 if (deltatc > 0.1) {
1045 throw CanteraError("RedlichKwongMFTP::calcCriticalConditions",
1046 "didn't converge");
1047 }
1048 }
1049
1050 pc = omega_b * GasConstant * tc / b;
1051 vc = omega_vc * GasConstant * tc / pc;
1052}
1053
1054int RedlichKwongMFTP::solveCubic(double T, double pres, double a, double b,
1055 span<double> Vroot) const
1056{
1057
1058 // Derive the coefficients of the cubic polynomial to solve.
1059 double an = 1.0;
1060 double bn = - GasConstant * T / pres;
1061 double sqt = sqrt(T);
1062 double cn = - (GasConstant * T * b / pres - a/(pres * sqt) + b * b);
1063 double dn = - (a * b / (pres * sqt));
1064
1065 double tmp = a * omega_b / (b * omega_a * GasConstant);
1066 double pp = 2./3.;
1067 double tc = pow(tmp, pp);
1068 double pc = omega_b * GasConstant * tc / b;
1069 double vc = omega_vc * GasConstant * tc / pc;
1070
1071 int nSolnValues = MixtureFugacityTP::solveCubic(T, pres, a, b, a, Vroot, an, bn, cn, dn, tc, vc);
1072
1073 return nSolnValues;
1074}
1075
1076}
Declaration for class Cantera::Species.
Headers for the factory class that can create known ThermoPhase objects (see Thermodynamic Properties...
A map of string keys to values whose type can vary at runtime.
Definition AnyMap.h:431
const UnitSystem & units() const
Return the default units that should be used to convert stored values.
Definition AnyMap.h:640
bool empty() const
Return boolean indicating whether AnyMap is empty.
Definition AnyMap.cpp:1468
double convert(const string &key, const string &units) const
Convert the item stored by the given key to the units specified in units.
Definition AnyMap.cpp:1595
static AnyMap fromYamlFile(const string &name, const string &parent_name="")
Create an AnyMap from a YAML file.
Definition AnyMap.cpp:1841
A wrapper for a variable whose type is determined at runtime.
Definition AnyMap.h:88
void getRow(size_t n, span< double > rw) const
Get the nth row and return it in a vector.
Definition Array.cpp:77
virtual void resize(size_t n, size_t m, double v=0.0)
Resize the array, and fill the new entries with 'v'.
Definition Array.cpp:52
Base class for exceptions thrown by Cantera classes.
Error thrown for problems processing information contained in an AnyMap or AnyValue.
Definition AnyMap.h:749
void getEntropy_R_ref(span< double > er) const override
Returns the vector of nondimensional entropies of the reference state at the current temperature of t...
double critPressure() const override
Critical pressure (Pa).
vector< double > m_h0_RT
Temporary storage for dimensionless reference state enthalpies.
void getStandardChemPotentials(span< double > mu) const override
Get the array of chemical potentials at unit activity.
double critTemperature() const override
Critical temperature (K).
virtual void _updateReferenceStateThermo() const
Updates the reference state thermodynamic functions at the current T of the solution.
void getEnthalpy_RT_ref(span< double > hrt) const override
Returns the vector of nondimensional enthalpies of the reference state at the current temperature of ...
vector< double > moleFractions_
Storage for the current values of the mole fractions of the species.
int solveCubic(double T, double pres, double a, double b, double aAlpha, span< double > Vroot, double an, double bn, double cn, double dn, double tc, double vc) const
Solve the cubic equation of state.
void setTemperature(const double temp) override
Set the temperature of the phase.
virtual double psatEst(double TKelvin) const
Estimate for the saturation pressure.
void getStandardVolumes(span< double > vol) const override
Get the molar volumes of each species in their standard states at the current T and P of the solution...
vector< double > m_cp0_R
Temporary storage for dimensionless reference state heat capacities.
bool addSpecies(shared_ptr< Species > spec) override
Add a Species to this Phase.
double z() const
Calculate the value of z.
vector< double > m_workS
Vector of size m_kk, used as a temporary holding area.
Definition Phase.h:899
size_t nSpecies() const
Returns the number of species in the phase.
Definition Phase.h:246
size_t m_kk
Number of species in the phase.
Definition Phase.h:875
size_t speciesIndex(const string &name, bool raise=true) const
Returns the index of a species named 'name' within the Phase object.
Definition Phase.cpp:127
double temperature() const
Temperature (K).
Definition Phase.h:585
double meanMolecularWeight() const
The mean molecular weight. Units: (kg/kmol)
Definition Phase.h:676
string speciesName(size_t k) const
Name of the species with index k.
Definition Phase.cpp:143
map< string, shared_ptr< Species > > m_species
Map of Species objects.
Definition Phase.h:888
double mean_X(span< const double > Q) const
Evaluate the mole-fraction-weighted mean of an array Q.
Definition Phase.cpp:627
double moleFraction(size_t k) const
Return the mole fraction of a single species.
Definition Phase.cpp:447
virtual double density() const
Density (kg/m^3).
Definition Phase.h:610
shared_ptr< Species > species(const string &name) const
Return the Species object for the named species.
Definition Phase.cpp:881
virtual double molarVolume() const
Molar volume (m^3/kmol).
Definition Phase.cpp:592
string name() const
Return the name of the phase.
Definition Phase.cpp:20
double dpdT_
The derivative of the pressure wrt the temperature.
void calculateAB(double temp, double &aCalc, double &bCalc) const
Calculate the a and the b parameters given the temperature.
double thermalExpansionCoeff() const override
Return the volumetric thermal expansion coefficient. Units: 1/K.
void setBinaryCoeffs(const string &species_i, const string &species_j, double a0, double a1)
Set values for the interaction parameter between two species.
vector< double > m_dAkdT
Temporary storage for dA_k/dT; length m_kk.
double sresid() const override
Calculate the deviation terms for the total entropy of the mixture from the ideal gas mixture.
void getPartialMolarEnthalpies(span< double > hbar) const override
Return partial molar enthalpies (J/kmol).
double soundSpeed() const override
Return the speed of sound. Units: m/s.
double pressure() const override
Return the thermodynamic pressure (Pa).
int m_formTempParam
Form of the temperature parameterization.
void getSpeciesParameters(const string &name, AnyMap &speciesNode) const override
Get phase-specific parameters of a Species object such that an identical one could be reconstructed a...
RedlichKwongMFTP(const string &infile="", const string &id="")
Construct a RedlichKwongMFTP object from an input file.
void getPartialMolarCp(span< double > cpbar) const override
Get the partial molar heat capacities at constant pressure [J/kmol/K].
static const double omega_b
Omega constant for b.
vector< double > m_pp
Temporary storage - length = m_kk.
void getActivityCoefficients(span< double > ac) const override
Get the array of non-dimensional activity coefficients at the current solution temperature,...
void initThermo() override
Initialize the ThermoPhase object after all species have been set up.
void getPartialMolarIntEnergies_TV(span< double > utilde) const override
Get the species molar internal energies associated with the derivatives of total internal energy at c...
double internalPressure() const override
Return the internal pressure [Pa].
double cv_mole() const override
Molar heat capacity at constant volume and composition [J/kmol/K].
void pressureDerivatives() const
Calculate dpdV and dpdT at the current conditions.
double isothermalCompressibility() const override
Returns the isothermal compressibility. Units: 1/Pa.
double hresid() const override
Calculate the deviation terms for the total enthalpy of the mixture from the ideal gas mixture.
static const double omega_a
Omega constant for a -> value of a in terms of critical properties.
double liquidVolEst(double TKelvin, double &pres) const override
Estimate for the molar volume of the liquid.
void getPartialMolarCv_TV(span< double > cvtilde) const override
Get the species molar heat capacities associated with the constant volume partial molar internal ener...
double dpdVCalc(double TKelvin, double molarVol, double &presCalc) const override
Calculate the pressure and the pressure derivative given the temperature and the molar volume.
static const double omega_vc
Omega constant for the critical molar volume.
vector< CoeffSource > m_coeffSource
For each species, specifies the source of the a and b coefficients.
void updateMixingExpressions() override
Update the a and b parameters.
double cp_mole() const override
Molar heat capacity at constant pressure and composition [J/kmol/K].
void getPartialMolarVolumes(span< double > vbar) const override
Return an array of partial molar volumes for the species in the mixture.
double m_a_current
Value of a in the equation of state.
void getPartialMolarEntropies(span< double > sbar) const override
Returns an array of partial molar entropies of the species in the solution.
double standardConcentration(size_t k=0) const override
Returns the standard concentration , which is used to normalize the generalized concentration.
double dpdV_
The derivative of the pressure wrt the volume.
bool addSpecies(shared_ptr< Species > spec) override
Add a Species to this Phase.
double m_b_current
Value of b in the equation of state.
void getChemPotentials(span< double > mu) const override
Get the species chemical potentials. Units: J/kmol.
int solveCubic(double T, double pres, double a, double b, span< double > Vroot) const
Prepare variables and call the function to solve the cubic equation of state.
double densityCalc(double T, double pressure, int phase, double rhoguess) override
Calculates the density given the temperature and the pressure and a guess at the density.
vector< double > dpdni_
Vector of derivatives of pressure wrt mole number.
map< string, map< string, pair< double, double > > > m_binaryParameters
Explicitly-specified binary interaction parameters.
void setSpeciesCoeffs(const string &species, double a0, double a1, double b)
Set the pure fluid interaction parameters for a species.
void getPartialMolarIntEnergies(span< double > ubar) const override
Return an array of partial molar internal energies for the species in the mixture.
double RT() const
Return the Gas Constant multiplied by the current temperature.
void initThermoFile(const string &inputFile, const string &id)
Initialize a ThermoPhase object using an input file.
shared_ptr< Solution > root() const
Get the Solution object containing this ThermoPhase object and linked Kinetics and Transport objects.
virtual void getSpeciesParameters(const string &name, AnyMap &speciesNode) const
Get phase-specific parameters of a Species object such that an identical one could be reconstructed a...
virtual double refPressure() const
Returns the reference pressure in Pa.
Unit conversion utility.
Definition Units.h:169
double convert(double value, const string &src, const string &dest) const
Convert value from the units of src to the units of dest.
Definition Units.cpp:539
void scale(InputIter begin, InputIter end, OutputIter out, S scale_factor)
Multiply elements of an array by a scale factor.
Definition utilities.h:118
const double GasConstant
Universal Gas Constant [J/kmol/K].
Definition ct_defs.h:123
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595
const double SmallNumber
smallest number to compare to zero.
Definition ct_defs.h:161
void checkArraySize(const char *procedure, size_t available, size_t required)
Wrapper for throwing ArraySizeError.
Contains declarations for string manipulation functions within Cantera.
Various templated functions that carry out common vector and polynomial operations (see Templated Arr...