Cantera  4.0.0a1
Loading...
Searching...
No Matches
PengRobinson.cpp
Go to the documentation of this file.
1//! @file PengRobinson.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 PengRobinson::omega_a = 4.5723552892138218E-01;
19const double PengRobinson::omega_b = 7.77960739038885E-02;
20const double PengRobinson::omega_vc = 3.07401308698703833E-01;
21
22PengRobinson::PengRobinson(const string& infile, const string& id_)
23{
24 initThermoFile(infile, id_);
25}
26
27void PengRobinson::setSpeciesCoeffs(const string& species, double a, double b, double w)
28{
29 size_t k = speciesIndex(species, true);
30
31 // Calculate value of kappa (independent of temperature)
32 // w is an acentric factor of species
33 if (w <= 0.491) {
34 m_kappa[k] = 0.37464 + 1.54226*w - 0.26992*w*w;
35 } else {
36 m_kappa[k] = 0.374642 + 1.487503*w - 0.164423*w*w + 0.016666*w*w*w;
37 }
38 m_acentric[k] = w; // store the original acentric factor to enable serialization
39
40 // Calculate alpha (temperature dependent interaction parameter)
41 double critTemp = speciesCritTemperature(a, b);
42 double sqt_T_r = sqrt(temperature() / critTemp);
43 double sqt_alpha = 1 + m_kappa[k] * (1 - sqt_T_r);
44 m_alpha[k] = sqt_alpha*sqt_alpha;
45
46 m_a_coeffs(k,k) = a;
47 double aAlpha_k = a*m_alpha[k];
48 m_aAlpha_binary(k,k) = aAlpha_k;
49
50 // standard mixing rule for cross-species interaction term
51 for (size_t j = 0; j < m_kk; j++) {
52 if (k == j) {
53 continue;
54 }
55 double a0kj = sqrt(m_a_coeffs(j,j) * a);
56 double aAlpha_j = a*m_alpha[j];
57 double a_Alpha = sqrt(aAlpha_j*aAlpha_k);
58 if (m_a_coeffs(j, k) == 0) {
59 m_a_coeffs(j, k) = a0kj;
60 m_aAlpha_binary(j, k) = a_Alpha;
61 m_a_coeffs(k, j) = a0kj;
62 m_aAlpha_binary(k, j) = a_Alpha;
63 }
64 }
65 m_b_coeffs[k] = b;
66}
67
68void PengRobinson::setBinaryCoeffs(const string& species_i,
69 const string& species_j, double a0)
70{
71 size_t ki = speciesIndex(species_i, true);
72 size_t kj = speciesIndex(species_j, true);
73
74 m_a_coeffs(ki, kj) = m_a_coeffs(kj, ki) = a0;
75 m_binaryParameters[species_i][species_j] = a0;
76 m_binaryParameters[species_j][species_i] = a0;
77 // Calculate alpha_ij
78 double alpha_ij = m_alpha[ki] * m_alpha[kj];
79 m_aAlpha_binary(ki, kj) = m_aAlpha_binary(kj, ki) = a0*alpha_ij;
80}
81
82// ------------Molar Thermodynamic Properties -------------------------
83
85{
87 double cpref = GasConstant * mean_X(m_cp0_R);
88 if (m_b == 0.0) {
89 return cpref;
90 }
91 double T = temperature();
92 double mv = molarVolume();
93 double vpb = mv + (1 + Sqrt2) * m_b;
94 double vmb = mv + (1 - Sqrt2) * m_b;
96 double dHdT_V = cpref + mv * m_dpdT - GasConstant
97 + 1.0 / (2.0 * Sqrt2 * m_b) * log(vpb / vmb) * T * d2aAlpha_dT2();
98 return dHdT_V - (mv + T * m_dpdT / m_dpdV) * m_dpdT;
99}
100
102{
104 double cvref = GasConstant * (mean_X(m_cp0_R) - 1.0);
105 if (m_b == 0.0) {
106 return cvref;
107 }
108 double T = temperature();
110 return (cp_mole() + T * m_dpdT * m_dpdT / m_dpdV);
111}
112
114{
116 // Get a copy of the private variables stored in the State object
117 double T = temperature();
118 double mv = molarVolume();
119 double denom = mv * mv + 2 * mv * m_b - m_b * m_b;
120 return GasConstant * T / (mv - m_b) - m_aAlpha_mix / denom;
121}
122
124{
126 return 1.0 / m_workS[k];
127}
128
129void PengRobinson::getActivityCoefficients(span<double> ac) const
130{
131 checkArraySize("PengRobinson::getActivityCoefficients", ac.size(), m_kk);
132 for (size_t k = 0; k < m_kk; k++) {
133 ac[k] = 1.0;
134 }
135 if (m_b == 0.0) {
136 return;
137 }
138 double mv = molarVolume();
139 double vpb2 = mv + (1 + Sqrt2) * m_b;
140 double vmb2 = mv + (1 - Sqrt2) * m_b;
141 double vmb = mv - m_b;
142 double pres = pressure();
143
144 for (size_t k = 0; k < m_kk; k++) {
145 m_pp[k] = 0.0;
146 for (size_t i = 0; i < m_kk; i++) {
147 m_pp[k] += moleFractions_[i] * m_aAlpha_binary(k, i);
148 }
149 }
150 double num = 0;
151 double denom = 2 * Sqrt2 * m_b * m_b;
152 double denom2 = m_b * (mv * mv + 2 * mv * m_b - m_b * m_b);
153 double RT_ = RT();
154 for (size_t k = 0; k < m_kk; k++) {
155 num = 2 * m_b * m_pp[k] - m_aAlpha_mix * m_b_coeffs[k];
156 ac[k] = (-RT_ * log(pres * mv/ RT_) + RT_ * log(mv / vmb)
157 + RT_ * m_b_coeffs[k] / vmb
158 - (num /denom) * log(vpb2/vmb2)
159 - m_aAlpha_mix * m_b_coeffs[k] * mv/denom2
160 );
161 }
162 for (size_t k = 0; k < m_kk; k++) {
163 ac[k] = exp(ac[k]/ RT_);
164 }
165}
166
167// ---- Partial Molar Properties of the Solution -----------------
168
169void PengRobinson::getChemPotentials(span<double> mu) const
170{
171 double RT_ = RT();
173 for (size_t k = 0; k < m_kk; k++) {
174 double xx = std::max(SmallNumber, moleFraction(k));
175 mu[k] += RT_ * log(xx);
176 }
177 if (m_b == 0.0) {
178 return;
179 }
180
181 double mv = molarVolume();
182 double vmb = mv - m_b;
183 double vpb2 = mv + (1 + Sqrt2) * m_b;
184 double vmb2 = mv + (1 - Sqrt2) * m_b;
185
186 for (size_t k = 0; k < m_kk; k++) {
187 m_pp[k] = 0.0;
188 for (size_t i = 0; i < m_kk; i++) {
189 m_pp[k] += moleFractions_[i] * m_aAlpha_binary(k, i);
190 }
191 }
192 double pres = pressure();
193 double denom = 2 * Sqrt2 * m_b * m_b;
194 double denom2 = m_b * (mv * mv + 2 * mv * m_b - m_b * m_b);
195
196 for (size_t k = 0; k < m_kk; k++) {
197 double num = 2 * m_b * m_pp[k] - m_aAlpha_mix * m_b_coeffs[k];
198
199 mu[k] += - RT_ * log(pres * mv / RT_)
200 + RT_ * log(mv / vmb) + RT_ * m_b_coeffs[k] / vmb
201 - (num /denom) * log(vpb2/vmb2)
202 - m_aAlpha_mix * m_b_coeffs[k] * mv/denom2;
203 }
204}
205
206void PengRobinson::getPartialMolarEnthalpies(span<double> hbar) const
207{
208 // First we get the reference state contributions
209 getEnthalpy_RT_ref(hbar);
210 scale(hbar.begin(), hbar.end(), hbar.begin(), RT());
211 if (m_b == 0.0) {
212 return;
213 }
214 vector<double> tmp;
215 tmp.resize(m_kk,0.0);
216
217 // We calculate m_dpdni
218 double T = temperature();
219 double mv = molarVolume();
220 double vmb = mv - m_b;
221 double vpb2 = mv + (1 + Sqrt2) * m_b;
222 double vmb2 = mv + (1 - Sqrt2) * m_b;
223 double daAlphadT = daAlpha_dT();
224
225 for (size_t k = 0; k < m_kk; k++) {
226 m_pp[k] = 0.0;
227 tmp[k] = 0.0;
228 for (size_t i = 0; i < m_kk; i++) {
229 double grad_aAlpha = m_dalphadT[i]/m_alpha[i] + m_dalphadT[k]/m_alpha[k];
230 m_pp[k] += moleFractions_[i] * m_aAlpha_binary(k, i);
231 tmp[k] +=moleFractions_[i] * m_aAlpha_binary(k, i) * grad_aAlpha;
232 }
233 }
234
235 double denom = mv * mv + 2 * mv * m_b - m_b * m_b;
236 double denom2 = denom * denom;
237 double RT_ = RT();
238 for (size_t k = 0; k < m_kk; k++) {
239 m_dpdni[k] = RT_ / vmb + RT_ * m_b_coeffs[k] / (vmb * vmb)
240 - 2.0 * m_pp[k] / denom + 2 * vmb * m_aAlpha_mix * m_b_coeffs[k] / denom2;
241 }
242
243 double fac = T * daAlphadT - m_aAlpha_mix;
245 double fac2 = mv + T * m_dpdT / m_dpdV;
246 double fac3 = 2 * Sqrt2 * m_b * m_b;
247 double fac4 = 0;
248 for (size_t k = 0; k < m_kk; k++) {
249 fac4 = T*tmp[k] -2 * m_pp[k];
250 double hE_v = mv * m_dpdni[k] - RT_
251 - m_b_coeffs[k] / fac3 * log(vpb2 / vmb2) * fac
252 + (mv * m_b_coeffs[k]) / (m_b * denom) * fac
253 + 1 / (2 * Sqrt2 * m_b) * log(vpb2 / vmb2) * fac4;
254 hbar[k] = hbar[k] + hE_v;
255 hbar[k] -= fac2 * m_dpdni[k];
256 }
257}
258
259void PengRobinson::getPartialMolarEntropies(span<double> sbar) const
260{
261 // Using the identity : (hk - T*sk) = gk
262 double T = temperature();
265 for (size_t k = 0; k < m_kk; k++) {
266 sbar[k] = (sbar[k] - m_workS[k])/T;
267 }
268}
269
270void PengRobinson::getPartialMolarIntEnergies(span<double> ubar) const
271{
272 // u_i = h_i - p*v_i
273 double p = pressure();
276 for (size_t k = 0; k < m_kk; k++) {
277 ubar[k] = ubar[k] - p*m_workS[k];
278 }
279}
280
281void PengRobinson::getPartialMolarCp(span<double> cpbar) const
282{
283 throw NotImplementedError("PengRobinson::getPartialMolarCp");
284}
285
286void PengRobinson::getPartialMolarVolumes(span<double> vbar) const
287{
288 checkArraySize("PengRobinson::getPartialMolarVolumes", vbar.size(), m_kk);
289 for (size_t k = 0; k < m_kk; k++) {
290 m_pp[k] = 0.0;
291 for (size_t i = 0; i < m_kk; i++) {
292 m_pp[k] += moleFractions_[i] * m_aAlpha_binary(k, i);
293 }
294 }
295
296 double mv = molarVolume();
297 double vmb = mv - m_b;
298 double vpb = mv + m_b;
299 double fac = mv * mv + 2 * mv * m_b - m_b * m_b;
300 double fac2 = fac * fac;
301 double RT_ = RT();
302
303 for (size_t k = 0; k < m_kk; k++) {
304 double num = RT_ + RT_ * m_b/ vmb + RT_ * m_b_coeffs[k] / vmb
305 + RT_ * m_b * m_b_coeffs[k] /(vmb * vmb) - 2 * mv * m_pp[k] / fac
306 + 2 * mv * vmb * m_aAlpha_mix * m_b_coeffs[k] / fac2;
307 double denom = pressure() + RT_ * m_b / (vmb * vmb) + m_aAlpha_mix/fac
308 - 2 * mv* vpb * m_aAlpha_mix / fac2;
309 vbar[k] = num / denom;
310 }
311}
312
313double PengRobinson::speciesCritTemperature(double a, double b) const
314{
315 if (b <= 0.0) {
316 return 1000000.;
317 } else if (a <= 0.0) {
318 return 0.0;
319 } else {
320 return a * omega_b / (b * omega_a * GasConstant);
321 }
322}
323
324bool PengRobinson::addSpecies(shared_ptr<Species> spec)
325{
326 bool added = MixtureFugacityTP::addSpecies(spec);
327 if (added) {
328 m_a_coeffs.resize(m_kk, m_kk, 0.0);
329 m_b_coeffs.push_back(0.0);
330 m_aAlpha_binary.resize(m_kk, m_kk, 0.0);
331 m_kappa.push_back(0.0);
332 m_acentric.push_back(0.0);
333 m_alpha.push_back(0.0);
334 m_dalphadT.push_back(0.0);
335 m_d2alphadT2.push_back(0.0);
336 m_pp.push_back(0.0);
337 m_partialMolarVolumes.push_back(0.0);
338 m_dpdni.push_back(0.0);
339 m_coeffSource.push_back(CoeffSource::EoS);
340 }
341 return added;
342}
343
345{
346 // Contents of 'critical-properties.yaml', loaded later if needed
347 AnyMap critPropsDb;
348 std::unordered_map<string, AnyMap*> dbSpecies;
349
350 for (auto& [name, species] : m_species) {
351 auto& data = species->input;
352 size_t k = speciesIndex(name, true);
353 if (m_a_coeffs(k, k) != 0.0) {
354 continue;
355 }
356 bool foundCoeffs = false;
357 if (data.hasKey("equation-of-state") &&
358 data["equation-of-state"].hasMapWhere("model", "Peng-Robinson"))
359 {
360 // Read a and b coefficients and acentric factor w_ac from species input
361 // information, specified in a YAML input file.
362 auto eos = data["equation-of-state"].getMapWhere(
363 "model", "Peng-Robinson");
364 if (eos.hasKey("a") && eos.hasKey("b") && eos.hasKey("acentric-factor")) {
365 double a0 = eos.convert("a", "Pa*m^6/kmol^2");
366 double b = eos.convert("b", "m^3/kmol");
367 // unitless acentric factor:
368 double w = eos["acentric-factor"].asDouble();
369 setSpeciesCoeffs(name, a0, b, w);
370 foundCoeffs = true;
371 }
372
373 if (eos.hasKey("binary-a")) {
374 AnyMap& binary_a = eos["binary-a"].as<AnyMap>();
375 const UnitSystem& units = binary_a.units();
376 for (auto& [name2, coeff] : binary_a) {
377 double a0 = units.convert(coeff, "Pa*m^6/kmol^2");
378 setBinaryCoeffs(name, name2, a0);
379 }
380 }
381 if (foundCoeffs) {
382 m_coeffSource[k] = CoeffSource::EoS;
383 continue;
384 }
385 }
386
387 // Coefficients have not been populated from model-specific input
388 double Tc = NAN, Pc = NAN, omega_ac = NAN;
389 if (data.hasKey("critical-parameters")) {
390 // Use critical state information stored in the species entry to
391 // calculate a, b, and the acentric factor.
392 auto& critProps = data["critical-parameters"].as<AnyMap>();
393 Tc = critProps.convert("critical-temperature", "K");
394 Pc = critProps.convert("critical-pressure", "Pa");
395 omega_ac = critProps["acentric-factor"].asDouble();
396 m_coeffSource[k] = CoeffSource::CritProps;
397 } else {
398 // Search 'crit-properties.yaml' to find Tc and Pc. Load data if needed.
399 if (critPropsDb.empty()) {
400 critPropsDb = AnyMap::fromYamlFile("critical-properties.yaml");
401 dbSpecies = critPropsDb["species"].asMap("name");
402 }
403
404 // All names in critical-properties.yaml are upper case
405 auto ucName = boost::algorithm::to_upper_copy(name);
406 if (dbSpecies.count(ucName)) {
407 auto& spec = *dbSpecies.at(ucName);
408 auto& critProps = spec["critical-parameters"].as<AnyMap>();
409 Tc = critProps.convert("critical-temperature", "K");
410 Pc = critProps.convert("critical-pressure", "Pa");
411 omega_ac = critProps["acentric-factor"].asDouble();
412 m_coeffSource[k] = CoeffSource::Database;
413 }
414 }
415
416 // Check if critical properties were found in either location
417 if (!isnan(Tc)) {
418 double a = omega_a * std::pow(GasConstant * Tc, 2) / Pc;
419 double b = omega_b * GasConstant * Tc / Pc;
420 setSpeciesCoeffs(name, a, b, omega_ac);
421 } else {
422 throw InputFileError("PengRobinson::initThermo", data,
423 "No Peng-Robinson model parameters or critical properties found for "
424 "species '{}'", name);
425 }
426 }
427}
428
429void PengRobinson::getSpeciesParameters(const string& name, AnyMap& speciesNode) const
430{
432 size_t k = speciesIndex(name, true);
433
434 // Pure species parameters
435 if (m_coeffSource[k] == CoeffSource::EoS) {
436 auto& eosNode = speciesNode["equation-of-state"].getMapWhere(
437 "model", "Peng-Robinson", true);
438 eosNode["a"].setQuantity(m_a_coeffs(k, k), "Pa*m^6/kmol^2");
439 eosNode["b"].setQuantity(m_b_coeffs[k], "m^3/kmol");
440 eosNode["acentric-factor"] = m_acentric[k];
441 } else if (m_coeffSource[k] == CoeffSource::CritProps) {
442 auto& critProps = speciesNode["critical-parameters"];
443 double Tc = speciesCritTemperature(m_a_coeffs(k, k), m_b_coeffs[k]);
444 double Pc = omega_b * GasConstant * Tc / m_b_coeffs[k];
445 critProps["critical-temperature"].setQuantity(Tc, "K");
446 critProps["critical-pressure"].setQuantity(Pc, "Pa");
447 critProps["acentric-factor"] = m_acentric[k];
448 }
449 // Nothing to do in the case where the parameters are from the database
450
451 if (m_binaryParameters.count(name)) {
452 // Include binary parameters regardless of where the pure species parameters
453 // were found
454 auto& eosNode = speciesNode["equation-of-state"].getMapWhere(
455 "model", "Peng-Robinson", true);
456 AnyMap bin_a;
457 for (const auto& [species, coeff] : m_binaryParameters.at(name)) {
458 bin_a[species].setQuantity(coeff, "Pa*m^6/kmol^2");
459 }
460 eosNode["binary-a"] = std::move(bin_a);
461 }
462}
463
465{
466 if (m_b == 0.0) {
467 return 0.0;
468 }
469 double molarV = molarVolume();
470 double hh = m_b / molarV;
471 double zz = z();
472 double alpha_1 = daAlpha_dT();
473 double vpb = molarV + (1.0 + Sqrt2) * m_b;
474 double vmb = molarV + (1.0 - Sqrt2) * m_b;
475 double fac = alpha_1 / (2.0 * Sqrt2 * m_b);
476 double sresid_mol_R = log(zz*(1.0 - hh)) + fac * log(vpb / vmb) / GasConstant;
477 return GasConstant * sresid_mol_R;
478}
479
481{
482 if (m_b == 0.0) {
483 return 0.0;
484 }
485 double molarV = molarVolume();
486 double zz = z();
487 double aAlpha_1 = daAlpha_dT();
488 double T = temperature();
489 double vpb = molarV + (1 + Sqrt2) * m_b;
490 double vmb = molarV + (1 - Sqrt2) * m_b;
491 double fac = 1 / (2.0 * Sqrt2 * m_b);
492 return GasConstant * T * (zz - 1.0)
493 + fac * log(vpb / vmb) * (T * aAlpha_1 - m_aAlpha_mix);
494}
495
496double PengRobinson::liquidVolEst(double T, double& presGuess) const
497{
498 double v = m_b * 1.1;
499 double atmp, btmp, aAlphatmp;
500 calculateAB(atmp, btmp, aAlphatmp);
501 double pres = std::max(psatEst(T), presGuess);
502 double Vroot[3];
503 bool foundLiq = false;
504 int m = 0;
505 while (m < 100 && !foundLiq) {
506 int nsol = solveCubic(T, pres, atmp, btmp, aAlphatmp, Vroot);
507 if (nsol == 1 || nsol == 2) {
508 double pc = critPressure();
509 if (pres > pc) {
510 foundLiq = true;
511 }
512 pres *= 1.04;
513 } else {
514 foundLiq = true;
515 }
516 }
517
518 if (foundLiq) {
519 v = Vroot[0];
520 presGuess = pres;
521 } else {
522 v = -1.0;
523 }
524 return v;
525}
526
527double PengRobinson::densityCalc(double T, double presPa, int phaseRequested,
528 double rhoGuess)
529{
530 // It's necessary to set the temperature so that m_aAlpha_mix is set correctly.
532 double tcrit = critTemperature();
533 double mmw = meanMolecularWeight();
534 if (rhoGuess == -1.0) {
535 if (phaseRequested >= FLUID_LIQUID_0) {
536 double lqvol = liquidVolEst(T, presPa);
537 rhoGuess = mmw / lqvol;
538 }
539 } else {
540 // Assume the Gas phase initial guess, if nothing is specified to the routine
541 rhoGuess = presPa * mmw / (GasConstant * T);
542 }
543
544 double volGuess = mmw / rhoGuess;
545 m_NSolns = solveCubic(T, presPa, m_a, m_b, m_aAlpha_mix, m_Vroot);
546
547 // Ensure root ordering used for branch selection only contains physical roots.
548 // This avoids selecting negative-volume roots in mixed-parameter mechanisms.
549 double vmin = std::max(0.0, m_b * (1.0 + 1e-12));
550 vector<double> physicalRoots;
551 for (double root : m_Vroot) {
552 if (std::isfinite(root) && root > vmin) {
553 physicalRoots.push_back(root);
554 }
555 }
556 std::sort(physicalRoots.begin(), physicalRoots.end());
557 if (physicalRoots.empty()) {
558 return -1.0;
559 } else if (physicalRoots.size() == 1) {
560 // Preserve branch-request semantics for single-root states:
561 // return -2 when the requested phase is inconsistent with the
562 // single branch identified by solveCubic() sign convention.
563 if ((phaseRequested == FLUID_GAS && m_NSolns < 0)
564 || (phaseRequested >= FLUID_LIQUID_0 && m_NSolns > 0))
565 {
566 return -2.0;
567 }
568 // Otherwise, accept the physically admissible root directly.
569 return mmw / physicalRoots[0];
570 } else if (physicalRoots.size() == 2) {
571 m_Vroot[0] = physicalRoots[0];
572 m_Vroot[1] = 0.5 * (physicalRoots[0] + physicalRoots[1]);
573 m_Vroot[2] = physicalRoots[1];
574 } else {
575 m_Vroot[0] = physicalRoots[0];
576 m_Vroot[1] = physicalRoots[1];
577 m_Vroot[2] = physicalRoots[2];
578 }
579
580 double molarVolLast = m_Vroot[0];
581 if (m_NSolns >= 2) {
582 if (phaseRequested >= FLUID_LIQUID_0) {
583 molarVolLast = m_Vroot[0];
584 } else if (phaseRequested == FLUID_GAS || phaseRequested == FLUID_SUPERCRIT) {
585 molarVolLast = m_Vroot[2];
586 } else {
587 if (volGuess > m_Vroot[1]) {
588 molarVolLast = m_Vroot[2];
589 } else {
590 molarVolLast = m_Vroot[0];
591 }
592 }
593 } else if (m_NSolns == 1) {
594 if (phaseRequested == FLUID_GAS || phaseRequested == FLUID_SUPERCRIT
595 || phaseRequested == FLUID_UNDEFINED)
596 {
597 molarVolLast = m_Vroot[0];
598 } else {
599 return -2.0;
600 }
601 } else if (m_NSolns == -1) {
602 if (phaseRequested >= FLUID_LIQUID_0 || phaseRequested == FLUID_UNDEFINED
603 || phaseRequested == FLUID_SUPERCRIT)
604 {
605 molarVolLast = m_Vroot[0];
606 } else if (T > tcrit) {
607 molarVolLast = m_Vroot[0];
608 } else {
609 return -2.0;
610 }
611 } else {
612 molarVolLast = m_Vroot[0];
613 return -1.0;
614 }
615 return mmw / molarVolLast;
616}
617
618double PengRobinson::dpdVCalc(double T, double molarVol, double& presCalc) const
619{
620 double denom = molarVol * molarVol + 2 * molarVol * m_b - m_b * m_b;
621 double vpb = molarVol + m_b;
622 double vmb = molarVol - m_b;
623 return -GasConstant * T / (vmb * vmb) + 2 * m_aAlpha_mix * vpb / (denom*denom);
624}
625
627{
629 return -1 / (molarVolume() * m_dpdV);
630}
631
633{
635 return -m_dpdT / (molarVolume() * m_dpdV);
636}
637
639{
641 return molarVolume() * sqrt(-cp_mole() / cv_mole() * m_dpdV / meanMolecularWeight());
642}
643
645{
646 double T = temperature();
647 double mv = molarVolume();
648 double pres;
649
650 m_dpdV = dpdVCalc(T, mv, pres);
651 double vmb = mv - m_b;
652 double denom = mv * mv + 2 * mv * m_b - m_b * m_b;
653 m_dpdT = (GasConstant / vmb - daAlpha_dT() / denom);
654}
655
657{
658 double temp = temperature();
659
660 // Update individual alpha
661 for (size_t j = 0; j < m_kk; j++) {
662 double critTemp_j = speciesCritTemperature(m_a_coeffs(j,j), m_b_coeffs[j]);
663 double sqt_alpha = 1 + m_kappa[j] * (1 - sqrt(temp / critTemp_j));
664 m_alpha[j] = sqt_alpha*sqt_alpha;
665 }
666
667 // Update aAlpha_i, j
668 for (size_t i = 0; i < m_kk; i++) {
669 for (size_t j = 0; j < m_kk; j++) {
670 m_aAlpha_binary(i, j) = sqrt(m_alpha[i] * m_alpha[j]) * m_a_coeffs(i,j);
671 }
672 }
674}
675
676void PengRobinson::calculateAB(double& aCalc, double& bCalc, double& aAlphaCalc) const
677{
678 bCalc = 0.0;
679 aCalc = 0.0;
680 aAlphaCalc = 0.0;
681 for (size_t i = 0; i < m_kk; i++) {
682 bCalc += moleFractions_[i] * m_b_coeffs[i];
683 for (size_t j = 0; j < m_kk; j++) {
684 aCalc += m_a_coeffs(i, j) * moleFractions_[i] * moleFractions_[j];
685 aAlphaCalc += m_aAlpha_binary(i, j) * moleFractions_[i] * moleFractions_[j];
686 }
687 }
688}
689
691{
692 double daAlphadT = 0.0, k, Tc, sqtTr, coeff1, coeff2;
693 for (size_t i = 0; i < m_kk; i++) {
694 // Calculate first derivative of alpha for individual species
695 Tc = speciesCritTemperature(m_a_coeffs(i,i), m_b_coeffs[i]);
696 sqtTr = sqrt(temperature() / Tc);
697 coeff1 = 1 / (Tc*sqtTr);
698 coeff2 = sqtTr - 1;
699 k = m_kappa[i];
700 m_dalphadT[i] = coeff1 * (k*k*coeff2 - k);
701 }
702 // Calculate mixture derivative
703 for (size_t i = 0; i < m_kk; i++) {
704 for (size_t j = 0; j < m_kk; j++) {
705 daAlphadT += moleFractions_[i] * moleFractions_[j] * 0.5
706 * m_aAlpha_binary(i, j)
707 * (m_dalphadT[i] / m_alpha[i] + m_dalphadT[j] / m_alpha[j]);
708 }
709 }
710 return daAlphadT;
711}
712
714{
715 for (size_t i = 0; i < m_kk; i++) {
716 double Tcrit_i = speciesCritTemperature(m_a_coeffs(i, i), m_b_coeffs[i]);
717 double sqt_Tr = sqrt(temperature() / Tcrit_i);
718 double coeff1 = 1 / (Tcrit_i*sqt_Tr);
719 double coeff2 = sqt_Tr - 1;
720 // Calculate first and second derivatives of alpha for individual species
721 double k = m_kappa[i];
722 m_dalphadT[i] = coeff1 * (k*k*coeff2 - k);
723 m_d2alphadT2[i] = (k*k + k) * coeff1 / (2*sqt_Tr*sqt_Tr*Tcrit_i);
724 }
725
726 // Calculate mixture derivative
727 double d2aAlphadT2 = 0.0;
728 for (size_t i = 0; i < m_kk; i++) {
729 double alphai = m_alpha[i];
730 for (size_t j = 0; j < m_kk; j++) {
731 double alphaj = m_alpha[j];
732 double alphaij = alphai * alphaj;
733 double term1 = m_d2alphadT2[i] / alphai + m_d2alphadT2[j] / alphaj;
734 double term2 = 2 * m_dalphadT[i] * m_dalphadT[j] / alphaij;
735 double term3 = m_dalphadT[i] / alphai + m_dalphadT[j] / alphaj;
736 d2aAlphadT2 += 0.5 * moleFractions_[i] * moleFractions_[j]
737 * m_aAlpha_binary(i, j)
738 * (term1 + term2 - 0.5 * term3 * term3);
739 }
740 }
741 return d2aAlphadT2;
742}
743
744void PengRobinson::calcCriticalConditions(double& pc, double& tc, double& vc) const
745{
746 if (m_b <= 0.0) {
747 tc = 1000000.;
748 pc = 1.0E13;
749 vc = omega_vc * GasConstant * tc / pc;
750 return;
751 }
752 if (m_a <= 0.0) {
753 tc = 0.0;
754 pc = 0.0;
755 vc = 2.0 * m_b;
756 return;
757 }
758 tc = m_a * omega_b / (m_b * omega_a * GasConstant);
759 pc = omega_b * GasConstant * tc / m_b;
760 vc = omega_vc * GasConstant * tc / pc;
761}
762
763int PengRobinson::solveCubic(double T, double pres, double a, double b, double aAlpha,
764 span<double> Vroot) const
765{
766 // Derive the coefficients of the cubic polynomial (in terms of molar volume v)
767 double bsqr = b * b;
768 double RT_p = GasConstant * T / pres;
769 double aAlpha_p = aAlpha / pres;
770 double an = 1.0;
771 double bn = (b - RT_p);
772 double cn = -(2 * RT_p * b - aAlpha_p + 3 * bsqr);
773 double dn = (bsqr * RT_p + bsqr * b - aAlpha_p * b);
774
775 double tc = a * omega_b / (b * omega_a * GasConstant);
776 double pc = omega_b * GasConstant * tc / b;
777 double vc = omega_vc * GasConstant * tc / pc;
778
779 return MixtureFugacityTP::solveCubic(T, pres, a, b, aAlpha, Vroot,
780 an, bn, cn, dn, tc, vc);
781}
782
784{
786
787 parameters["aAlpha_mix"] = m_aAlpha_mix;
788 parameters["b_mix"] = m_b;
789 parameters["daAlpha_dT"] = daAlpha_dT();
790 parameters["d2aAlpha_dT2"] = d2aAlpha_dT2();
791
792 return parameters;
793}
794
795}
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
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
Error thrown for problems processing information contained in an AnyMap or AnyValue.
Definition AnyMap.h:749
double critPressure() const override
Critical pressure (Pa).
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.
An error indicating that an unimplemented function has been called.
double thermalExpansionCoeff() const override
Return the volumetric thermal expansion coefficient. Units: 1/K.
void setSpeciesCoeffs(const string &species, double a, double b, double w)
Set the pure fluid interaction parameters for a species.
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
Returns an array of partial molar enthalpies for the species in the mixture.
double soundSpeed() const override
Return the speed of sound. Units: m/s.
double pressure() const override
Return the thermodynamic pressure (Pa).
map< string, map< string, double > > m_binaryParameters
Explicitly-specified binary interaction parameters, to enable serialization.
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...
void getPartialMolarCp(span< double > cpbar) const override
Calculate species-specific molar specific heats.
static const double omega_b
Omega constant: b0 (= omega_b) used in Peng-Robinson equation of state.
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 setBinaryCoeffs(const string &species_i, const string &species_j, double a)
Set values for the interaction parameter between two species.
vector< double > m_dpdni
Vector of derivatives of pressure with respect to mole number.
PengRobinson(const string &infile="", const string &id="")
Construct and initialize a PengRobinson object directly from an input file.
void initThermo() override
Initialize the ThermoPhase object after all species have been set up.
double m_dpdT
The derivative of the pressure with respect to the temperature.
double cv_mole() const override
Molar heat capacity at constant volume and composition [J/kmol/K].
double liquidVolEst(double T, double &pres) const override
Estimate for the molar volume of the liquid.
double daAlpha_dT() const
Calculate temperature derivative .
void calculatePressureDerivatives() const
Calculate and 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.
double m_aAlpha_mix
Value of in the equation of state.
double m_a
Value of in the equation of state.
static const double omega_a
Omega constant: a0 (= omega_a) used in Peng-Robinson equation of state.
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, b, and omega coefficients.
void calculateAB(double &aCalc, double &bCalc, double &aAlpha) const
Calculate the , , and parameters given the temperature.
void updateMixingExpressions() override
Update the , , and parameters.
double cp_mole() const override
Molar heat capacity at constant pressure and composition [J/kmol/K].
double m_b
Value of in the equation of state.
void getPartialMolarVolumes(span< double > vbar) const override
Return an array of partial molar volumes for the species in the mixture.
int solveCubic(double T, double pres, double a, double b, double aAlpha, span< double > Vroot) const
Prepare variables and call the function to solve the cubic 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.
AnyMap getAuxiliaryData() override
Return parameters used by the Peng-Robinson equation of state.
bool addSpecies(shared_ptr< Species > spec) override
Add a Species to this Phase.
double dpdVCalc(double T, double molarVol, double &presCalc) const override
Calculate the pressure and the pressure derivative given the temperature and the molar volume.
void getChemPotentials(span< double > mu) const override
Get the species chemical potentials. Units: J/kmol.
double d2aAlpha_dT2() const
Calculate second derivative .
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.
double m_dpdV
The derivative of the pressure with respect to the volume.
double speciesCritTemperature(double a, double b) const
Calculate species-specific critical temperature.
void getPartialMolarIntEnergies(span< double > ubar) const override
Return an array of partial molar internal energies for the species in the mixture.
vector< double > m_acentric
acentric factor for each species, length m_kk
vector< double > m_workS
Vector of size m_kk, used as a temporary holding area.
Definition Phase.h:899
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
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
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 RT() const
Return the Gas Constant multiplied by the current temperature.
AnyMap parameters(bool withInput=true) const
Returns the parameters of a ThermoPhase object such that an identical one could be reconstructed usin...
virtual AnyMap getAuxiliaryData()
Return intermediate or model-specific parameters used by particular derived classes.
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...
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
const double Sqrt2
Sqrt(2)
Definition ct_defs.h:74
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...