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
306void RedlichKwongMFTP::getPartialMolarVolumes(span<double> vbar) const
307{
308 checkArraySize("RedlichKwongMFTP::getPartialMolarVolumes", vbar.size(), m_kk);
309 for (size_t k = 0; k < m_kk; k++) {
310 m_pp[k] = 0.0;
311 for (size_t i = 0; i < m_kk; i++) {
312 size_t counter = k + m_kk*i;
313 m_pp[k] += moleFractions_[i] * a_vec_Curr_[counter];
314 }
315 }
316 for (size_t k = 0; k < m_kk; k++) {
317 m_workS[k] = 0.0;
318 for (size_t i = 0; i < m_kk; i++) {
319 size_t counter = k + m_kk*i;
320 m_workS[k] += moleFractions_[i] * a_coeff_vec(1,counter);
321 }
322 }
323
324 double sqt = sqrt(temperature());
325 double mv = molarVolume();
326 double vmb = mv - m_b_current;
327 double vpb = mv + m_b_current;
328 for (size_t k = 0; k < m_kk; k++) {
329 double num = (RT() + RT() * m_b_current/ vmb + RT() * b_vec_Curr_[k] / vmb
330 + RT() * m_b_current * b_vec_Curr_[k] /(vmb * vmb)
331 - 2.0 * m_pp[k] / (sqt * vpb)
332 + m_a_current * b_vec_Curr_[k] / (sqt * vpb * vpb)
333 );
334 double denom = (pressure() + RT() * m_b_current/(vmb * vmb) - m_a_current / (sqt * vpb * vpb)
335 );
336 vbar[k] = num / denom;
337 }
338}
339
340bool RedlichKwongMFTP::addSpecies(shared_ptr<Species> spec)
341{
342 bool added = MixtureFugacityTP::addSpecies(spec);
343 if (added) {
344 a_vec_Curr_.resize(m_kk * m_kk, 0.0);
345
346 // Initialize a_vec and b_vec to NaN, to screen for species with
347 // pureFluidParameters which are undefined in the input file:
348 b_vec_Curr_.push_back(NAN);
349 a_coeff_vec.resize(2, m_kk * m_kk, NAN);
350
351 m_pp.push_back(0.0);
352 m_coeffSource.push_back(CoeffSource::EoS);
353 m_partialMolarVolumes.push_back(0.0);
354 dpdni_.push_back(0.0);
355 }
356 return added;
357}
358
360{
361 // Contents of 'critical-properties.yaml', loaded later if needed
362 AnyMap critPropsDb;
363 std::unordered_map<string, AnyMap*> dbSpecies;
364
365 for (auto& [name, species] : m_species) {
366 auto& data = species->input;
367 size_t k = speciesIndex(name, true);
368 if (!isnan(a_coeff_vec(0, k + m_kk * k))) {
369 continue;
370 }
371 bool foundCoeffs = false;
372
373 if (data.hasKey("equation-of-state") &&
374 data["equation-of-state"].hasMapWhere("model", "Redlich-Kwong"))
375 {
376 // Read a and b coefficients from species 'input' information (that is, as
377 // specified in a YAML input file) specific to the Redlich-Kwong model
378 auto eos = data["equation-of-state"].getMapWhere(
379 "model", "Redlich-Kwong");
380
381 if (eos.hasKey("a") && eos.hasKey("b")) {
382 double a0 = 0, a1 = 0;
383 if (eos["a"].isScalar()) {
384 a0 = eos.convert("a", "Pa*m^6/kmol^2*K^0.5");
385 } else {
386 auto avec = eos["a"].asVector<AnyValue>(2);
387 a0 = eos.units().convert(avec[0], "Pa*m^6/kmol^2*K^0.5");
388 a1 = eos.units().convert(avec[1], "Pa*m^6/kmol^2/K^0.5");
389 }
390 double b = eos.convert("b", "m^3/kmol");
391 foundCoeffs = true;
392 setSpeciesCoeffs(name, a0, a1, b);
393 m_coeffSource[k] = CoeffSource::EoS;
394 }
395
396 if (eos.hasKey("binary-a")) {
397 AnyMap& binary_a = eos["binary-a"].as<AnyMap>();
398 const UnitSystem& units = binary_a.units();
399 for (auto& [name2, coeff] : binary_a) {
400 double a0 = 0, a1 = 0;
401 if (coeff.isScalar()) {
402 a0 = units.convert(coeff, "Pa*m^6/kmol^2*K^0.5");
403 } else {
404 auto avec = coeff.asVector<AnyValue>(2);
405 a0 = units.convert(avec[0], "Pa*m^6/kmol^2*K^0.5");
406 a1 = units.convert(avec[1], "Pa*m^6/kmol^2/K^0.5");
407 }
408 setBinaryCoeffs(name, name2, a0, a1);
409 }
410 }
411
412 if (foundCoeffs) {
413 continue;
414 }
415 }
416
417 // Coefficients have not been populated from model-specific input
418 double Tc = NAN, Pc = NAN;
419 if (data.hasKey("critical-parameters")) {
420 // Use critical state information stored in the species entry to
421 // calculate a and b
422 auto& critProps = data["critical-parameters"].as<AnyMap>();
423 Tc = critProps.convert("critical-temperature", "K");
424 Pc = critProps.convert("critical-pressure", "Pa");
425 m_coeffSource[k] = CoeffSource::CritProps;
426 } else {
427 // Search 'crit-properties.yaml' to find Tc and Pc. Load data if needed
428 if (critPropsDb.empty()) {
429 critPropsDb = AnyMap::fromYamlFile("critical-properties.yaml");
430 dbSpecies = critPropsDb["species"].asMap("name");
431 }
432
433 // All names in critical-properties.yaml are upper case
434 auto ucName = boost::algorithm::to_upper_copy(name);
435 if (dbSpecies.count(ucName)) {
436 auto& spec = *dbSpecies.at(ucName);
437 auto& critProps = spec["critical-parameters"].as<AnyMap>();
438 Tc = critProps.convert("critical-temperature", "K");
439 Pc = critProps.convert("critical-pressure", "Pa");
440 m_coeffSource[k] = CoeffSource::Database;
441 }
442 }
443
444 // Check if critical properties were found in either location
445 if (!isnan(Tc)) {
446 // Assuming no temperature dependence (that is, a1 = 0)
447 double a = omega_a * pow(GasConstant, 2) * pow(Tc, 2.5) / Pc;
448 double b = omega_b * GasConstant * Tc / Pc;
449 setSpeciesCoeffs(name, a, 0.0, b);
450 } else {
451 throw InputFileError("RedlichKwongMFTP::initThermo", data,
452 "No critical property or Redlich-Kwong parameters found "
453 "for species {}.", name);
454 }
455 }
456}
457
459 AnyMap& speciesNode) const
460{
462 size_t k = speciesIndex(name, true);
463 if (m_coeffSource[k] == CoeffSource::EoS) {
464 auto& eosNode = speciesNode["equation-of-state"].getMapWhere(
465 "model", "Redlich-Kwong", true);
466
467 size_t counter = k + m_kk * k;
468 if (a_coeff_vec(1, counter) != 0.0) {
469 vector<AnyValue> coeffs(2);
470 coeffs[0].setQuantity(a_coeff_vec(0, counter), "Pa*m^6/kmol^2*K^0.5");
471 coeffs[1].setQuantity(a_coeff_vec(1, counter), "Pa*m^6/kmol^2/K^0.5");
472 eosNode["a"] = std::move(coeffs);
473 } else {
474 eosNode["a"].setQuantity(a_coeff_vec(0, counter),
475 "Pa*m^6/kmol^2*K^0.5");
476 }
477 eosNode["b"].setQuantity(b_vec_Curr_[k], "m^3/kmol");
478 } else if (m_coeffSource[k] == CoeffSource::CritProps) {
479 auto& critProps = speciesNode["critical-parameters"];
480 double a = a_coeff_vec(0, k + m_kk * k);
481 double b = b_vec_Curr_[k];
482 double Tc = pow(a * omega_b / (b * omega_a * GasConstant), 2.0/3.0);
483 double Pc = omega_b * GasConstant * Tc / b;
484 critProps["critical-temperature"].setQuantity(Tc, "K");
485 critProps["critical-pressure"].setQuantity(Pc, "Pa");
486 }
487
488 if (m_binaryParameters.count(name)) {
489 auto& eosNode = speciesNode["equation-of-state"].getMapWhere(
490 "model", "Redlich-Kwong", true);
491 AnyMap bin_a;
492 for (const auto& [name2, coeffs] : m_binaryParameters.at(name)) {
493 if (coeffs.second == 0) {
494 bin_a[name2].setQuantity(coeffs.first, "Pa*m^6/kmol^2*K^0.5");
495 } else {
496 vector<AnyValue> C(2);
497 C[0].setQuantity(coeffs.first, "Pa*m^6/kmol^2*K^0.5");
498 C[1].setQuantity(coeffs.second, "Pa*m^6/kmol^2/K^0.5");
499 bin_a[name2] = std::move(C);
500 }
501 }
502 eosNode["binary-a"] = std::move(bin_a);
503 }
504}
505
507{
508 // note this agrees with tpx
509 double rho = density();
510 double mmw = meanMolecularWeight();
511 double molarV = mmw / rho;
512 double hh = m_b_current / molarV;
513 double zz = z();
514 double dadt = da_dt();
515 double T = temperature();
516 double sqT = sqrt(T);
517 double fac = dadt - m_a_current / (2.0 * T);
518 double sresid_mol_R = log(zz*(1.0 - hh)) + log(1.0 + hh) * fac / (sqT * GasConstant * m_b_current);
519 return GasConstant * sresid_mol_R;
520}
521
523{
524 // note this agrees with tpx
525 double rho = density();
526 double mmw = meanMolecularWeight();
527 double molarV = mmw / rho;
528 double hh = m_b_current / molarV;
529 double zz = z();
530 double dadt = da_dt();
531 double T = temperature();
532 double sqT = sqrt(T);
533 double fac = T * dadt - 3.0 *m_a_current / (2.0);
534 return GasConstant * T * (zz - 1.0) + fac * log(1.0 + hh) / (sqT * m_b_current);
535}
536
537double RedlichKwongMFTP::liquidVolEst(double TKelvin, double& presGuess) const
538{
539 double v = m_b_current * 1.1;
540 double atmp;
541 double btmp;
542 calculateAB(TKelvin, atmp, btmp);
543 double pres = std::max(psatEst(TKelvin), presGuess);
544 double Vroot[3];
545 bool foundLiq = false;
546 int m = 0;
547 while (m < 100 && !foundLiq) {
548 int nsol = solveCubic(TKelvin, pres, atmp, btmp, Vroot);
549 if (nsol == 1 || nsol == 2) {
550 double pc = critPressure();
551 if (pres > pc) {
552 foundLiq = true;
553 }
554 pres *= 1.04;
555 } else {
556 foundLiq = true;
557 }
558 }
559
560 if (foundLiq) {
561 v = Vroot[0];
562 presGuess = pres;
563 } else {
564 v = -1.0;
565 }
566 return v;
567}
568
569double RedlichKwongMFTP::densityCalc(double TKelvin, double presPa, int phaseRequested,
570 double rhoguess)
571{
572 // It's necessary to set the temperature so that m_a_current is set correctly.
573 setTemperature(TKelvin);
574 double tcrit = critTemperature();
575 double mmw = meanMolecularWeight();
576 if (rhoguess == -1.0) {
577 if (phaseRequested >= FLUID_LIQUID_0) {
578 double lqvol = liquidVolEst(TKelvin, presPa);
579 rhoguess = mmw / lqvol;
580 }
581 } else {
582 // Assume the Gas phase initial guess, if nothing is specified to the routine
583 rhoguess = presPa * mmw / (GasConstant * TKelvin);
584 }
585
586
587 double volguess = mmw / rhoguess;
588 NSolns_ = solveCubic(TKelvin, presPa, m_a_current, m_b_current, Vroot_);
589
590 double molarVolLast = Vroot_[0];
591 if (NSolns_ >= 2) {
592 if (phaseRequested >= FLUID_LIQUID_0) {
593 molarVolLast = Vroot_[0];
594 } else if (phaseRequested == FLUID_GAS || phaseRequested == FLUID_SUPERCRIT) {
595 molarVolLast = Vroot_[2];
596 } else {
597 if (volguess > Vroot_[1]) {
598 molarVolLast = Vroot_[2];
599 } else {
600 molarVolLast = Vroot_[0];
601 }
602 }
603 } else if (NSolns_ == 1) {
604 if (phaseRequested == FLUID_GAS || phaseRequested == FLUID_SUPERCRIT || phaseRequested == FLUID_UNDEFINED) {
605 molarVolLast = Vroot_[0];
606 } else {
607 return -2.0;
608 }
609 } else if (NSolns_ == -1) {
610 if (phaseRequested >= FLUID_LIQUID_0 || phaseRequested == FLUID_UNDEFINED || phaseRequested == FLUID_SUPERCRIT) {
611 molarVolLast = Vroot_[0];
612 } else if (TKelvin > tcrit) {
613 molarVolLast = Vroot_[0];
614 } else {
615 return -2.0;
616 }
617 } else {
618 molarVolLast = Vroot_[0];
619 return -1.0;
620 }
621 return mmw / molarVolLast;
622}
623
624double RedlichKwongMFTP::dpdVCalc(double TKelvin, double molarVol, double& presCalc) const
625{
626 double sqt = sqrt(TKelvin);
627 presCalc = GasConstant * TKelvin / (molarVol - m_b_current)
628 - m_a_current / (sqt * molarVol * (molarVol + m_b_current));
629
630 double vpb = molarVol + m_b_current;
631 double vmb = molarVol - m_b_current;
632 double dpdv = (- GasConstant * TKelvin / (vmb * vmb)
633 + m_a_current * (2 * molarVol + m_b_current) / (sqt * molarVol * molarVol * vpb * vpb));
634 return dpdv;
635}
636
638{
640 return -1 / (molarVolume() * dpdV_);
641}
642
644{
646 return -dpdT_ / (molarVolume() * dpdV_);
647}
648
650{
652 return molarVolume() * sqrt(-cp_mole() / cv_mole() * dpdV_ / meanMolecularWeight());
653}
654
656{
657 double TKelvin = temperature();
658 double mv = molarVolume();
659 double pres;
660
661 dpdV_ = dpdVCalc(TKelvin, mv, pres);
662 double sqt = sqrt(TKelvin);
663 double vpb = mv + m_b_current;
664 double vmb = mv - m_b_current;
665 double dadt = da_dt();
666 double fac = dadt - m_a_current/(2.0 * TKelvin);
667 dpdT_ = (GasConstant / vmb - fac / (sqt * mv * vpb));
668}
669
671{
672 double temp = temperature();
673 if (m_formTempParam == 1) {
674 for (size_t i = 0; i < m_kk; i++) {
675 for (size_t j = 0; j < m_kk; j++) {
676 size_t counter = i * m_kk + j;
677 a_vec_Curr_[counter] = a_coeff_vec(0,counter) + a_coeff_vec(1,counter) * temp;
678 }
679 }
680 }
681
682 m_b_current = 0.0;
683 m_a_current = 0.0;
684 for (size_t i = 0; i < m_kk; i++) {
685 m_b_current += moleFractions_[i] * b_vec_Curr_[i];
686 for (size_t j = 0; j < m_kk; j++) {
687 m_a_current += a_vec_Curr_[i * m_kk + j] * moleFractions_[i] * moleFractions_[j];
688 }
689 }
690 if (isnan(m_b_current)) {
691 // One or more species do not have specified coefficients.
692 fmt::memory_buffer b;
693 for (size_t k = 0; k < m_kk; k++) {
694 if (isnan(b_vec_Curr_[k])) {
695 if (b.size() > 0) {
696 fmt_append(b, ", {}", speciesName(k));
697 } else {
698 fmt_append(b, "{}", speciesName(k));
699 }
700 }
701 }
702 throw CanteraError("RedlichKwongMFTP::updateMixingExpressions",
703 "Missing Redlich-Kwong coefficients for species: {}", to_string(b));
704 }
705}
706
707void RedlichKwongMFTP::calculateAB(double temp, double& aCalc, double& bCalc) const
708{
709 bCalc = 0.0;
710 aCalc = 0.0;
711 if (m_formTempParam == 1) {
712 for (size_t i = 0; i < m_kk; i++) {
713 bCalc += moleFractions_[i] * b_vec_Curr_[i];
714 for (size_t j = 0; j < m_kk; j++) {
715 size_t counter = i * m_kk + j;
716 double a_vec_Curr = a_coeff_vec(0,counter) + a_coeff_vec(1,counter) * temp;
717 aCalc += a_vec_Curr * moleFractions_[i] * moleFractions_[j];
718 }
719 }
720 } else {
721 for (size_t i = 0; i < m_kk; i++) {
722 bCalc += moleFractions_[i] * b_vec_Curr_[i];
723 for (size_t j = 0; j < m_kk; j++) {
724 size_t counter = i * m_kk + j;
725 double a_vec_Curr = a_coeff_vec(0,counter);
726 aCalc += a_vec_Curr * moleFractions_[i] * moleFractions_[j];
727 }
728 }
729 }
730}
731
732double RedlichKwongMFTP::da_dt() const
733{
734 double dadT = 0.0;
735 if (m_formTempParam == 1) {
736 for (size_t i = 0; i < m_kk; i++) {
737 for (size_t j = 0; j < m_kk; j++) {
738 size_t counter = i * m_kk + j;
739 dadT+= a_coeff_vec(1,counter) * moleFractions_[i] * moleFractions_[j];
740 }
741 }
742 }
743 return dadT;
744}
745
746void RedlichKwongMFTP::calcCriticalConditions(double& pc, double& tc, double& vc) const
747{
748 double a0 = 0.0;
749 double aT = 0.0;
750 for (size_t i = 0; i < m_kk; i++) {
751 for (size_t j = 0; j <m_kk; j++) {
752 size_t counter = i + m_kk * j;
753 a0 += moleFractions_[i] * moleFractions_[j] * a_coeff_vec(0, counter);
754 aT += moleFractions_[i] * moleFractions_[j] * a_coeff_vec(1, counter);
755 }
756 }
757 double a = m_a_current;
758 double b = m_b_current;
759 if (m_formTempParam != 0) {
760 a = a0;
761 }
762 if (b <= 0.0) {
763 tc = 1000000.;
764 pc = 1.0E13;
765 vc = omega_vc * GasConstant * tc / pc;
766 return;
767 }
768 if (a <= 0.0) {
769 tc = 0.0;
770 pc = 0.0;
771 vc = 2.0 * b;
772 return;
773 }
774 double tmp = a * omega_b / (b * omega_a * GasConstant);
775 double pp = 2./3.;
776 double sqrttc, f, dfdt, deltatc;
777
778 if (m_formTempParam == 0) {
779 tc = pow(tmp, pp);
780 } else {
781 tc = pow(tmp, pp);
782 for (int j = 0; j < 10; j++) {
783 sqrttc = sqrt(tc);
784 f = omega_a * b * GasConstant * tc * sqrttc / omega_b - aT * tc - a0;
785 dfdt = 1.5 * omega_a * b * GasConstant * sqrttc / omega_b - aT;
786 deltatc = - f / dfdt;
787 tc += deltatc;
788 }
789 if (deltatc > 0.1) {
790 throw CanteraError("RedlichKwongMFTP::calcCriticalConditions",
791 "didn't converge");
792 }
793 }
794
795 pc = omega_b * GasConstant * tc / b;
796 vc = omega_vc * GasConstant * tc / pc;
797}
798
799int RedlichKwongMFTP::solveCubic(double T, double pres, double a, double b,
800 span<double> Vroot) const
801{
802
803 // Derive the coefficients of the cubic polynomial to solve.
804 double an = 1.0;
805 double bn = - GasConstant * T / pres;
806 double sqt = sqrt(T);
807 double cn = - (GasConstant * T * b / pres - a/(pres * sqt) + b * b);
808 double dn = - (a * b / (pres * sqt));
809
810 double tmp = a * omega_b / (b * omega_a * GasConstant);
811 double pp = 2./3.;
812 double tc = pow(tmp, pp);
813 double pc = omega_b * GasConstant * tc / b;
814 double vc = omega_vc * GasConstant * tc / pc;
815
816 int nSolnValues = MixtureFugacityTP::solveCubic(T, pres, a, b, a, Vroot, an, bn, cn, dn, tc, vc);
817
818 return nSolnValues;
819}
820
821}
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).
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.
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).
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.
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.
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.
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...