Cantera  2.5.1
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 
9 #include "cantera/base/ctml.h"
10 
11 #include <boost/math/tools/roots.hpp>
12 
13 #include <algorithm>
14 
15 using namespace std;
16 namespace bmt = boost::math::tools;
17 
18 namespace Cantera
19 {
20 
21 const doublereal RedlichKwongMFTP::omega_a = 4.27480233540E-01;
22 const doublereal RedlichKwongMFTP::omega_b = 8.66403499650E-02;
23 const doublereal RedlichKwongMFTP::omega_vc = 3.33333333333333E-01;
24 
25 RedlichKwongMFTP::RedlichKwongMFTP() :
26  m_formTempParam(0),
27  m_b_current(0.0),
28  m_a_current(0.0),
29  NSolns_(0),
30  dpdV_(0.0),
31  dpdT_(0.0)
32 {
33  fill_n(Vroot_, 3, 0.0);
34 }
35 
36 RedlichKwongMFTP::RedlichKwongMFTP(const std::string& infile, const std::string& id_) :
37  m_formTempParam(0),
38  m_b_current(0.0),
39  m_a_current(0.0),
40  NSolns_(0),
41  dpdV_(0.0),
42  dpdT_(0.0)
43 {
44  fill_n(Vroot_, 3, 0.0);
45  initThermoFile(infile, id_);
46 }
47 
48 RedlichKwongMFTP::RedlichKwongMFTP(XML_Node& phaseRefRoot, const std::string& id_) :
49  m_formTempParam(0),
50  m_b_current(0.0),
51  m_a_current(0.0),
52  NSolns_(0),
53  dpdV_(0.0),
54  dpdT_(0.0)
55 {
56  fill_n(Vroot_, 3, 0.0);
57  importPhase(phaseRefRoot, this);
58 }
59 
60 void RedlichKwongMFTP::setSpeciesCoeffs(const std::string& species,
61  double a0, double a1, double b)
62 {
63  size_t k = speciesIndex(species);
64  if (k == npos) {
65  throw CanteraError("RedlichKwongMFTP::setSpeciesCoeffs",
66  "Unknown species '{}'.", species);
67  }
68 
69  if (a1 != 0.0) {
70  m_formTempParam = 1; // expression is temperature-dependent
71  }
72 
73  size_t counter = k + m_kk * k;
74  a_coeff_vec(0, counter) = a0;
75  a_coeff_vec(1, counter) = a1;
76 
77  // standard mixing rule for cross-species interaction term
78  for (size_t j = 0; j < m_kk; j++) {
79  if (k == j) {
80  continue;
81  }
82 
83  // a_coeff_vec(0) is initialized to NaN to mark uninitialized species
84  if (isnan(a_coeff_vec(0, j + m_kk * j))) {
85  // The diagonal element of the jth species has not yet been defined.
86  continue;
87  } else if (isnan(a_coeff_vec(0, j + m_kk * k))) {
88  // Only use the mixing rules if the off-diagonal element has not already been defined by a
89  // user-specified crossFluidParameters entry:
90  double a0kj = sqrt(a_coeff_vec(0, j + m_kk * j) * a0);
91  double a1kj = sqrt(a_coeff_vec(1, j + m_kk * j) * a1);
92  a_coeff_vec(0, j + m_kk * k) = a0kj;
93  a_coeff_vec(1, j + m_kk * k) = a1kj;
94  a_coeff_vec(0, k + m_kk * j) = a0kj;
95  a_coeff_vec(1, k + m_kk * j) = a1kj;
96  }
97  }
98  a_coeff_vec.getRow(0, a_vec_Curr_.data());
99  b_vec_Curr_[k] = b;
100 }
101 
102 void RedlichKwongMFTP::setBinaryCoeffs(const std::string& species_i,
103  const std::string& species_j, double a0, double a1)
104 {
105  size_t ki = speciesIndex(species_i);
106  if (ki == npos) {
107  throw CanteraError("RedlichKwongMFTP::setBinaryCoeffs",
108  "Unknown species '{}'.", species_i);
109  }
110  size_t kj = speciesIndex(species_j);
111  if (kj == npos) {
112  throw CanteraError("RedlichKwongMFTP::setBinaryCoeffs",
113  "Unknown species '{}'.", species_j);
114  }
115 
116  if (a1 != 0.0) {
117  m_formTempParam = 1; // expression is temperature-dependent
118  }
119  size_t counter1 = ki + m_kk * kj;
120  size_t counter2 = kj + m_kk * ki;
121  a_coeff_vec(0, counter1) = a_coeff_vec(0, counter2) = a0;
122  a_coeff_vec(1, counter1) = a_coeff_vec(1, counter2) = a1;
123  a_vec_Curr_[counter1] = a_vec_Curr_[counter2] = a0;
124 }
125 
126 // ------------Molar Thermodynamic Properties -------------------------
127 
129 {
131  doublereal h_ideal = RT() * mean_X(m_h0_RT);
132  doublereal h_nonideal = hresid();
133  return h_ideal + h_nonideal;
134 }
135 
137 {
139  doublereal sr_ideal = GasConstant * (mean_X(m_s0_R)
140  - sum_xlogx() - std::log(pressure()/refPressure()));
141  doublereal sr_nonideal = sresid();
142  return sr_ideal + sr_nonideal;
143 }
144 
145 doublereal RedlichKwongMFTP::cp_mole() const
146 {
148  doublereal TKelvin = temperature();
149  doublereal sqt = sqrt(TKelvin);
150  doublereal mv = molarVolume();
151  doublereal vpb = mv + m_b_current;
153  doublereal cpref = GasConstant * mean_X(m_cp0_R);
154  doublereal dadt = da_dt();
155  doublereal fac = TKelvin * dadt - 3.0 * m_a_current / 2.0;
156  doublereal dHdT_V = (cpref + mv * dpdT_ - GasConstant - 1.0 / (2.0 * m_b_current * TKelvin * sqt) * log(vpb/mv) * fac
157  +1.0/(m_b_current * sqt) * log(vpb/mv) * (-0.5 * dadt));
158  return dHdT_V - (mv + TKelvin * dpdT_ / dpdV_) * dpdT_;
159 }
160 
161 doublereal RedlichKwongMFTP::cv_mole() const
162 {
164  doublereal TKelvin = temperature();
165  doublereal sqt = sqrt(TKelvin);
166  doublereal mv = molarVolume();
167  doublereal vpb = mv + m_b_current;
168  doublereal cvref = GasConstant * (mean_X(m_cp0_R) - 1.0);
169  doublereal dadt = da_dt();
170  doublereal fac = TKelvin * dadt - 3.0 * m_a_current / 2.0;
171  return (cvref - 1.0/(2.0 * m_b_current * TKelvin * sqt) * log(vpb/mv)*fac
172  +1.0/(m_b_current * sqt) * log(vpb/mv)*(-0.5*dadt));
173 }
174 
175 doublereal RedlichKwongMFTP::pressure() const
176 {
178 
179  // Get a copy of the private variables stored in the State object
180  doublereal T = temperature();
181  double molarV = meanMolecularWeight() / density();
182  double pp = GasConstant * T/(molarV - m_b_current) - m_a_current/(sqrt(T) * molarV * (molarV + m_b_current));
183  return pp;
184 }
185 
187 {
188  // Calculate the molarVolume of the solution (m**3 kmol-1)
189  const doublereal* const dtmp = moleFractdivMMW();
191  double invDens = dot(m_tmpV.begin(), m_tmpV.end(), dtmp);
192 
193  // Set the density in the parent State object directly, by calling the
194  // Phase::setDensity() function.
195  Phase::setDensity(1.0/invDens);
196 }
197 
198 void RedlichKwongMFTP::setTemperature(const doublereal temp)
199 {
200  Phase::setTemperature(temp);
202  updateAB();
203 }
204 
206 {
208  updateAB();
209 }
210 
212 {
214  for (size_t k = 0; k < m_kk; k++) {
215  c[k] *= moleFraction(k)*pressure()/RT();
216  }
217 }
218 
219 doublereal RedlichKwongMFTP::standardConcentration(size_t k) const
220 {
221  getStandardVolumes(m_tmpV.data());
222  return 1.0 / m_tmpV[k];
223 }
224 
226 {
227  doublereal mv = molarVolume();
228  doublereal sqt = sqrt(temperature());
229  doublereal vpb = mv + m_b_current;
230  doublereal vmb = mv - m_b_current;
231 
232  for (size_t k = 0; k < m_kk; k++) {
233  m_pp[k] = 0.0;
234  for (size_t i = 0; i < m_kk; i++) {
235  size_t counter = k + m_kk*i;
236  m_pp[k] += moleFractions_[i] * a_vec_Curr_[counter];
237  }
238  }
239  doublereal pres = pressure();
240 
241  for (size_t k = 0; k < m_kk; k++) {
242  ac[k] = (- RT() * log(pres * mv / RT())
243  + RT() * log(mv / vmb)
244  + RT() * b_vec_Curr_[k] / vmb
245  - 2.0 * m_pp[k] / (m_b_current * sqt) * log(vpb/mv)
246  + m_a_current * b_vec_Curr_[k] / (m_b_current * m_b_current * sqt) * log(vpb/mv)
247  - m_a_current / (m_b_current * sqt) * (b_vec_Curr_[k]/vpb)
248  );
249  }
250  for (size_t k = 0; k < m_kk; k++) {
251  ac[k] = exp(ac[k]/RT());
252  }
253 }
254 
255 // ---- Partial Molar Properties of the Solution -----------------
256 
257 void RedlichKwongMFTP::getChemPotentials_RT(doublereal* muRT) const
258 {
259  getChemPotentials(muRT);
260  for (size_t k = 0; k < m_kk; k++) {
261  muRT[k] *= 1.0 / RT();
262  }
263 }
264 
265 void RedlichKwongMFTP::getChemPotentials(doublereal* mu) const
266 {
267  getGibbs_ref(mu);
268  for (size_t k = 0; k < m_kk; k++) {
269  double xx = std::max(SmallNumber, moleFraction(k));
270  mu[k] += RT()*(log(xx));
271  }
272 
273  doublereal mv = molarVolume();
274  doublereal sqt = sqrt(temperature());
275  doublereal vpb = mv + m_b_current;
276  doublereal vmb = mv - m_b_current;
277 
278  for (size_t k = 0; k < m_kk; k++) {
279  m_pp[k] = 0.0;
280  for (size_t i = 0; i < m_kk; i++) {
281  size_t counter = k + m_kk*i;
282  m_pp[k] += moleFractions_[i] * a_vec_Curr_[counter];
283  }
284  }
285  doublereal pres = pressure();
286  doublereal refP = refPressure();
287 
288  for (size_t k = 0; k < m_kk; k++) {
289  mu[k] += (RT() * log(pres/refP) - RT() * log(pres * mv / RT())
290  + RT() * log(mv / vmb)
291  + RT() * b_vec_Curr_[k] / vmb
292  - 2.0 * m_pp[k] / (m_b_current * sqt) * log(vpb/mv)
293  + m_a_current * b_vec_Curr_[k] / (m_b_current * m_b_current * sqt) * log(vpb/mv)
294  - m_a_current / (m_b_current * sqt) * (b_vec_Curr_[k]/vpb)
295  );
296  }
297 }
298 
300 {
301  // First we get the reference state contributions
302  getEnthalpy_RT_ref(hbar);
303  scale(hbar, hbar+m_kk, hbar, RT());
304 
305  // We calculate dpdni_
306  doublereal TKelvin = temperature();
307  doublereal mv = molarVolume();
308  doublereal sqt = sqrt(TKelvin);
309  doublereal vpb = mv + m_b_current;
310  doublereal vmb = mv - m_b_current;
311  for (size_t k = 0; k < m_kk; k++) {
312  m_pp[k] = 0.0;
313  for (size_t i = 0; i < m_kk; i++) {
314  size_t counter = k + m_kk*i;
315  m_pp[k] += moleFractions_[i] * a_vec_Curr_[counter];
316  }
317  }
318  for (size_t k = 0; k < m_kk; k++) {
319  dpdni_[k] = RT()/vmb + RT() * b_vec_Curr_[k] / (vmb * vmb) - 2.0 * m_pp[k] / (sqt * mv * vpb)
320  + m_a_current * b_vec_Curr_[k]/(sqt * mv * vpb * vpb);
321  }
322  doublereal dadt = da_dt();
323  doublereal fac = TKelvin * dadt - 3.0 * m_a_current / 2.0;
324 
325  for (size_t k = 0; k < m_kk; k++) {
326  m_tmpV[k] = 0.0;
327  for (size_t i = 0; i < m_kk; i++) {
328  size_t counter = k + m_kk*i;
329  m_tmpV[k] += 2.0 * moleFractions_[i] * TKelvin * a_coeff_vec(1,counter) - 3.0 * moleFractions_[i] * a_vec_Curr_[counter];
330  }
331  }
332 
334  doublereal fac2 = mv + TKelvin * dpdT_ / dpdV_;
335  for (size_t k = 0; k < m_kk; k++) {
336  double hE_v = (mv * dpdni_[k] - RT() - b_vec_Curr_[k]/ (m_b_current * m_b_current * sqt) * log(vpb/mv)*fac
337  + 1.0 / (m_b_current * sqt) * log(vpb/mv) * m_tmpV[k]
338  + b_vec_Curr_[k] / vpb / (m_b_current * sqt) * fac);
339  hbar[k] = hbar[k] + hE_v;
340  hbar[k] -= fac2 * dpdni_[k];
341  }
342 }
343 
344 void RedlichKwongMFTP::getPartialMolarEntropies(doublereal* sbar) const
345 {
346  getEntropy_R_ref(sbar);
347  scale(sbar, sbar+m_kk, sbar, GasConstant);
348  doublereal TKelvin = temperature();
349  doublereal sqt = sqrt(TKelvin);
350  doublereal mv = molarVolume();
351  doublereal refP = refPressure();
352 
353  for (size_t k = 0; k < m_kk; k++) {
354  doublereal xx = std::max(SmallNumber, moleFraction(k));
355  sbar[k] += GasConstant * (- log(xx));
356  }
357  for (size_t k = 0; k < m_kk; k++) {
358  m_pp[k] = 0.0;
359  for (size_t i = 0; i < m_kk; i++) {
360  size_t counter = k + m_kk*i;
361  m_pp[k] += moleFractions_[i] * a_vec_Curr_[counter];
362  }
363  }
364  for (size_t k = 0; k < m_kk; k++) {
365  m_tmpV[k] = 0.0;
366  for (size_t i = 0; i < m_kk; i++) {
367  size_t counter = k + m_kk*i;
368  m_tmpV[k] += moleFractions_[i] * a_coeff_vec(1,counter);
369  }
370  }
371 
372  doublereal dadt = da_dt();
373  doublereal fac = dadt - m_a_current / (2.0 * TKelvin);
374  doublereal vmb = mv - m_b_current;
375  doublereal vpb = mv + m_b_current;
376  for (size_t k = 0; k < m_kk; k++) {
377  sbar[k] -=(GasConstant * log(GasConstant * TKelvin / (refP * mv))
378  + GasConstant
379  + GasConstant * log(mv/vmb)
380  + GasConstant * b_vec_Curr_[k]/vmb
381  + m_pp[k]/(m_b_current * TKelvin * sqt) * log(vpb/mv)
382  - 2.0 * m_tmpV[k]/(m_b_current * sqt) * log(vpb/mv)
383  + b_vec_Curr_[k] / (m_b_current * m_b_current * sqt) * log(vpb/mv) * fac
384  - 1.0 / (m_b_current * sqt) * b_vec_Curr_[k] / vpb * fac
385  );
386  }
387 
389  getPartialMolarVolumes(m_partialMolarVolumes.data());
390  for (size_t k = 0; k < m_kk; k++) {
391  sbar[k] -= -m_partialMolarVolumes[k] * dpdT_;
392  }
393 }
394 
396 {
397  getIntEnergy_RT(ubar);
398  scale(ubar, ubar+m_kk, ubar, RT());
399 }
400 
401 void RedlichKwongMFTP::getPartialMolarCp(doublereal* cpbar) const
402 {
403  getCp_R(cpbar);
404  scale(cpbar, cpbar+m_kk, cpbar, GasConstant);
405 }
406 
407 void RedlichKwongMFTP::getPartialMolarVolumes(doublereal* vbar) const
408 {
409  for (size_t k = 0; k < m_kk; k++) {
410  m_pp[k] = 0.0;
411  for (size_t i = 0; i < m_kk; i++) {
412  size_t counter = k + m_kk*i;
413  m_pp[k] += moleFractions_[i] * a_vec_Curr_[counter];
414  }
415  }
416  for (size_t k = 0; k < m_kk; k++) {
417  m_tmpV[k] = 0.0;
418  for (size_t i = 0; i < m_kk; i++) {
419  size_t counter = k + m_kk*i;
420  m_tmpV[k] += moleFractions_[i] * a_coeff_vec(1,counter);
421  }
422  }
423 
424  doublereal sqt = sqrt(temperature());
425  doublereal mv = molarVolume();
426  doublereal vmb = mv - m_b_current;
427  doublereal vpb = mv + m_b_current;
428  for (size_t k = 0; k < m_kk; k++) {
429  doublereal num = (RT() + RT() * m_b_current/ vmb + RT() * b_vec_Curr_[k] / vmb
430  + RT() * m_b_current * b_vec_Curr_[k] /(vmb * vmb)
431  - 2.0 * m_pp[k] / (sqt * vpb)
432  + m_a_current * b_vec_Curr_[k] / (sqt * vpb * vpb)
433  );
434  doublereal denom = (pressure() + RT() * m_b_current/(vmb * vmb) - m_a_current / (sqt * vpb * vpb)
435  );
436  vbar[k] = num / denom;
437  }
438 }
439 
441 {
442  double pc, tc, vc;
443  double a0 = 0.0;
444  double aT = 0.0;
445  for (size_t i = 0; i < m_kk; i++) {
446  for (size_t j = 0; j <m_kk; j++) {
447  size_t counter = i + m_kk * j;
448  a0 += moleFractions_[i] * moleFractions_[j] * a_coeff_vec(0, counter);
449  aT += moleFractions_[i] * moleFractions_[j] * a_coeff_vec(1, counter);
450  }
451  }
452  calcCriticalConditions(m_a_current, m_b_current, a0, aT, pc, tc, vc);
453  return tc;
454 }
455 
457 {
458  double pc, tc, vc;
459  double a0 = 0.0;
460  double aT = 0.0;
461  for (size_t i = 0; i < m_kk; i++) {
462  for (size_t j = 0; j <m_kk; j++) {
463  size_t counter = i + m_kk * j;
464  a0 += moleFractions_[i] * moleFractions_[j] * a_coeff_vec(0, counter);
465  aT += moleFractions_[i] * moleFractions_[j] * a_coeff_vec(1, counter);
466  }
467  }
468  calcCriticalConditions(m_a_current, m_b_current, a0, aT, pc, tc, vc);
469  return pc;
470 }
471 
473 {
474  double pc, tc, vc;
475  double a0 = 0.0;
476  double aT = 0.0;
477  for (size_t i = 0; i < m_kk; i++) {
478  for (size_t j = 0; j <m_kk; j++) {
479  size_t counter = i + m_kk * j;
480  a0 += moleFractions_[i] * moleFractions_[j] * a_coeff_vec(0, counter);
481  aT += moleFractions_[i] * moleFractions_[j] * a_coeff_vec(1, counter);
482  }
483  }
484  calcCriticalConditions(m_a_current, m_b_current, a0, aT, pc, tc, vc);
485  return vc;
486 }
487 
489 {
490  double pc, tc, vc;
491  double a0 = 0.0;
492  double aT = 0.0;
493  for (size_t i = 0; i < m_kk; i++) {
494  for (size_t j = 0; j <m_kk; j++) {
495  size_t counter = i + m_kk * j;
496  a0 += moleFractions_[i] * moleFractions_[j] * a_coeff_vec(0, counter);
497  aT += moleFractions_[i] * moleFractions_[j] * a_coeff_vec(1, counter);
498  }
499  }
500  calcCriticalConditions(m_a_current, m_b_current, a0, aT, pc, tc, vc);
501  return pc*vc/tc/GasConstant;
502 }
503 
505 {
506  double pc, tc, vc;
507  double a0 = 0.0;
508  double aT = 0.0;
509  for (size_t i = 0; i < m_kk; i++) {
510  for (size_t j = 0; j <m_kk; j++) {
511  size_t counter = i + m_kk * j;
512  a0 += moleFractions_[i] * moleFractions_[j] *a_coeff_vec(0, counter);
513  aT += moleFractions_[i] * moleFractions_[j] *a_coeff_vec(1, counter);
514  }
515  }
516  calcCriticalConditions(m_a_current, m_b_current, a0, aT, pc, tc, vc);
517  double mmw = meanMolecularWeight();
518  return mmw / vc;
519 }
520 
521 bool RedlichKwongMFTP::addSpecies(shared_ptr<Species> spec)
522 {
523  bool added = MixtureFugacityTP::addSpecies(spec);
524  if (added) {
525  a_vec_Curr_.resize(m_kk * m_kk, 0.0);
526 
527  // Initialize a_vec and b_vec to NaN, to screen for species with
528  // pureFluidParameters which are undefined in the input file:
529  b_vec_Curr_.push_back(NAN);
530  a_coeff_vec.resize(2, m_kk * m_kk, NAN);
531 
532  m_pp.push_back(0.0);
533  m_tmpV.push_back(0.0);
534  m_partialMolarVolumes.push_back(0.0);
535  dpdni_.push_back(0.0);
536  }
537  return added;
538 }
539 
540 void RedlichKwongMFTP::initThermoXML(XML_Node& phaseNode, const std::string& id)
541 {
542  if (phaseNode.hasChild("thermo")) {
543  XML_Node& thermoNode = phaseNode.child("thermo");
544  std::string model = thermoNode["model"];
545  if (model != "RedlichKwong" && model != "RedlichKwongMFTP") {
546  throw CanteraError("RedlichKwongMFTP::initThermoXML",
547  "Unknown thermo model : " + model);
548  }
549 
550  // Reset any coefficients which may have been set using values from
551  // 'critProperties.xml' as part of non-XML initialization, so that
552  // off-diagonal elements can be correctly initialized
553  a_coeff_vec.data().assign(a_coeff_vec.data().size(), NAN);
554 
555  // Go get all of the coefficients and factors in the
556  // activityCoefficients XML block
557  if (thermoNode.hasChild("activityCoefficients")) {
558  XML_Node& acNode = thermoNode.child("activityCoefficients");
559 
560  // Loop through the children and read out fluid parameters. Process
561  // all the pureFluidParameters, first:
562  // Loop back through the "activityCoefficients" children and process the
563  // crossFluidParameters in the XML tree:
564  for (size_t i = 0; i < acNode.nChildren(); i++) {
565  XML_Node& xmlACChild = acNode.child(i);
566  if (caseInsensitiveEquals(xmlACChild.name(), "purefluidparameters")) {
567  readXMLPureFluid(xmlACChild);
568  } else if (caseInsensitiveEquals(xmlACChild.name(), "crossfluidparameters")) {
569  readXMLCrossFluid(xmlACChild);
570  }
571  }
572  }
573  // If any species exist which have undefined pureFluidParameters,
574  // search the database in 'critProperties.xml' to find critical
575  // temperature and pressure to calculate a and b.
576 
577  // Loop through all species in the CTI file
578  size_t iSpecies = 0;
579 
580  for (size_t i = 0; i < m_kk; i++) {
581  string iName = speciesName(i);
582 
583  // Get the index of the species
584  iSpecies = speciesIndex(iName);
585 
586  // Check if a and b are already populated (only the diagonal elements of a).
587  size_t counter = iSpecies + m_kk * iSpecies;
588 
589  // If not, then search the database:
590  if (isnan(a_coeff_vec(0, counter))) {
591 
592  vector<double> coeffArray;
593 
594  // Search the database for the species name and calculate
595  // coefficients a and b, from critical properties:
596  // coeffArray[0] = a0, coeffArray[1] = b;
597  coeffArray = getCoeff(iName);
598 
599  // Check if species was found in the database of critical properties,
600  // and assign the results
601  if (!isnan(coeffArray[0])) {
602  //Assuming no temperature dependence (i,e a1 = 0)
603  setSpeciesCoeffs(iName, coeffArray[0], 0.0, coeffArray[1]);
604  }
605  }
606  }
607  }
608 
609  MixtureFugacityTP::initThermoXML(phaseNode, id);
610 }
611 
613 {
614  for (auto& item : m_species) {
615  // Read a and b coefficients from species 'input' information (i.e. as
616  // specified in a YAML input file)
617  if (item.second->input.hasKey("equation-of-state")) {
618  auto eos = item.second->input["equation-of-state"].getMapWhere(
619  "model", "Redlich-Kwong");
620  double a0 = 0, a1 = 0;
621  if (eos["a"].isScalar()) {
622  a0 = eos.convert("a", "Pa*m^6/kmol^2*K^0.5");
623  } else {
624  auto avec = eos["a"].asVector<AnyValue>(2);
625  a0 = eos.units().convert(avec[0], "Pa*m^6/kmol^2*K^0.5");
626  a1 = eos.units().convert(avec[1], "Pa*m^6/kmol^2/K^0.5");
627  }
628  double b = eos.convert("b", "m^3/kmol");
629  setSpeciesCoeffs(item.first, a0, a1, b);
630  if (eos.hasKey("binary-a")) {
631  AnyMap& binary_a = eos["binary-a"].as<AnyMap>();
632  const UnitSystem& units = binary_a.units();
633  for (auto& item2 : binary_a) {
634  double a0 = 0, a1 = 0;
635  if (item2.second.isScalar()) {
636  a0 = units.convert(item2.second, "Pa*m^6/kmol^2*K^0.5");
637  } else {
638  auto avec = item2.second.asVector<AnyValue>(2);
639  a0 = units.convert(avec[0], "Pa*m^6/kmol^2*K^0.5");
640  a1 = units.convert(avec[1], "Pa*m^6/kmol^2/K^0.5");
641  }
642  setBinaryCoeffs(item.first, item2.first, a0, a1);
643  }
644  }
645  } else {
646  // Check if a and b are already populated for this species (only the
647  // diagonal elements of a). If not, then search 'critProperties.xml'
648  // to find critical temperature and pressure to calculate a and b.
649  size_t k = speciesIndex(item.first);
650  if (isnan(a_coeff_vec(0, k + m_kk * k))) {
651  // coeffs[0] = a0, coeffs[1] = b;
652  vector<double> coeffs = getCoeff(item.first);
653 
654  // Check if species was found in the database of critical
655  // properties, and assign the results
656  if (!isnan(coeffs[0])) {
657  // Assuming no temperature dependence (i.e. a1 = 0)
658  setSpeciesCoeffs(item.first, coeffs[0], 0.0, coeffs[1]);
659  }
660  }
661  }
662  }
663 }
664 
665 vector<double> RedlichKwongMFTP::getCoeff(const std::string& iName)
666 {
667  vector_fp spCoeff{NAN, NAN};
668 
669  // Get number of species in the database
670  // open xml file critProperties.xml
671  XML_Node* doc = get_XML_File("critProperties.xml");
672  size_t nDatabase = doc->nChildren();
673 
674  // Loop through all species in the database and attempt to match supplied
675  // species to each. If present, calculate pureFluidParameters a_k and b_k
676  // based on crit properties T_c and P_c:
677  for (size_t isp = 0; isp < nDatabase; isp++) {
678  XML_Node& acNodeDoc = doc->child(isp);
679  std::string iNameLower = toLowerCopy(iName);
680  std::string dbName = toLowerCopy(acNodeDoc.attrib("name"));
681 
682  // Attempt to match provided specie iName to current database species
683  // dbName:
684  if (iNameLower == dbName) {
685  // Read from database and calculate a and b coefficients
686  double vParams;
687  double T_crit=0.;
688  double P_crit=0.;
689 
690  if (acNodeDoc.hasChild("Tc")) {
691  vParams = 0.0;
692  XML_Node& xmlChildCoeff = acNodeDoc.child("Tc");
693  if (xmlChildCoeff.hasAttrib("value"))
694  {
695  std::string critTemp = xmlChildCoeff.attrib("value");
696  vParams = strSItoDbl(critTemp);
697  }
698  if (vParams <= 0.0) //Assuming that Pc and Tc are non zero.
699  {
700  throw CanteraError("RedlichKwongMFTP::getCoeff",
701  "Critical Temperature must be positive ");
702  }
703  T_crit = vParams;
704  } else {
705  throw CanteraError("RedlichKwongMFTP::getCoeff",
706  "Critical Temperature not in database ");
707  }
708  if (acNodeDoc.hasChild("Pc")) {
709  vParams = 0.0;
710  XML_Node& xmlChildCoeff = acNodeDoc.child("Pc");
711  if (xmlChildCoeff.hasAttrib("value"))
712  {
713  std::string critPressure = xmlChildCoeff.attrib("value");
714  vParams = strSItoDbl(critPressure);
715  }
716  if (vParams <= 0.0) //Assuming that Pc and Tc are non zero.
717  {
718  throw CanteraError("RedlichKwongMFTP::getCoeff",
719  "Critical Pressure must be positive ");
720  }
721  P_crit = vParams;
722  } else {
723  throw CanteraError("RedlichKwongMFTP::getCoeff",
724  "Critical Pressure not in database ");
725  }
726 
727  //Assuming no temperature dependence
728  spCoeff[0] = omega_a * pow(GasConstant, 2) * pow(T_crit, 2.5) / P_crit; //coeff a
729  spCoeff[1] = omega_b * GasConstant * T_crit / P_crit; // coeff b
730  break;
731  }
732  }
733  return spCoeff;
734 }
735 
737 {
738  string xname = pureFluidParam.name();
739  if (xname != "pureFluidParameters") {
740  throw CanteraError("RedlichKwongMFTP::readXMLPureFluid",
741  "Incorrect name for processing this routine: " + xname);
742  }
743 
744  double a0 = 0.0;
745  double a1 = 0.0;
746  double b = 0.0;
747  for (size_t iChild = 0; iChild < pureFluidParam.nChildren(); iChild++) {
748  XML_Node& xmlChild = pureFluidParam.child(iChild);
749  string nodeName = toLowerCopy(xmlChild.name());
750 
751  if (nodeName == "a_coeff") {
752  vector_fp vParams;
753  string iModel = toLowerCopy(xmlChild.attrib("model"));
754  getFloatArray(xmlChild, vParams, true, "Pascal-m6/kmol2", "a_coeff");
755 
756  if (iModel == "constant" && vParams.size() == 1) {
757  a0 = vParams[0];
758  a1 = 0;
759  } else if (iModel == "linear_a" && vParams.size() == 2) {
760  a0 = vParams[0];
761  a1 = vParams[1];
762  } else {
763  throw CanteraError("RedlichKwongMFTP::readXMLPureFluid",
764  "unknown model or incorrect number of parameters");
765  }
766 
767  } else if (nodeName == "b_coeff") {
768  b = getFloatCurrent(xmlChild, "toSI");
769  }
770  }
771  setSpeciesCoeffs(pureFluidParam.attrib("species"), a0, a1, b);
772 }
773 
775 {
776  string xname = CrossFluidParam.name();
777  if (xname != "crossFluidParameters") {
778  throw CanteraError("RedlichKwongMFTP::readXMLCrossFluid",
779  "Incorrect name for processing this routine: " + xname);
780  }
781 
782  string iName = CrossFluidParam.attrib("species1");
783  string jName = CrossFluidParam.attrib("species2");
784 
785  size_t num = CrossFluidParam.nChildren();
786  for (size_t iChild = 0; iChild < num; iChild++) {
787  XML_Node& xmlChild = CrossFluidParam.child(iChild);
788  string nodeName = toLowerCopy(xmlChild.name());
789 
790  if (nodeName == "a_coeff") {
791  vector_fp vParams;
792  getFloatArray(xmlChild, vParams, true, "Pascal-m6/kmol2", "a_coeff");
793  string iModel = toLowerCopy(xmlChild.attrib("model"));
794  if (iModel == "constant" && vParams.size() == 1) {
795  setBinaryCoeffs(iName, jName, vParams[0], 0.0);
796  } else if (iModel == "linear_a") {
797  setBinaryCoeffs(iName, jName, vParams[0], vParams[1]);
798  } else {
799  throw CanteraError("RedlichKwongMFTP::readXMLCrossFluid",
800  "unknown model ({}) or wrong number of parameters ({})",
801  iModel, vParams.size());
802  }
803  }
804  }
805 }
806 
808 {
810  std::string model = thermoNode["model"];
811 }
812 
813 doublereal RedlichKwongMFTP::sresid() const
814 {
815  // note this agrees with tpx
816  doublereal rho = density();
817  doublereal mmw = meanMolecularWeight();
818  doublereal molarV = mmw / rho;
819  double hh = m_b_current / molarV;
820  doublereal zz = z();
821  doublereal dadt = da_dt();
822  doublereal T = temperature();
823  doublereal sqT = sqrt(T);
824  doublereal fac = dadt - m_a_current / (2.0 * T);
825  double sresid_mol_R = log(zz*(1.0 - hh)) + log(1.0 + hh) * fac / (sqT * GasConstant * m_b_current);
826  return GasConstant * sresid_mol_R;
827 }
828 
829 doublereal RedlichKwongMFTP::hresid() const
830 {
831  // note this agrees with tpx
832  doublereal rho = density();
833  doublereal mmw = meanMolecularWeight();
834  doublereal molarV = mmw / rho;
835  double hh = m_b_current / molarV;
836  doublereal zz = z();
837  doublereal dadt = da_dt();
838  doublereal T = temperature();
839  doublereal sqT = sqrt(T);
840  doublereal fac = T * dadt - 3.0 *m_a_current / (2.0);
841  return GasConstant * T * (zz - 1.0) + fac * log(1.0 + hh) / (sqT * m_b_current);
842 }
843 
844 doublereal RedlichKwongMFTP::liquidVolEst(doublereal TKelvin, doublereal& presGuess) const
845 {
846  double v = m_b_current * 1.1;
847  double atmp;
848  double btmp;
849  calculateAB(TKelvin, atmp, btmp);
850  doublereal pres = std::max(psatEst(TKelvin), presGuess);
851  double Vroot[3];
852  bool foundLiq = false;
853  int m = 0;
854  while (m < 100 && !foundLiq) {
855  int nsol = NicholsSolve(TKelvin, pres, atmp, btmp, Vroot);
856  if (nsol == 1 || nsol == 2) {
857  double pc = critPressure();
858  if (pres > pc) {
859  foundLiq = true;
860  }
861  pres *= 1.04;
862  } else {
863  foundLiq = true;
864  }
865  }
866 
867  if (foundLiq) {
868  v = Vroot[0];
869  presGuess = pres;
870  } else {
871  v = -1.0;
872  }
873  return v;
874 }
875 
876 doublereal RedlichKwongMFTP::densityCalc(doublereal TKelvin, doublereal presPa, int phaseRequested, doublereal rhoguess)
877 {
878  // It's necessary to set the temperature so that m_a_current is set correctly.
879  setTemperature(TKelvin);
880  double tcrit = critTemperature();
881  doublereal mmw = meanMolecularWeight();
882  if (rhoguess == -1.0) {
883  if (phaseRequested != FLUID_GAS) {
884  if (TKelvin > tcrit) {
885  rhoguess = presPa * mmw / (GasConstant * TKelvin);
886  } else {
887  if (phaseRequested == FLUID_GAS || phaseRequested == FLUID_SUPERCRIT) {
888  rhoguess = presPa * mmw / (GasConstant * TKelvin);
889  } else if (phaseRequested >= FLUID_LIQUID_0) {
890  double lqvol = liquidVolEst(TKelvin, presPa);
891  rhoguess = mmw / lqvol;
892  }
893  }
894  } else {
895  // Assume the Gas phase initial guess, if nothing is specified to
896  // the routine
897  rhoguess = presPa * mmw / (GasConstant * TKelvin);
898  }
899  }
900 
901  doublereal volguess = mmw / rhoguess;
902  NSolns_ = NicholsSolve(TKelvin, presPa, m_a_current, m_b_current, Vroot_);
903 
904  doublereal molarVolLast = Vroot_[0];
905  if (NSolns_ >= 2) {
906  if (phaseRequested >= FLUID_LIQUID_0) {
907  molarVolLast = Vroot_[0];
908  } else if (phaseRequested == FLUID_GAS || phaseRequested == FLUID_SUPERCRIT) {
909  molarVolLast = Vroot_[2];
910  } else {
911  if (volguess > Vroot_[1]) {
912  molarVolLast = Vroot_[2];
913  } else {
914  molarVolLast = Vroot_[0];
915  }
916  }
917  } else if (NSolns_ == 1) {
918  if (phaseRequested == FLUID_GAS || phaseRequested == FLUID_SUPERCRIT || phaseRequested == FLUID_UNDEFINED) {
919  molarVolLast = Vroot_[0];
920  } else {
921  return -2.0;
922  }
923  } else if (NSolns_ == -1) {
924  if (phaseRequested >= FLUID_LIQUID_0 || phaseRequested == FLUID_UNDEFINED || phaseRequested == FLUID_SUPERCRIT) {
925  molarVolLast = Vroot_[0];
926  } else if (TKelvin > tcrit) {
927  molarVolLast = Vroot_[0];
928  } else {
929  return -2.0;
930  }
931  } else {
932  molarVolLast = Vroot_[0];
933  return -1.0;
934  }
935  return mmw / molarVolLast;
936 }
937 
939 {
940  double Vroot[3];
941  double T = temperature();
942  int nsol = NicholsSolve(T, pressure(), m_a_current, m_b_current, Vroot);
943  if (nsol != 3) {
944  return critDensity();
945  }
946 
947  auto resid = [this, T](double v) {
948  double pp;
949  return dpdVCalc(T, v, pp);
950  };
951 
952  boost::uintmax_t maxiter = 100;
953  std::pair<double, double> vv = bmt::toms748_solve(
954  resid, Vroot[0], Vroot[1], bmt::eps_tolerance<double>(48), maxiter);
955 
956  doublereal mmw = meanMolecularWeight();
957  return mmw / (0.5 * (vv.first + vv.second));
958 }
959 
961 {
962  double Vroot[3];
963  double T = temperature();
964  int nsol = NicholsSolve(T, pressure(), m_a_current, m_b_current, Vroot);
965  if (nsol != 3) {
966  return critDensity();
967  }
968 
969  auto resid = [this, T](double v) {
970  double pp;
971  return dpdVCalc(T, v, pp);
972  };
973 
974  boost::uintmax_t maxiter = 100;
975  std::pair<double, double> vv = bmt::toms748_solve(
976  resid, Vroot[1], Vroot[2], bmt::eps_tolerance<double>(48), maxiter);
977 
978  doublereal mmw = meanMolecularWeight();
979  return mmw / (0.5 * (vv.first + vv.second));
980 }
981 
982 doublereal RedlichKwongMFTP::pressureCalc(doublereal TKelvin, doublereal molarVol) const
983 {
984  doublereal sqt = sqrt(TKelvin);
985  double pres = GasConstant * TKelvin / (molarVol - m_b_current)
986  - m_a_current / (sqt * molarVol * (molarVol + m_b_current));
987  return pres;
988 }
989 
990 doublereal RedlichKwongMFTP::dpdVCalc(doublereal TKelvin, doublereal molarVol, doublereal& presCalc) const
991 {
992  doublereal sqt = sqrt(TKelvin);
993  presCalc = GasConstant * TKelvin / (molarVol - m_b_current)
994  - m_a_current / (sqt * molarVol * (molarVol + m_b_current));
995 
996  doublereal vpb = molarVol + m_b_current;
997  doublereal vmb = molarVol - m_b_current;
998  doublereal dpdv = (- GasConstant * TKelvin / (vmb * vmb)
999  + m_a_current * (2 * molarVol + m_b_current) / (sqt * molarVol * molarVol * vpb * vpb));
1000  return dpdv;
1001 }
1002 
1004 {
1005  doublereal TKelvin = temperature();
1006  doublereal mv = molarVolume();
1007  doublereal pres;
1008 
1009  dpdV_ = dpdVCalc(TKelvin, mv, pres);
1010  doublereal sqt = sqrt(TKelvin);
1011  doublereal vpb = mv + m_b_current;
1012  doublereal vmb = mv - m_b_current;
1013  doublereal dadt = da_dt();
1014  doublereal fac = dadt - m_a_current/(2.0 * TKelvin);
1015  dpdT_ = (GasConstant / vmb - fac / (sqt * mv * vpb));
1016 }
1017 
1018 void RedlichKwongMFTP::updateMixingExpressions()
1019 {
1020  updateAB();
1021 }
1022 
1024 {
1025  double temp = temperature();
1026  if (m_formTempParam == 1) {
1027  for (size_t i = 0; i < m_kk; i++) {
1028  for (size_t j = 0; j < m_kk; j++) {
1029  size_t counter = i * m_kk + j;
1030  a_vec_Curr_[counter] = a_coeff_vec(0,counter) + a_coeff_vec(1,counter) * temp;
1031  }
1032  }
1033  }
1034 
1035  m_b_current = 0.0;
1036  m_a_current = 0.0;
1037  for (size_t i = 0; i < m_kk; i++) {
1038  m_b_current += moleFractions_[i] * b_vec_Curr_[i];
1039  for (size_t j = 0; j < m_kk; j++) {
1040  m_a_current += a_vec_Curr_[i * m_kk + j] * moleFractions_[i] * moleFractions_[j];
1041  }
1042  }
1043  if (isnan(m_b_current)) {
1044  // One or more species do not have specified coefficients.
1045  fmt::memory_buffer b;
1046  for (size_t k = 0; k < m_kk; k++) {
1047  if (isnan(b_vec_Curr_[k])) {
1048  if (b.size() > 0) {
1049  format_to(b, ", {}", speciesName(k));
1050  } else {
1051  format_to(b, "{}", speciesName(k));
1052  }
1053  }
1054  }
1055  throw CanteraError("RedlichKwongMFTP::updateAB",
1056  "Missing Redlich-Kwong coefficients for species: {}", to_string(b));
1057  }
1058 }
1059 
1060 void RedlichKwongMFTP::calculateAB(doublereal temp, doublereal& aCalc, doublereal& bCalc) const
1061 {
1062  bCalc = 0.0;
1063  aCalc = 0.0;
1064  if (m_formTempParam == 1) {
1065  for (size_t i = 0; i < m_kk; i++) {
1066  bCalc += moleFractions_[i] * b_vec_Curr_[i];
1067  for (size_t j = 0; j < m_kk; j++) {
1068  size_t counter = i * m_kk + j;
1069  doublereal a_vec_Curr = a_coeff_vec(0,counter) + a_coeff_vec(1,counter) * temp;
1070  aCalc += a_vec_Curr * moleFractions_[i] * moleFractions_[j];
1071  }
1072  }
1073  } else {
1074  for (size_t i = 0; i < m_kk; i++) {
1075  bCalc += moleFractions_[i] * b_vec_Curr_[i];
1076  for (size_t j = 0; j < m_kk; j++) {
1077  size_t counter = i * m_kk + j;
1078  doublereal a_vec_Curr = a_coeff_vec(0,counter);
1079  aCalc += a_vec_Curr * moleFractions_[i] * moleFractions_[j];
1080  }
1081  }
1082  }
1083 }
1084 
1085 doublereal RedlichKwongMFTP::da_dt() const
1086 {
1087  doublereal dadT = 0.0;
1088  if (m_formTempParam == 1) {
1089  for (size_t i = 0; i < m_kk; i++) {
1090  for (size_t j = 0; j < m_kk; j++) {
1091  size_t counter = i * m_kk + j;
1092  dadT+= a_coeff_vec(1,counter) * moleFractions_[i] * moleFractions_[j];
1093  }
1094  }
1095  }
1096  return dadT;
1097 }
1098 
1099 void RedlichKwongMFTP::calcCriticalConditions(doublereal a, doublereal b, doublereal a0_coeff, doublereal aT_coeff,
1100  doublereal& pc, doublereal& tc, doublereal& vc) const
1101 {
1102  if (m_formTempParam != 0) {
1103  a = a0_coeff;
1104  }
1105  if (b <= 0.0) {
1106  tc = 1000000.;
1107  pc = 1.0E13;
1108  vc = omega_vc * GasConstant * tc / pc;
1109  return;
1110  }
1111  if (a <= 0.0) {
1112  tc = 0.0;
1113  pc = 0.0;
1114  vc = 2.0 * b;
1115  return;
1116  }
1117  double tmp = a * omega_b / (b * omega_a * GasConstant);
1118  double pp = 2./3.;
1119  doublereal sqrttc, f, dfdt, deltatc;
1120 
1121  if (m_formTempParam == 0) {
1122  tc = pow(tmp, pp);
1123  } else {
1124  tc = pow(tmp, pp);
1125  for (int j = 0; j < 10; j++) {
1126  sqrttc = sqrt(tc);
1127  f = omega_a * b * GasConstant * tc * sqrttc / omega_b - aT_coeff * tc - a0_coeff;
1128  dfdt = 1.5 * omega_a * b * GasConstant * sqrttc / omega_b - aT_coeff;
1129  deltatc = - f / dfdt;
1130  tc += deltatc;
1131  }
1132  if (deltatc > 0.1) {
1133  throw CanteraError("RedlichKwongMFTP::calcCriticalConditions", "didn't converge");
1134  }
1135  }
1136 
1137  pc = omega_b * GasConstant * tc / b;
1138  vc = omega_vc * GasConstant * tc / pc;
1139 }
1140 
1141 int RedlichKwongMFTP::NicholsSolve(double TKelvin, double pres, doublereal a, doublereal b,
1142  doublereal Vroot[3]) const
1143 {
1144  Vroot[0] = 0.0;
1145  Vroot[1] = 0.0;
1146  Vroot[2] = 0.0;
1147  if (TKelvin <= 0.0) {
1148  throw CanteraError("RedlichKwongMFTP::NicholsSolve", "neg temperature");
1149  }
1150 
1151  // Derive the coefficients of the cubic polynomial to solve.
1152  doublereal an = 1.0;
1153  doublereal bn = - GasConstant * TKelvin / pres;
1154  doublereal sqt = sqrt(TKelvin);
1155  doublereal cn = - (GasConstant * TKelvin * b / pres - a/(pres * sqt) + b * b);
1156  doublereal dn = - (a * b / (pres * sqt));
1157 
1158  double tmp = a * omega_b / (b * omega_a * GasConstant);
1159  double pp = 2./3.;
1160  double tc = pow(tmp, pp);
1161  double pc = omega_b * GasConstant * tc / b;
1162  double vc = omega_vc * GasConstant * tc / pc;
1163  // Derive the center of the cubic, x_N
1164  doublereal xN = - bn /(3 * an);
1165 
1166  // Derive the value of delta**2. This is a key quantity that determines the
1167  // number of turning points
1168  doublereal delta2 = (bn * bn - 3 * an * cn) / (9 * an * an);
1169  doublereal delta = 0.0;
1170 
1171  // Calculate a couple of ratios
1172  doublereal ratio1 = 3.0 * an * cn / (bn * bn);
1173  doublereal ratio2 = pres * b / (GasConstant * TKelvin);
1174  if (fabs(ratio1) < 1.0E-7) {
1175  doublereal ratio3 = a / (GasConstant * sqt) * pres / (GasConstant * TKelvin);
1176  if (fabs(ratio2) < 1.0E-5 && fabs(ratio3) < 1.0E-5) {
1177  doublereal zz = 1.0;
1178  for (int i = 0; i < 10; i++) {
1179  doublereal znew = zz / (zz - ratio2) - ratio3 / (zz + ratio1);
1180  doublereal deltaz = znew - zz;
1181  zz = znew;
1182  if (fabs(deltaz) < 1.0E-14) {
1183  break;
1184  }
1185  }
1186  doublereal v = zz * GasConstant * TKelvin / pres;
1187  Vroot[0] = v;
1188  return 1;
1189  }
1190  }
1191 
1192  int nSolnValues;
1193  double h2 = 4. * an * an * delta2 * delta2 * delta2;
1194  if (delta2 > 0.0) {
1195  delta = sqrt(delta2);
1196  }
1197 
1198  doublereal h = 2.0 * an * delta * delta2;
1199  doublereal yN = 2.0 * bn * bn * bn / (27.0 * an * an) - bn * cn / (3.0 * an) + dn;
1200  doublereal desc = yN * yN - h2;
1201 
1202  if (fabs(fabs(h) - fabs(yN)) < 1.0E-10) {
1203  if (desc != 0.0) {
1204  // this is for getting to other cases
1205  throw CanteraError("RedlichKwongMFTP::NicholsSolve", "numerical issues");
1206  }
1207  desc = 0.0;
1208  }
1209 
1210  if (desc < 0.0) {
1211  nSolnValues = 3;
1212  } else if (desc == 0.0) {
1213  nSolnValues = 2;
1214  // We are here as p goes to zero.
1215  } else if (desc > 0.0) {
1216  nSolnValues = 1;
1217  }
1218 
1219  // One real root -> have to determine whether gas or liquid is the root
1220  if (desc > 0.0) {
1221  doublereal tmpD = sqrt(desc);
1222  doublereal tmp1 = (- yN + tmpD) / (2.0 * an);
1223  doublereal sgn1 = 1.0;
1224  if (tmp1 < 0.0) {
1225  sgn1 = -1.0;
1226  tmp1 = -tmp1;
1227  }
1228  doublereal tmp2 = (- yN - tmpD) / (2.0 * an);
1229  doublereal sgn2 = 1.0;
1230  if (tmp2 < 0.0) {
1231  sgn2 = -1.0;
1232  tmp2 = -tmp2;
1233  }
1234  doublereal p1 = pow(tmp1, 1./3.);
1235  doublereal p2 = pow(tmp2, 1./3.);
1236  doublereal alpha = xN + sgn1 * p1 + sgn2 * p2;
1237  Vroot[0] = alpha;
1238  Vroot[1] = 0.0;
1239  Vroot[2] = 0.0;
1240  tmp = an * Vroot[0] * Vroot[0] * Vroot[0] + bn * Vroot[0] * Vroot[0] + cn * Vroot[0] + dn;
1241  } else if (desc < 0.0) {
1242  doublereal tmp = - yN/h;
1243  doublereal val = acos(tmp);
1244  doublereal theta = val / 3.0;
1245  doublereal oo = 2. * Pi / 3.;
1246  doublereal alpha = xN + 2. * delta * cos(theta);
1247  doublereal beta = xN + 2. * delta * cos(theta + oo);
1248  doublereal gamma = xN + 2. * delta * cos(theta + 2.0 * oo);
1249  Vroot[0] = beta;
1250  Vroot[1] = gamma;
1251  Vroot[2] = alpha;
1252 
1253  for (int i = 0; i < 3; i++) {
1254  tmp = an * Vroot[i] * Vroot[i] * Vroot[i] + bn * Vroot[i] * Vroot[i] + cn * Vroot[i] + dn;
1255  if (fabs(tmp) > 1.0E-4) {
1256  for (int j = 0; j < 3; j++) {
1257  if (j != i && fabs(Vroot[i] - Vroot[j]) < 1.0E-4 * (fabs(Vroot[i]) + fabs(Vroot[j]))) {
1258  warn_user("RedlichKwongMFTP::NicholsSolve",
1259  "roots have merged: {}, {} (T = {}, p = {})",
1260  Vroot[i], Vroot[j], TKelvin, pres);
1261  }
1262  }
1263  }
1264  }
1265  } else if (desc == 0.0) {
1266  if (yN == 0.0 && h == 0.0) {
1267  Vroot[0] = xN;
1268  Vroot[1] = xN;
1269  Vroot[2] = xN;
1270  } else {
1271  // need to figure out whether delta is pos or neg
1272  if (yN > 0.0) {
1273  tmp = pow(yN/(2*an), 1./3.);
1274  if (fabs(tmp - delta) > 1.0E-9) {
1275  throw CanteraError("RedlichKwongMFTP::NicholsSolve", "unexpected");
1276  }
1277  Vroot[1] = xN + delta;
1278  Vroot[0] = xN - 2.0*delta; // liquid phase root
1279  } else {
1280  tmp = pow(yN/(2*an), 1./3.);
1281  if (fabs(tmp - delta) > 1.0E-9) {
1282  throw CanteraError("RedlichKwongMFTP::NicholsSolve", "unexpected");
1283  }
1284  delta = -delta;
1285  Vroot[0] = xN + delta;
1286  Vroot[1] = xN - 2.0*delta; // gas phase root
1287  }
1288  }
1289  for (int i = 0; i < 2; i++) {
1290  tmp = an * Vroot[i] * Vroot[i] * Vroot[i] + bn * Vroot[i] * Vroot[i] + cn * Vroot[i] + dn;
1291  }
1292  }
1293 
1294  // Unfortunately, there is a heavy amount of roundoff error due to bad
1295  // conditioning in this
1296  double res, dresdV = 0.0;
1297  for (int i = 0; i < nSolnValues; i++) {
1298  for (int n = 0; n < 20; n++) {
1299  res = an * Vroot[i] * Vroot[i] * Vroot[i] + bn * Vroot[i] * Vroot[i] + cn * Vroot[i] + dn;
1300  if (fabs(res) < 1.0E-14) {
1301  break;
1302  }
1303  dresdV = 3.0 * an * Vroot[i] * Vroot[i] + 2.0 * bn * Vroot[i] + cn;
1304  double del = - res / dresdV;
1305  Vroot[i] += del;
1306  if (fabs(del) / (fabs(Vroot[i]) + fabs(del)) < 1.0E-14) {
1307  break;
1308  }
1309  double res2 = an * Vroot[i] * Vroot[i] * Vroot[i] + bn * Vroot[i] * Vroot[i] + cn * Vroot[i] + dn;
1310  if (fabs(res2) < fabs(res)) {
1311  continue;
1312  } else {
1313  Vroot[i] -= del;
1314  Vroot[i] += 0.1 * del;
1315  }
1316  }
1317  if ((fabs(res) > 1.0E-14) && (fabs(res) > 1.0E-14 * fabs(dresdV) * fabs(Vroot[i]))) {
1318  warn_user("RedlichKwongMFTP::NicholsSolve",
1319  "root did not converge: V = {} (T = {}, p = {})",
1320  Vroot[i], TKelvin, pres);
1321  }
1322  }
1323 
1324  if (nSolnValues == 1) {
1325  if (TKelvin > tc) {
1326  if (Vroot[0] < vc) {
1327  nSolnValues = -1;
1328  }
1329  } else {
1330  if (Vroot[0] < xN) {
1331  nSolnValues = -1;
1332  }
1333  }
1334  } else {
1335  if (nSolnValues == 2 && delta > 0.0) {
1336  nSolnValues = -2;
1337  }
1338  }
1339  return nSolnValues;
1340 }
1341 
1342 }
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:360
const UnitSystem & units() const
Return the default units that should be used to convert stored values.
Definition: AnyMap.h:492
A wrapper for a variable whose type is determined at runtime.
Definition: AnyMap.h:77
void resize(size_t n, size_t m, doublereal v=0.0)
Resize the array, and fill the new entries with 'v'.
Definition: Array.h:112
vector_fp & data()
Return a reference to the data vector.
Definition: Array.h:277
void getRow(size_t n, doublereal *const rw)
Get the nth row and return it in a vector.
Definition: Array.h:165
Base class for exceptions thrown by Cantera classes.
Definition: ctexceptions.h:61
vector_fp m_cp0_R
Temporary storage for dimensionless reference state heat capacities.
virtual bool addSpecies(shared_ptr< Species > spec)
doublereal z() const
Calculate the value of z.
vector_fp m_h0_RT
Temporary storage for dimensionless reference state enthalpies.
virtual void _updateReferenceStateThermo() const
Updates the reference state thermodynamic functions at the current T of the solution.
virtual void getStandardVolumes(doublereal *vol) const
Get the molar volumes of each species in their standard states at the current T and P of the solution...
vector_fp m_s0_R
Temporary storage for dimensionless reference state entropies.
virtual void getCp_R(doublereal *cpr) const
Get the nondimensional Heat Capacities at constant pressure for the standard state of the species at ...
virtual doublereal psatEst(doublereal TKelvin) const
Estimate for the saturation pressure.
virtual void getEntropy_R_ref(doublereal *er) const
Returns the vector of nondimensional entropies of the reference state at the current temperature of t...
virtual void getIntEnergy_RT(doublereal *urt) const
Returns the vector of nondimensional internal Energies of the standard state at the current temperatu...
vector_fp moleFractions_
Storage for the current values of the mole fractions of the species.
virtual void compositionChanged()
Apply changes to the state which are needed after the composition changes.
virtual void getGibbs_ref(doublereal *g) const
Returns the vector of the Gibbs function of the reference state at the current temperature of the sol...
virtual void getEnthalpy_RT_ref(doublereal *hrt) const
doublereal mean_X(const doublereal *const Q) const
Evaluate the mole-fraction-weighted mean of an array Q.
Definition: Phase.cpp:746
const double * moleFractdivMMW() const
Returns a const pointer to the start of the moleFraction/MW array.
Definition: Phase.cpp:593
size_t m_kk
Number of species in the phase.
Definition: Phase.h:942
std::string speciesName(size_t k) const
Name of the species with index k.
Definition: Phase.cpp:229
virtual void setDensity(const double density_)
Set the internally stored density (kg/m^3) of the phase.
Definition: Phase.cpp:716
doublereal meanMolecularWeight() const
The mean molecular weight. Units: (kg/kmol)
Definition: Phase.h:748
double moleFraction(size_t k) const
Return the mole fraction of a single species.
Definition: Phase.cpp:577
virtual double density() const
Density (kg/m^3).
Definition: Phase.h:685
doublereal temperature() const
Temperature (K).
Definition: Phase.h:667
size_t speciesIndex(const std::string &name) const
Returns the index of a species named 'name' within the Phase object.
Definition: Phase.cpp:201
virtual void setTemperature(const doublereal temp)
Set the internally stored temperature of the phase (K).
Definition: Phase.h:724
double molarVolume() const
Molar volume (m^3/kmol).
Definition: Phase.cpp:711
shared_ptr< Species > species(const std::string &name) const
Return the Species object for the named species.
Definition: Phase.cpp:980
doublereal sum_xlogx() const
Evaluate .
Definition: Phase.cpp:756
virtual bool addSpecies(shared_ptr< Species > spec)
virtual void getActivityCoefficients(doublereal *ac) const
Get the array of non-dimensional activity coefficients at the current solution temperature,...
void updateAB()
Update the a and b parameters.
static const doublereal omega_b
Omega constant for b.
int m_formTempParam
Form of the temperature parameterization.
virtual doublereal hresid() const
Calculate the deviation terms for the total enthalpy of the mixture from the ideal gas mixture.
doublereal m_b_current
Value of b in the equation of state.
doublereal dpdT_
The derivative of the pressure wrt the temperature.
virtual doublereal pressure() const
Return the thermodynamic pressure (Pa).
virtual void getPartialMolarIntEnergies(doublereal *ubar) const
Return an array of partial molar internal energies for the species in the mixture.
virtual doublereal critVolume() const
Critical volume (m3/kmol).
virtual doublereal cp_mole() const
Molar heat capacity at constant pressure. Units: J/kmol/K.
virtual std::vector< double > getCoeff(const std::string &iName)
Retrieve a and b coefficients by looking up tabulated critical parameters.
virtual void getPartialMolarEnthalpies(doublereal *hbar) const
Returns an array of partial molar enthalpies for the species in the mixture.
virtual void getPartialMolarEntropies(doublereal *sbar) const
Returns an array of partial molar entropies of the species in the solution.
virtual void initThermoXML(XML_Node &phaseNode, const std::string &id)
Import and initialize a ThermoPhase object using an XML tree.
virtual doublereal critPressure() const
Critical pressure (Pa).
virtual doublereal enthalpy_mole() const
Molar enthalpy. Units: J/kmol.
virtual doublereal densityCalc(doublereal TKelvin, doublereal pressure, int phase, doublereal rhoguess)
Calculates the density given the temperature and the pressure and a guess at the density.
void readXMLCrossFluid(XML_Node &pureFluidParam)
Read the cross species RedlichKwong input parameters.
virtual void getPartialMolarVolumes(doublereal *vbar) const
Return an array of partial molar volumes for the species in the mixture.
vector_fp m_pp
Temporary storage - length = m_kk.
virtual doublereal cv_mole() const
Molar heat capacity at constant volume. Units: J/kmol/K.
virtual void getActivityConcentrations(doublereal *c) const
This method returns an array of generalized concentrations.
void readXMLPureFluid(XML_Node &pureFluidParam)
Read the pure species RedlichKwong input parameters.
virtual doublereal dpdVCalc(doublereal TKelvin, doublereal molarVol, doublereal &presCalc) const
Calculate the pressure and the pressure derivative given the temperature and the molar volume.
virtual doublereal sresid() const
Calculate the deviation terms for the total entropy of the mixture from the ideal gas mixture.
virtual void getPartialMolarCp(doublereal *cpbar) const
Return an array of partial molar heat capacities for the species in the mixture.
void setBinaryCoeffs(const std::string &species_i, const std::string &species_j, double a0, double a1)
Set values for the interaction parameter between two species.
void pressureDerivatives() const
Calculate dpdV and dpdT at the current conditions.
virtual doublereal entropy_mole() const
Molar entropy. Units: J/kmol/K.
static const doublereal omega_vc
Omega constant for the critical molar volume.
virtual doublereal liquidVolEst(doublereal TKelvin, doublereal &pres) const
Estimate for the molar volume of the liquid.
virtual void initThermo()
Initialize the ThermoPhase object after all species have been set up.
virtual doublereal critTemperature() const
Critical temperature (K).
doublereal m_a_current
Value of a in the equation of state.
virtual doublereal densSpinodalGas() const
Return the value of the density at the gas spinodal point (on the gas side) for the current temperatu...
void setSpeciesCoeffs(const std::string &species, double a0, double a1, double b)
Set the pure fluid interaction parameters for a species.
RedlichKwongMFTP()
Base constructor.
virtual void setTemperature(const doublereal temp)
Set the temperature of the phase.
int NicholsSolve(double TKelvin, double pres, doublereal a, doublereal b, doublereal Vroot[3]) const
Solve the cubic equation of state.
virtual doublereal pressureCalc(doublereal TKelvin, doublereal molarVol) const
Calculate the pressure given the temperature and the molar volume.
virtual void compositionChanged()
Apply changes to the state which are needed after the composition changes.
void calculateAB(doublereal temp, doublereal &aCalc, doublereal &bCalc) const
Calculate the a and the b parameters given the temperature.
virtual void getChemPotentials_RT(doublereal *mu) const
Get the array of non-dimensional species chemical potentials.
doublereal dpdV_
The derivative of the pressure wrt the volume.
virtual doublereal critCompressibility() const
Critical compressibility (unitless).
virtual doublereal densSpinodalLiquid() const
Return the value of the density at the liquid spinodal point (on the liquid side) for the current tem...
virtual void getChemPotentials(doublereal *mu) const
Get the species chemical potentials. Units: J/kmol.
virtual doublereal critDensity() const
Critical density (kg/m3).
vector_fp m_tmpV
Temporary storage - length = m_kk.
static const doublereal omega_a
Omega constant for a -> value of a in terms of critical properties.
vector_fp dpdni_
Vector of derivatives of pressure wrt mole number.
virtual void setParametersFromXML(const XML_Node &thermoNode)
Set equation of state parameter values from XML entries.
virtual void calcDensity()
Calculate the density of the mixture using the partial molar volumes and mole fractions as input.
virtual doublereal standardConcentration(size_t k=0) const
Returns the standard concentration , which is used to normalize the generalized concentration.
doublereal RT() const
Return the Gas Constant multiplied by the current temperature.
Definition: ThermoPhase.h:776
virtual void initThermoFile(const std::string &inputFile, const std::string &id)
virtual void initThermoXML(XML_Node &phaseNode, const std::string &id)
Import and initialize a ThermoPhase object using an XML tree.
virtual doublereal refPressure() const
Returns the reference pressure in Pa.
Definition: ThermoPhase.h:145
virtual void setParametersFromXML(const XML_Node &eosdata)
Set equation of state parameter values from XML entries.
Definition: ThermoPhase.h:1728
Unit conversion utility.
Definition: Units.h:101
double convert(double value, const std::string &src, const std::string &dest) const
Convert value from the units of src to the units of dest.
Definition: Units.cpp:322
Class XML_Node is a tree-based representation of the contents of an XML file.
Definition: xml.h:104
std::string attrib(const std::string &attr) const
Function returns the value of an attribute.
Definition: xml.cpp:492
std::string name() const
Returns the name of the XML node.
Definition: xml.h:372
bool hasChild(const std::string &ch) const
Tests whether the current node has a child node with a particular name.
Definition: xml.cpp:528
size_t nChildren(bool discardComments=false) const
Return the number of children.
Definition: xml.cpp:556
XML_Node & child(const size_t n) const
Return a changeable reference to the n'th child of the current node.
Definition: xml.cpp:546
bool hasAttrib(const std::string &a) const
Tests whether the current node has an attribute with a particular name.
Definition: xml.cpp:533
CTML ("Cantera Markup Language") is the variant of XML that Cantera uses to store data.
XML_Node * get_XML_File(const std::string &file, int debug)
Return a pointer to the XML tree for a Cantera input file.
Definition: global.cpp:110
void warn_user(const std::string &method, const std::string &msg, const Args &... args)
Definition: global.h:206
const size_t npos
index returned by functions to indicate "no position"
Definition: ct_defs.h:188
const double SmallNumber
smallest number to compare to zero.
Definition: ct_defs.h:149
std::vector< double > vector_fp
Turn on the use of stl vectors for the basic array type within cantera Vector of doubles.
Definition: ct_defs.h:180
const double GasConstant
Universal Gas Constant [J/kmol/K].
Definition: ct_defs.h:109
const double Pi
Pi.
Definition: ct_defs.h:53
void importPhase(XML_Node &phase, ThermoPhase *th)
Import a phase information into an empty ThermoPhase object.
Namespace for the Cantera kernel.
Definition: AnyMap.cpp:264
doublereal getFloatCurrent(const XML_Node &node, const std::string &type)
Get a floating-point value from the current XML element.
Definition: ctml.cpp:177
doublereal dot(InputIter x_begin, InputIter x_end, InputIter2 y_begin)
Function that calculates a templated inner product.
Definition: utilities.h:112
bool caseInsensitiveEquals(const std::string &input, const std::string &test)
Case insensitive equality predicate.
std::string toLowerCopy(const std::string &input)
Convert to lower case.
void scale(InputIter begin, InputIter end, OutputIter out, S scale_factor)
Multiply elements of an array by a scale factor.
Definition: utilities.h:135
size_t getFloatArray(const XML_Node &node, vector_fp &v, const bool convert, const std::string &unitsString, const std::string &nodeName)
This function reads the current node or a child node of the current node with the default name,...
Definition: ctml.cpp:256
doublereal strSItoDbl(const std::string &strSI)
Interpret one or two token string as a single double.
Contains declarations for string manipulation functions within Cantera.