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
14namespace Cantera
15{
16
17const double RedlichKwongMFTP::omega_a = 4.27480233540E-01;
18const double RedlichKwongMFTP::omega_b = 8.66403499650E-02;
19const double RedlichKwongMFTP::omega_vc = 3.33333333333333E-01;
20
21RedlichKwongMFTP::RedlichKwongMFTP(const string& infile, const string& id_)
22{
23 initThermoFile(infile, id_);
24}
25
26void RedlichKwongMFTP::setSpeciesCoeffs(const string& species,
27 double a0, double a1, double b)
28{
29 size_t k = speciesIndex(species, true);
30
31 if (a1 != 0.0) {
32 m_formTempParam = 1; // expression is temperature-dependent
33 }
34
35 size_t counter = k + m_kk * k;
36 a_coeff_vec(0, counter) = a0;
37 a_coeff_vec(1, counter) = a1;
38
39 // standard mixing rule for cross-species interaction term
40 for (size_t j = 0; j < m_kk; j++) {
41 if (k == j) {
42 continue;
43 }
44
45 // a_coeff_vec(0) is initialized to NaN to mark uninitialized species
46 if (isnan(a_coeff_vec(0, j + m_kk * j))) {
47 // The diagonal element of the jth species has not yet been defined.
48 continue;
49 } else if (isnan(a_coeff_vec(0, j + m_kk * k))) {
50 // Only use the mixing rules if the off-diagonal element has not already been defined by a
51 // user-specified crossFluidParameters entry:
52 double a0kj = sqrt(a_coeff_vec(0, j + m_kk * j) * a0);
53 double a1kj = sqrt(a_coeff_vec(1, j + m_kk * j) * a1);
54 a_coeff_vec(0, j + m_kk * k) = a0kj;
55 a_coeff_vec(1, j + m_kk * k) = a1kj;
56 a_coeff_vec(0, k + m_kk * j) = a0kj;
57 a_coeff_vec(1, k + m_kk * j) = a1kj;
58 }
59 }
60 a_coeff_vec.getRow(0, a_vec_Curr_);
61 b_vec_Curr_[k] = b;
62}
63
64void RedlichKwongMFTP::setBinaryCoeffs(const string& species_i, const string& species_j,
65 double a0, double a1)
66{
67 size_t ki = speciesIndex(species_i, true);
68 size_t kj = speciesIndex(species_j, true);
69
70 if (a1 != 0.0) {
71 m_formTempParam = 1; // expression is temperature-dependent
72 }
73
74 m_binaryParameters[species_i][species_j] = {a0, a1};
75 m_binaryParameters[species_j][species_i] = {a0, a1};
76 size_t counter1 = ki + m_kk * kj;
77 size_t counter2 = kj + m_kk * ki;
78 a_coeff_vec(0, counter1) = a_coeff_vec(0, counter2) = a0;
79 a_coeff_vec(1, counter1) = a_coeff_vec(1, counter2) = a1;
80 a_vec_Curr_[counter1] = a_vec_Curr_[counter2] = a0;
81}
82
83// ------------Molar Thermodynamic Properties -------------------------
84
86{
88 double TKelvin = temperature();
89 double sqt = sqrt(TKelvin);
90 double mv = molarVolume();
91 double vpb = mv + m_b_current;
93 double cpref = GasConstant * mean_X(m_cp0_R);
94 double dadt = da_dt();
95 double fac = TKelvin * dadt - 3.0 * m_a_current / 2.0;
96 double dHdT_V = (cpref + mv * dpdT_ - GasConstant - 1.0 / (2.0 * m_b_current * TKelvin * sqt) * log(vpb/mv) * fac
97 +1.0/(m_b_current * sqt) * log(vpb/mv) * (-0.5 * dadt));
98 return dHdT_V - (mv + TKelvin * dpdT_ / dpdV_) * dpdT_;
99}
100
102{
104 double TKelvin = temperature();
105 double sqt = sqrt(TKelvin);
106 double mv = molarVolume();
107 double vpb = mv + m_b_current;
108 double cvref = GasConstant * (mean_X(m_cp0_R) - 1.0);
109 double dadt = da_dt();
110 double fac = TKelvin * dadt - 3.0 * m_a_current / 2.0;
111 return (cvref - 1.0/(2.0 * m_b_current * TKelvin * sqt) * log(vpb/mv)*fac
112 +1.0/(m_b_current * sqt) * log(vpb/mv)*(-0.5*dadt));
113}
114
116{
118
119 // Get a copy of the private variables stored in the State object
120 double T = temperature();
121 double molarV = meanMolecularWeight() / density();
122 double pp = GasConstant * T/(molarV - m_b_current) - m_a_current/(sqrt(T) * molarV * (molarV + m_b_current));
123 return pp;
124}
125
127{
129 return 1.0 / m_workS[k];
130}
131
133{
134 checkArraySize("RedlichKwongMFTP::getActivityCoefficients", ac.size(), m_kk);
135 double mv = molarVolume();
136 double sqt = sqrt(temperature());
137 double vpb = mv + m_b_current;
138 double vmb = mv - m_b_current;
139
140 for (size_t k = 0; k < m_kk; k++) {
141 m_pp[k] = 0.0;
142 for (size_t i = 0; i < m_kk; i++) {
143 size_t counter = k + m_kk*i;
144 m_pp[k] += moleFractions_[i] * a_vec_Curr_[counter];
145 }
146 }
147 double pres = pressure();
148
149 for (size_t k = 0; k < m_kk; k++) {
150 ac[k] = (- RT() * log(pres * mv / RT())
151 + RT() * log(mv / vmb)
152 + RT() * b_vec_Curr_[k] / vmb
153 - 2.0 * m_pp[k] / (m_b_current * sqt) * log(vpb/mv)
154 + m_a_current * b_vec_Curr_[k] / (m_b_current * m_b_current * sqt) * log(vpb/mv)
155 - m_a_current / (m_b_current * sqt) * (b_vec_Curr_[k]/vpb)
156 );
157 }
158 for (size_t k = 0; k < m_kk; k++) {
159 ac[k] = exp(ac[k]/RT());
160 }
161}
162
163// ---- Partial Molar Properties of the Solution -----------------
164
165void RedlichKwongMFTP::getChemPotentials(span<double> mu) const
166{
167 getGibbs_ref(mu);
168 for (size_t k = 0; k < m_kk; k++) {
169 double xx = std::max(SmallNumber, moleFraction(k));
170 mu[k] += RT()*(log(xx));
171 }
172
173 double mv = molarVolume();
174 double sqt = sqrt(temperature());
175 double vpb = mv + m_b_current;
176 double vmb = mv - m_b_current;
177
178 for (size_t k = 0; k < m_kk; k++) {
179 m_pp[k] = 0.0;
180 for (size_t i = 0; i < m_kk; i++) {
181 size_t counter = k + m_kk*i;
182 m_pp[k] += moleFractions_[i] * a_vec_Curr_[counter];
183 }
184 }
185 double pres = pressure();
186 double refP = refPressure();
187
188 for (size_t k = 0; k < m_kk; k++) {
189 mu[k] += (RT() * log(pres/refP) - RT() * log(pres * mv / RT())
190 + RT() * log(mv / vmb)
191 + RT() * b_vec_Curr_[k] / vmb
192 - 2.0 * m_pp[k] / (m_b_current * sqt) * log(vpb/mv)
193 + m_a_current * b_vec_Curr_[k] / (m_b_current * m_b_current * sqt) * log(vpb/mv)
194 - m_a_current / (m_b_current * sqt) * (b_vec_Curr_[k]/vpb)
195 );
196 }
197}
198
200{
201 // First we get the reference state contributions
202 getEnthalpy_RT_ref(hbar);
203 scale(hbar.begin(), hbar.end(), hbar.begin(), RT());
204
205 // We calculate dpdni_
206 double TKelvin = temperature();
207 double mv = molarVolume();
208 double sqt = sqrt(TKelvin);
209 double vpb = mv + m_b_current;
210 double vmb = mv - m_b_current;
211 for (size_t k = 0; k < m_kk; k++) {
212 m_pp[k] = 0.0;
213 for (size_t i = 0; i < m_kk; i++) {
214 size_t counter = k + m_kk*i;
215 m_pp[k] += moleFractions_[i] * a_vec_Curr_[counter];
216 }
217 }
218 for (size_t k = 0; k < m_kk; k++) {
219 dpdni_[k] = RT()/vmb + RT() * b_vec_Curr_[k] / (vmb * vmb) - 2.0 * m_pp[k] / (sqt * mv * vpb)
220 + m_a_current * b_vec_Curr_[k]/(sqt * mv * vpb * vpb);
221 }
222 double dadt = da_dt();
223 double fac = TKelvin * dadt - 3.0 * m_a_current / 2.0;
224
225 for (size_t k = 0; k < m_kk; k++) {
226 m_workS[k] = 0.0;
227 for (size_t i = 0; i < m_kk; i++) {
228 size_t counter = k + m_kk*i;
229 m_workS[k] += 2.0 * moleFractions_[i] * TKelvin * a_coeff_vec(1,counter) - 3.0 * moleFractions_[i] * a_vec_Curr_[counter];
230 }
231 }
232
234 double fac2 = mv + TKelvin * dpdT_ / dpdV_;
235 for (size_t k = 0; k < m_kk; k++) {
236 double hE_v = (mv * dpdni_[k] - RT() - b_vec_Curr_[k]/ (m_b_current * m_b_current * sqt) * log(vpb/mv)*fac
237 + 1.0 / (m_b_current * sqt) * log(vpb/mv) * m_workS[k]
238 + b_vec_Curr_[k] / vpb / (m_b_current * sqt) * fac);
239 hbar[k] = hbar[k] + hE_v;
240 hbar[k] -= fac2 * dpdni_[k];
241 }
242}
243
245{
246 getEntropy_R_ref(sbar);
247 scale(sbar.begin(), sbar.end(), sbar.begin(), GasConstant);
248 double TKelvin = temperature();
249 double sqt = sqrt(TKelvin);
250 double mv = molarVolume();
251 double refP = refPressure();
252
253 for (size_t k = 0; k < m_kk; k++) {
254 double xx = std::max(SmallNumber, moleFraction(k));
255 sbar[k] += GasConstant * (- log(xx));
256 }
257 for (size_t k = 0; k < m_kk; k++) {
258 m_pp[k] = 0.0;
259 for (size_t i = 0; i < m_kk; i++) {
260 size_t counter = k + m_kk*i;
261 m_pp[k] += moleFractions_[i] * a_vec_Curr_[counter];
262 }
263 }
264 for (size_t k = 0; k < m_kk; k++) {
265 m_workS[k] = 0.0;
266 for (size_t i = 0; i < m_kk; i++) {
267 size_t counter = k + m_kk*i;
268 m_workS[k] += moleFractions_[i] * a_coeff_vec(1,counter);
269 }
270 }
271
272 double dadt = da_dt();
273 double fac = dadt - m_a_current / (2.0 * TKelvin);
274 double vmb = mv - m_b_current;
275 double vpb = mv + m_b_current;
276 for (size_t k = 0; k < m_kk; k++) {
277 sbar[k] -=(GasConstant * log(GasConstant * TKelvin / (refP * mv))
279 + GasConstant * log(mv/vmb)
280 + GasConstant * b_vec_Curr_[k]/vmb
281 + m_pp[k]/(m_b_current * TKelvin * sqt) * log(vpb/mv)
282 - 2.0 * m_workS[k]/(m_b_current * sqt) * log(vpb/mv)
283 + b_vec_Curr_[k] / (m_b_current * m_b_current * sqt) * log(vpb/mv) * fac
284 - 1.0 / (m_b_current * sqt) * b_vec_Curr_[k] / vpb * fac
285 );
286 }
287
289 getPartialMolarVolumes(m_partialMolarVolumes);
290 for (size_t k = 0; k < m_kk; k++) {
291 sbar[k] -= -m_partialMolarVolumes[k] * dpdT_;
292 }
293}
294
296{
297 // u_k = h_k - P * v_k
298 getPartialMolarVolumes(m_partialMolarVolumes);
300 double p = pressure();
301 for (size_t k = 0; k < nSpecies(); k++) {
302 ubar[k] -= p * m_partialMolarVolumes[k];
303 }
304}
305
307{
308 checkArraySize("RedlichKwongMFTP::getPartialMolarIntEnergies_TV", utilde.size(), m_kk);
310
311 double T = temperature();
312 double sqt = sqrt(T);
313 double mv = molarVolume();
314 double b = m_b_current;
315 double a = m_a_current;
316 double dadt = da_dt();
317
318 for (size_t k = 0; k < m_kk; k++) {
319 utilde[k] = RT() * (m_h0_RT[k] - 1.0);
320 }
321
322 if (fabs(b) < SmallNumber) {
323 return;
324 }
325
326 // A_k = sum_i x_i a_ki and A'_k = dA_k/dT = sum_i x_i a_ki,1
327 for (size_t k = 0; k < m_kk; k++) {
328 m_pp[k] = 0.0;
329 m_dAkdT[k] = 0.0;
330 for (size_t i = 0; i < m_kk; i++) {
331 size_t counter = k + m_kk * i;
332 m_pp[k] += moleFractions_[i] * a_vec_Curr_[counter];
333 m_dAkdT[k] += moleFractions_[i] * a_coeff_vec(1, counter);
334 }
335 }
336
337 double vpb = mv + b;
338 double logv = log(vpb / mv);
339 double F = T * dadt - 1.5 * a;
340 double C = -logv / b + 1.0 / vpb;
341 double pref = 1.0 / (b * sqt);
342
343 for (size_t k = 0; k < m_kk; k++) {
344 double Sk = 2.0 * T * m_dAkdT[k] - 3.0 * m_pp[k];
345 double ures = pref * (logv * Sk + b_vec_Curr_[k] * F * C);
346 utilde[k] += ures;
347 }
348}
349
350void RedlichKwongMFTP::getPartialMolarCp(span<double> cpbar) const
351{
352 checkArraySize("RedlichKwongMFTP::getPartialMolarCp", cpbar.size(), m_kk);
354
355 double T = temperature();
356 double sqt = sqrt(T);
357 double v = molarVolume();
358 double b = m_b_current;
359 double a = m_a_current;
360 double dadt = da_dt();
361 double d2adt2 = 0.0; // current RK model uses a(T)=a0+a1*T
362
363 for (size_t k = 0; k < m_kk; k++) {
364 cpbar[k] = GasConstant * m_cp0_R[k];
365 }
366
367 if (fabs(b) < SmallNumber) {
368 return;
369 }
370
371 // Mixture pressure derivatives
373 double pT = dpdT_;
374 double pV = dpdV_;
375 double dvdT_P = -pT / pV;
376
377 double vpb = v + b;
378 double vmb = v - b;
379
380 // P_TT, P_TV and P_VV for v''(T)|P from EOS:
381 // v'' = -(P_TT + 2 P_TV v' + P_VV v'^2) / P_V
382 double pTT = (-d2adt2 + dadt / T - 3.0 * a / (4.0 * T * T)) / (sqt * v * vpb);
383 double pTV = -GasConstant / (vmb * vmb)
384 + (dadt - a / (2.0 * T)) * (2.0 * v + b) / (sqt * v * v * vpb * vpb);
385 double pVV = 2.0 * GasConstant * T / (vmb * vmb * vmb)
386 - 2.0 * a * (3.0 * v * v + 3.0 * b * v + b * b)
387 / (sqt * v * v * v * vpb * vpb * vpb);
388 double d2vdT2_P = -(pTT + 2.0 * pTV * dvdT_P + pVV * dvdT_P * dvdT_P) / pV;
389
390 // Species-dependent A_k and dA_k/dT
391 for (size_t k = 0; k < m_kk; k++) {
392 m_pp[k] = 0.0;
393 m_dAkdT[k] = 0.0;
394 for (size_t i = 0; i < m_kk; i++) {
395 size_t counter = k + m_kk * i;
396 m_pp[k] += moleFractions_[i] * a_vec_Curr_[counter];
397 m_dAkdT[k] += moleFractions_[i] * a_coeff_vec(1, counter);
398 }
399 }
400
401 double F = T * dadt - 1.5 * a;
402 double dF_dT = -0.5 * dadt + T * d2adt2;
403 double logvpv = log(vpb / v);
404 double termLPrime = -b * dvdT_P / (v * vpb);
405
406 for (size_t k = 0; k < m_kk; k++) {
407 double bk = b_vec_Curr_[k];
408 double Ak = m_pp[k];
409 double dAk = m_dAkdT[k];
410 double Sk = 2.0 * T * dAk - 3.0 * Ak;
411 double dSk_dT = -dAk; // for linear a(T)
412
413 // dp/dn_k at constant T,V,n_j
414 double dpdni = RT() / vmb + RT() * bk / (vmb * vmb)
415 - 2.0 * Ak / (sqt * v * vpb)
416 + a * bk / (sqt * v * vpb * vpb);
417
418 // d/dT|P [dp/dn_k]
419 double d_dpdni_dT = GasConstant / vmb
420 - GasConstant * T * dvdT_P / (vmb * vmb)
421 + GasConstant * bk / (vmb * vmb)
422 - 2.0 * GasConstant * T * bk * dvdT_P / (vmb * vmb * vmb)
423 - 2.0 * dAk / (sqt * v * vpb)
424 + Ak / (T * sqt * v * vpb)
425 + 2.0 * Ak * (2.0 * v + b) * dvdT_P / (sqt * v * v * vpb * vpb)
426 + bk * dadt / (sqt * v * vpb * vpb)
427 - bk * a / (2.0 * T * sqt * v * vpb * vpb)
428 - bk * a * (3.0 * v + b) * dvdT_P / (sqt * v * v * vpb * vpb * vpb);
429
430 // d/dT|P of departure terms used in getPartialMolarEnthalpies
431 double dterm1 = -bk / (b * b * sqt)
432 * (termLPrime * F + logvpv * dF_dT - logvpv * F / (2.0 * T));
433 double dterm2 = 1.0 / (b * sqt)
434 * (termLPrime * Sk + logvpv * dSk_dT - logvpv * Sk / (2.0 * T));
435 double dterm3 = bk / (b * sqt)
436 * (dF_dT / vpb - F * dvdT_P / (vpb * vpb) - F / (2.0 * T * vpb));
437
438 // cp_k = d hbar_k / dT at constant P and composition
439 cpbar[k] += -GasConstant
440 + (dvdT_P + T * d2vdT2_P) * dpdni
441 + T * dvdT_P * d_dpdni_dT
442 + dterm1 + dterm2 + dterm3;
443 }
444}
445
446void RedlichKwongMFTP::getPartialMolarCv_TV(span<double> cvtilde) const
447{
448 checkArraySize("RedlichKwongMFTP::getPartialMolarCv_TV", cvtilde.size(), m_kk);
450
451 double T = temperature();
452 double sqt = sqrt(T);
453 double mv = molarVolume();
454 double b = m_b_current;
455 double a = m_a_current;
456 double dadt = da_dt();
457
458 for (size_t k = 0; k < m_kk; k++) {
459 cvtilde[k] = GasConstant * (m_cp0_R[k] - 1.0);
460 }
461
462 if (fabs(b) < SmallNumber) {
463 return;
464 }
465
466 // A_k = sum_i x_i a_ki and A'_k = dA_k/dT = sum_i x_i a_ki,1
467 for (size_t k = 0; k < m_kk; k++) {
468 m_pp[k] = 0.0;
469 m_dAkdT[k] = 0.0;
470 for (size_t i = 0; i < m_kk; i++) {
471 size_t counter = k + m_kk * i;
472 m_pp[k] += moleFractions_[i] * a_vec_Curr_[counter];
473 m_dAkdT[k] += moleFractions_[i] * a_coeff_vec(1, counter);
474 }
475 }
476
477 // Residual contribution for
478 // \tilde{u}_k = (\partial U / \partial n_k)_{T,V,n_j}
479 // with linear a(T): a_ij = a_ij,0 + a_ij,1 T
480 double vpb = mv + b;
481 double logv = log(vpb / mv);
482 double F = T * dadt - 1.5 * a;
483 double C = -logv / b + 1.0 / vpb;
484 double pref = 1.0 / (b * sqt);
485
486 for (size_t k = 0; k < m_kk; k++) {
487 double Sk = 2.0 * T * m_dAkdT[k] - 3.0 * m_pp[k];
488 double ures = pref * (logv * Sk + b_vec_Curr_[k] * F * C);
489 double dresdT = pref * (-logv * m_dAkdT[k] - 0.5 * b_vec_Curr_[k] * dadt * C)
490 - 0.5 * ures / T;
491 cvtilde[k] += dresdT;
492 }
493}
494
495void RedlichKwongMFTP::getPartialMolarVolumes(span<double> vbar) const
496{
497 checkArraySize("RedlichKwongMFTP::getPartialMolarVolumes", vbar.size(), m_kk);
498 for (size_t k = 0; k < m_kk; k++) {
499 m_pp[k] = 0.0;
500 for (size_t i = 0; i < m_kk; i++) {
501 size_t counter = k + m_kk*i;
502 m_pp[k] += moleFractions_[i] * a_vec_Curr_[counter];
503 }
504 }
505 for (size_t k = 0; k < m_kk; k++) {
506 m_workS[k] = 0.0;
507 for (size_t i = 0; i < m_kk; i++) {
508 size_t counter = k + m_kk*i;
509 m_workS[k] += moleFractions_[i] * a_coeff_vec(1,counter);
510 }
511 }
512
513 double sqt = sqrt(temperature());
514 double mv = molarVolume();
515 double vmb = mv - m_b_current;
516 double vpb = mv + m_b_current;
517 for (size_t k = 0; k < m_kk; k++) {
518 double num = (RT() + RT() * m_b_current/ vmb + RT() * b_vec_Curr_[k] / vmb
519 + RT() * m_b_current * b_vec_Curr_[k] /(vmb * vmb)
520 - 2.0 * m_pp[k] / (sqt * vpb)
521 + m_a_current * b_vec_Curr_[k] / (sqt * vpb * vpb)
522 );
523 double denom = (pressure() + RT() * m_b_current/(vmb * vmb) - m_a_current / (sqt * vpb * vpb)
524 );
525 vbar[k] = num / denom;
526 }
527}
528
529bool RedlichKwongMFTP::addSpecies(shared_ptr<Species> spec)
530{
531 bool added = MixtureFugacityTP::addSpecies(spec);
532 if (added) {
533 a_vec_Curr_.resize(m_kk * m_kk, 0.0);
534
535 // Initialize a_vec and b_vec to NaN, to screen for species with
536 // pureFluidParameters which are undefined in the input file:
537 b_vec_Curr_.push_back(NAN);
538 a_coeff_vec.resize(2, m_kk * m_kk, NAN);
539
540 m_pp.push_back(0.0);
541 m_dAkdT.push_back(0.0);
542 m_coeffSource.push_back(CoeffSource::EoS);
543 m_partialMolarVolumes.push_back(0.0);
544 dpdni_.push_back(0.0);
545 }
546 return added;
547}
548
550{
551 // Contents of 'critical-properties.yaml', loaded later if needed
552 AnyMap critPropsDb;
553 std::unordered_map<string, AnyMap*> dbSpecies;
554
555 for (auto& [name, species] : m_species) {
556 auto& data = species->input;
557 size_t k = speciesIndex(name, true);
558 if (!isnan(a_coeff_vec(0, k + m_kk * k))) {
559 continue;
560 }
561 bool foundCoeffs = false;
562
563 if (data.hasKey("equation-of-state") &&
564 data["equation-of-state"].hasMapWhere("model", "Redlich-Kwong"))
565 {
566 // Read a and b coefficients from species 'input' information (that is, as
567 // specified in a YAML input file) specific to the Redlich-Kwong model
568 auto eos = data["equation-of-state"].getMapWhere(
569 "model", "Redlich-Kwong");
570
571 if (eos.hasKey("a") && eos.hasKey("b")) {
572 double a0 = 0, a1 = 0;
573 if (eos["a"].isScalar()) {
574 a0 = eos.convert("a", "Pa*m^6/kmol^2*K^0.5");
575 } else {
576 auto avec = eos["a"].asVector<AnyValue>(2);
577 a0 = eos.units().convert(avec[0], "Pa*m^6/kmol^2*K^0.5");
578 a1 = eos.units().convert(avec[1], "Pa*m^6/kmol^2/K^0.5");
579 }
580 double b = eos.convert("b", "m^3/kmol");
581 foundCoeffs = true;
582 setSpeciesCoeffs(name, a0, a1, b);
583 m_coeffSource[k] = CoeffSource::EoS;
584 }
585
586 if (eos.hasKey("binary-a")) {
587 AnyMap& binary_a = eos["binary-a"].as<AnyMap>();
588 const UnitSystem& units = binary_a.units();
589 for (auto& [name2, coeff] : binary_a) {
590 double a0 = 0, a1 = 0;
591 if (coeff.isScalar()) {
592 a0 = units.convert(coeff, "Pa*m^6/kmol^2*K^0.5");
593 } else {
594 auto avec = coeff.asVector<AnyValue>(2);
595 a0 = units.convert(avec[0], "Pa*m^6/kmol^2*K^0.5");
596 a1 = units.convert(avec[1], "Pa*m^6/kmol^2/K^0.5");
597 }
598 setBinaryCoeffs(name, name2, a0, a1);
599 }
600 }
601
602 if (foundCoeffs) {
603 continue;
604 }
605 }
606
607 // Coefficients have not been populated from model-specific input
608 double Tc = NAN, Pc = NAN;
609 if (data.hasKey("critical-parameters")) {
610 // Use critical state information stored in the species entry to
611 // calculate a and b
612 auto& critProps = data["critical-parameters"].as<AnyMap>();
613 Tc = critProps.convert("critical-temperature", "K");
614 Pc = critProps.convert("critical-pressure", "Pa");
615 m_coeffSource[k] = CoeffSource::CritProps;
616 } else {
617 // Search 'crit-properties.yaml' to find Tc and Pc. Load data if needed
618 if (critPropsDb.empty()) {
619 critPropsDb = AnyMap::fromYamlFile("critical-properties.yaml");
620 dbSpecies = critPropsDb["species"].asMap("name");
621 }
622
623 // All names in critical-properties.yaml are upper case
624 auto ucName = boost::algorithm::to_upper_copy(name);
625 if (dbSpecies.count(ucName)) {
626 auto& spec = *dbSpecies.at(ucName);
627 auto& critProps = spec["critical-parameters"].as<AnyMap>();
628 Tc = critProps.convert("critical-temperature", "K");
629 Pc = critProps.convert("critical-pressure", "Pa");
630 m_coeffSource[k] = CoeffSource::Database;
631 }
632 }
633
634 // Check if critical properties were found in either location
635 if (!isnan(Tc)) {
636 // Assuming no temperature dependence (that is, a1 = 0)
637 double a = omega_a * pow(GasConstant, 2) * pow(Tc, 2.5) / Pc;
638 double b = omega_b * GasConstant * Tc / Pc;
639 setSpeciesCoeffs(name, a, 0.0, b);
640 } else {
641 throw InputFileError("RedlichKwongMFTP::initThermo", data,
642 "No critical property or Redlich-Kwong parameters found "
643 "for species {}.", name);
644 }
645 }
646}
647
649 AnyMap& speciesNode) const
650{
652 size_t k = speciesIndex(name, true);
653 if (m_coeffSource[k] == CoeffSource::EoS) {
654 auto& eosNode = speciesNode["equation-of-state"].getMapWhere(
655 "model", "Redlich-Kwong", true);
656
657 size_t counter = k + m_kk * k;
658 if (a_coeff_vec(1, counter) != 0.0) {
659 vector<AnyValue> coeffs(2);
660 coeffs[0].setQuantity(a_coeff_vec(0, counter), "Pa*m^6/kmol^2*K^0.5");
661 coeffs[1].setQuantity(a_coeff_vec(1, counter), "Pa*m^6/kmol^2/K^0.5");
662 eosNode["a"] = std::move(coeffs);
663 } else {
664 eosNode["a"].setQuantity(a_coeff_vec(0, counter),
665 "Pa*m^6/kmol^2*K^0.5");
666 }
667 eosNode["b"].setQuantity(b_vec_Curr_[k], "m^3/kmol");
668 } else if (m_coeffSource[k] == CoeffSource::CritProps) {
669 auto& critProps = speciesNode["critical-parameters"];
670 double a = a_coeff_vec(0, k + m_kk * k);
671 double b = b_vec_Curr_[k];
672 double Tc = pow(a * omega_b / (b * omega_a * GasConstant), 2.0/3.0);
673 double Pc = omega_b * GasConstant * Tc / b;
674 critProps["critical-temperature"].setQuantity(Tc, "K");
675 critProps["critical-pressure"].setQuantity(Pc, "Pa");
676 }
677
678 if (m_binaryParameters.count(name)) {
679 auto& eosNode = speciesNode["equation-of-state"].getMapWhere(
680 "model", "Redlich-Kwong", true);
681 AnyMap bin_a;
682 for (const auto& [name2, coeffs] : m_binaryParameters.at(name)) {
683 if (coeffs.second == 0) {
684 bin_a[name2].setQuantity(coeffs.first, "Pa*m^6/kmol^2*K^0.5");
685 } else {
686 vector<AnyValue> C(2);
687 C[0].setQuantity(coeffs.first, "Pa*m^6/kmol^2*K^0.5");
688 C[1].setQuantity(coeffs.second, "Pa*m^6/kmol^2/K^0.5");
689 bin_a[name2] = std::move(C);
690 }
691 }
692 eosNode["binary-a"] = std::move(bin_a);
693 }
694}
695
697{
698 // note this agrees with tpx
699 double rho = density();
700 double mmw = meanMolecularWeight();
701 double molarV = mmw / rho;
702 double hh = m_b_current / molarV;
703 double zz = z();
704 double dadt = da_dt();
705 double T = temperature();
706 double sqT = sqrt(T);
707 double fac = dadt - m_a_current / (2.0 * T);
708 double sresid_mol_R = log(zz*(1.0 - hh)) + log(1.0 + hh) * fac / (sqT * GasConstant * m_b_current);
709 return GasConstant * sresid_mol_R;
710}
711
713{
714 // note this agrees with tpx
715 double rho = density();
716 double mmw = meanMolecularWeight();
717 double molarV = mmw / rho;
718 double hh = m_b_current / molarV;
719 double zz = z();
720 double dadt = da_dt();
721 double T = temperature();
722 double sqT = sqrt(T);
723 double fac = T * dadt - 3.0 *m_a_current / (2.0);
724 return GasConstant * T * (zz - 1.0) + fac * log(1.0 + hh) / (sqT * m_b_current);
725}
726
727double RedlichKwongMFTP::liquidVolEst(double TKelvin, double& presGuess) const
728{
729 double v = m_b_current * 1.1;
730 double atmp;
731 double btmp;
732 calculateAB(TKelvin, atmp, btmp);
733 double pres = std::max(psatEst(TKelvin), presGuess);
734 double Vroot[3];
735 bool foundLiq = false;
736 int m = 0;
737 while (m < 100 && !foundLiq) {
738 int nsol = solveCubic(TKelvin, pres, atmp, btmp, Vroot);
739 if (nsol == 1 || nsol == 2) {
740 double pc = critPressure();
741 if (pres > pc) {
742 foundLiq = true;
743 }
744 pres *= 1.04;
745 } else {
746 foundLiq = true;
747 }
748 }
749
750 if (foundLiq) {
751 v = Vroot[0];
752 presGuess = pres;
753 } else {
754 v = -1.0;
755 }
756 return v;
757}
758
759double RedlichKwongMFTP::densityCalc(double TKelvin, double presPa, int phaseRequested,
760 double rhoguess)
761{
762 // It's necessary to set the temperature so that m_a_current is set correctly.
763 setTemperature(TKelvin);
764 double tcrit = critTemperature();
765 double mmw = meanMolecularWeight();
766 if (rhoguess == -1.0) {
767 if (phaseRequested >= FLUID_LIQUID_0) {
768 double lqvol = liquidVolEst(TKelvin, presPa);
769 rhoguess = mmw / lqvol;
770 }
771 } else {
772 // Assume the Gas phase initial guess, if nothing is specified to the routine
773 rhoguess = presPa * mmw / (GasConstant * TKelvin);
774 }
775
776
777 double volguess = mmw / rhoguess;
778 NSolns_ = solveCubic(TKelvin, presPa, m_a_current, m_b_current, Vroot_);
779
780 double molarVolLast = Vroot_[0];
781 if (NSolns_ >= 2) {
782 if (phaseRequested >= FLUID_LIQUID_0) {
783 molarVolLast = Vroot_[0];
784 } else if (phaseRequested == FLUID_GAS || phaseRequested == FLUID_SUPERCRIT) {
785 molarVolLast = Vroot_[2];
786 } else {
787 if (volguess > Vroot_[1]) {
788 molarVolLast = Vroot_[2];
789 } else {
790 molarVolLast = Vroot_[0];
791 }
792 }
793 } else if (NSolns_ == 1) {
794 if (phaseRequested == FLUID_GAS || phaseRequested == FLUID_SUPERCRIT || phaseRequested == FLUID_UNDEFINED) {
795 molarVolLast = Vroot_[0];
796 } else {
797 return -2.0;
798 }
799 } else if (NSolns_ == -1) {
800 if (phaseRequested >= FLUID_LIQUID_0 || phaseRequested == FLUID_UNDEFINED || phaseRequested == FLUID_SUPERCRIT) {
801 molarVolLast = Vroot_[0];
802 } else if (TKelvin > tcrit) {
803 molarVolLast = Vroot_[0];
804 } else {
805 return -2.0;
806 }
807 } else {
808 molarVolLast = Vroot_[0];
809 return -1.0;
810 }
811 return mmw / molarVolLast;
812}
813
814double RedlichKwongMFTP::dpdVCalc(double TKelvin, double molarVol, double& presCalc) const
815{
816 double sqt = sqrt(TKelvin);
817 presCalc = GasConstant * TKelvin / (molarVol - m_b_current)
818 - m_a_current / (sqt * molarVol * (molarVol + m_b_current));
819
820 double vpb = molarVol + m_b_current;
821 double vmb = molarVol - m_b_current;
822 double dpdv = (- GasConstant * TKelvin / (vmb * vmb)
823 + m_a_current * (2 * molarVol + m_b_current) / (sqt * molarVol * molarVol * vpb * vpb));
824 return dpdv;
825}
826
828{
830 return -1 / (molarVolume() * dpdV_);
831}
832
834{
836 return -dpdT_ / (molarVolume() * dpdV_);
837}
838
840{
842 return temperature() * dpdT_ - pressure();
843}
844
846{
848 return molarVolume() * sqrt(-cp_mole() / cv_mole() * dpdV_ / meanMolecularWeight());
849}
850
852{
853 double TKelvin = temperature();
854 double mv = molarVolume();
855 double pres;
856
857 dpdV_ = dpdVCalc(TKelvin, mv, pres);
858 double sqt = sqrt(TKelvin);
859 double vpb = mv + m_b_current;
860 double vmb = mv - m_b_current;
861 double dadt = da_dt();
862 double fac = dadt - m_a_current/(2.0 * TKelvin);
863 dpdT_ = (GasConstant / vmb - fac / (sqt * mv * vpb));
864}
865
867{
868 double temp = temperature();
869 if (m_formTempParam == 1) {
870 for (size_t i = 0; i < m_kk; i++) {
871 for (size_t j = 0; j < m_kk; j++) {
872 size_t counter = i * m_kk + j;
873 a_vec_Curr_[counter] = a_coeff_vec(0,counter) + a_coeff_vec(1,counter) * temp;
874 }
875 }
876 }
877
878 m_b_current = 0.0;
879 m_a_current = 0.0;
880 for (size_t i = 0; i < m_kk; i++) {
881 m_b_current += moleFractions_[i] * b_vec_Curr_[i];
882 for (size_t j = 0; j < m_kk; j++) {
883 m_a_current += a_vec_Curr_[i * m_kk + j] * moleFractions_[i] * moleFractions_[j];
884 }
885 }
886 if (isnan(m_b_current)) {
887 // One or more species do not have specified coefficients.
888 fmt::memory_buffer b;
889 for (size_t k = 0; k < m_kk; k++) {
890 if (isnan(b_vec_Curr_[k])) {
891 if (b.size() > 0) {
892 fmt_append(b, ", {}", speciesName(k));
893 } else {
894 fmt_append(b, "{}", speciesName(k));
895 }
896 }
897 }
898 throw CanteraError("RedlichKwongMFTP::updateMixingExpressions",
899 "Missing Redlich-Kwong coefficients for species: {}", to_string(b));
900 }
901}
902
903void RedlichKwongMFTP::calculateAB(double temp, double& aCalc, double& bCalc) const
904{
905 bCalc = 0.0;
906 aCalc = 0.0;
907 if (m_formTempParam == 1) {
908 for (size_t i = 0; i < m_kk; i++) {
909 bCalc += moleFractions_[i] * b_vec_Curr_[i];
910 for (size_t j = 0; j < m_kk; j++) {
911 size_t counter = i * m_kk + j;
912 double a_vec_Curr = a_coeff_vec(0,counter) + a_coeff_vec(1,counter) * temp;
913 aCalc += a_vec_Curr * moleFractions_[i] * moleFractions_[j];
914 }
915 }
916 } else {
917 for (size_t i = 0; i < m_kk; i++) {
918 bCalc += moleFractions_[i] * b_vec_Curr_[i];
919 for (size_t j = 0; j < m_kk; j++) {
920 size_t counter = i * m_kk + j;
921 double a_vec_Curr = a_coeff_vec(0,counter);
922 aCalc += a_vec_Curr * moleFractions_[i] * moleFractions_[j];
923 }
924 }
925 }
926}
927
928double RedlichKwongMFTP::da_dt() const
929{
930 double dadT = 0.0;
931 if (m_formTempParam == 1) {
932 for (size_t i = 0; i < m_kk; i++) {
933 for (size_t j = 0; j < m_kk; j++) {
934 size_t counter = i * m_kk + j;
935 dadT+= a_coeff_vec(1,counter) * moleFractions_[i] * moleFractions_[j];
936 }
937 }
938 }
939 return dadT;
940}
941
942void RedlichKwongMFTP::calcCriticalConditions(double& pc, double& tc, double& vc) const
943{
944 double a0 = 0.0;
945 double aT = 0.0;
946 for (size_t i = 0; i < m_kk; i++) {
947 for (size_t j = 0; j <m_kk; j++) {
948 size_t counter = i + m_kk * j;
949 a0 += moleFractions_[i] * moleFractions_[j] * a_coeff_vec(0, counter);
950 aT += moleFractions_[i] * moleFractions_[j] * a_coeff_vec(1, counter);
951 }
952 }
953 double a = m_a_current;
954 double b = m_b_current;
955 if (m_formTempParam != 0) {
956 a = a0;
957 }
958 if (b <= 0.0) {
959 tc = 1000000.;
960 pc = 1.0E13;
961 vc = omega_vc * GasConstant * tc / pc;
962 return;
963 }
964 if (a <= 0.0) {
965 tc = 0.0;
966 pc = 0.0;
967 vc = 2.0 * b;
968 return;
969 }
970 double tmp = a * omega_b / (b * omega_a * GasConstant);
971 double pp = 2./3.;
972 double sqrttc, f, dfdt, deltatc;
973
974 if (m_formTempParam == 0) {
975 tc = pow(tmp, pp);
976 } else {
977 tc = pow(tmp, pp);
978 for (int j = 0; j < 10; j++) {
979 sqrttc = sqrt(tc);
980 f = omega_a * b * GasConstant * tc * sqrttc / omega_b - aT * tc - a0;
981 dfdt = 1.5 * omega_a * b * GasConstant * sqrttc / omega_b - aT;
982 deltatc = - f / dfdt;
983 tc += deltatc;
984 }
985 if (deltatc > 0.1) {
986 throw CanteraError("RedlichKwongMFTP::calcCriticalConditions",
987 "didn't converge");
988 }
989 }
990
991 pc = omega_b * GasConstant * tc / b;
992 vc = omega_vc * GasConstant * tc / pc;
993}
994
995int RedlichKwongMFTP::solveCubic(double T, double pres, double a, double b,
996 span<double> Vroot) const
997{
998
999 // Derive the coefficients of the cubic polynomial to solve.
1000 double an = 1.0;
1001 double bn = - GasConstant * T / pres;
1002 double sqt = sqrt(T);
1003 double cn = - (GasConstant * T * b / pres - a/(pres * sqt) + b * b);
1004 double dn = - (a * b / (pres * sqt));
1005
1006 double tmp = a * omega_b / (b * omega_a * GasConstant);
1007 double pp = 2./3.;
1008 double tc = pow(tmp, pp);
1009 double pc = omega_b * GasConstant * tc / b;
1010 double vc = omega_vc * GasConstant * tc / pc;
1011
1012 int nSolnValues = MixtureFugacityTP::solveCubic(T, pres, a, b, a, Vroot, an, bn, cn, dn, tc, vc);
1013
1014 return nSolnValues;
1015}
1016
1017}
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 getGibbs_ref(span< double > g) const override
Returns the vector of the Gibbs function of the reference state at the current temperature of the sol...
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.
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.
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...