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