Cantera  2.1.2
ChemEquil.cpp
Go to the documentation of this file.
1 /**
2  * @file ChemEquil.cpp
3  * Chemical equilibrium. Implementation file for class
4  * ChemEquil.
5  */
6 
7 // Copyright 2001 California Institute of Technology
8 
11 
12 #include "cantera/base/ct_defs.h"
13 #include "cantera/base/global.h"
14 #include "PropertyCalculator.h"
19 
20 using namespace std;
21 
22 #include <cstdio>
23 
24 int Cantera::ChemEquil_print_lvl = 0;
25 
26 namespace Cantera
27 {
28 
29 int _equilflag(const char* xy)
30 {
31  string flag = string(xy);
32  if (flag == "TP") {
33  return TP;
34  } else if (flag == "TV") {
35  return TV;
36  } else if (flag == "HP") {
37  return HP;
38  } else if (flag == "UV") {
39  return UV;
40  } else if (flag == "SP") {
41  return SP;
42  } else if (flag == "SV") {
43  return SV;
44  } else if (flag == "UP") {
45  return UP;
46  } else {
47  throw CanteraError("_equilflag","unknown property pair "+flag);
48  }
49  return -1;
50 }
51 
52 ChemEquil::ChemEquil() : m_skip(npos), m_elementTotalSum(1.0),
53  m_p0(OneAtm), m_eloc(npos),
54  m_elemFracCutoff(1.0E-100),
55  m_doResPerturb(false)
56 {}
57 
58 ChemEquil::ChemEquil(thermo_t& s) :
59  m_skip(npos),
60  m_elementTotalSum(1.0),
61  m_p0(OneAtm), m_eloc(npos),
62  m_elemFracCutoff(1.0E-100),
63  m_doResPerturb(false)
64 {
65  initialize(s);
66 }
67 
68 ChemEquil::~ChemEquil()
69 {
70 }
71 
73 {
74  // store a pointer to s and some of its properties locally.
75  m_phase = &s;
76 
77  m_p0 = s.refPressure();
78  m_kk = s.nSpecies();
79  m_mm = s.nElements();
81 
82  // allocate space in internal work arrays within the ChemEquil object
83  m_molefractions.resize(m_kk);
84  m_lambda.resize(m_mm, -100.0);
85  m_elementmolefracs.resize(m_mm);
86  m_comp.resize(m_mm * m_kk);
87  m_jwork1.resize(m_mm+2);
88  m_jwork2.resize(m_mm+2);
89  m_startSoln.resize(m_mm+1);
90  m_grt.resize(m_kk);
91  m_mu_RT.resize(m_kk);
92  m_muSS_RT.resize(m_kk);
93  m_component.resize(m_mm,npos);
94  m_orderVectorElements.resize(m_mm);
95 
96  for (size_t m = 0; m < m_mm; m++) {
97  m_orderVectorElements[m] = m;
98  }
99  m_orderVectorSpecies.resize(m_kk);
100  for (size_t k = 0; k < m_kk; k++) {
101  m_orderVectorSpecies[k] = k;
102  }
103 
104  // set up elemental composition matrix
105  size_t mneg = npos;
106  doublereal na, ewt;
107  for (size_t m = 0; m < m_mm; m++) {
108  for (size_t k = 0; k < m_kk; k++) {
109  na = s.nAtoms(k,m);
110 
111  // handle the case of negative atom numbers (used to
112  // represent positive ions, where the 'element' is an
113  // electron
114  if (na < 0.0) {
115 
116  // if negative atom numbers have already been specified
117  // for some element other than this one, throw
118  // an exception
119  if (mneg != npos && mneg != m)
120  throw CanteraError("ChemEquil::initialize",
121  "negative atom numbers allowed for only one element");
122  mneg = m;
123  ewt = s.atomicWeight(m);
124 
125  // the element should be an electron... if it isn't
126  // print a warning.
127  if (ewt > 1.0e-3)
128  writelog(string("WARNING: species "
129  +s.speciesName(k)
130  +" has "+fp2str(s.nAtoms(k,m))
131  +" atoms of element "
132  +s.elementName(m)+
133  ", but this element is not an electron.\n"));
134  }
135  }
136  }
137  m_eloc = mneg;
138 
139  // set up the elemental composition matrix
140  for (size_t k = 0; k < m_kk; k++) {
141  for (size_t m = 0; m < m_mm; m++) {
142  m_comp[k*m_mm + m] = s.nAtoms(k,m);
143  }
144  }
145 }
146 
148  const vector_fp& lambda_RT, doublereal t)
149 {
150  // Construct the chemical potentials by summing element potentials
151  fill(m_mu_RT.begin(), m_mu_RT.end(), 0.0);
152  for (size_t k = 0; k < m_kk; k++)
153  for (size_t m = 0; m < m_mm; m++) {
154  m_mu_RT[k] += lambda_RT[m]*nAtoms(k,m);
155  }
156 
157  // Set the temperature
158  s.setTemperature(t);
159 
160  // Call the phase-specific method to set the phase to the
161  // equilibrium state with the specified species chemical
162  // potentials.
163  s.setToEquilState(DATA_PTR(m_mu_RT));
164  update(s);
165 }
166 
168 {
169 
170  // get the mole fractions, temperature, and density
172  m_temp = s.temperature();
173  m_dens = s.density();
174 
175  // compute the elemental mole fractions
176  double sum = 0.0;
177  for (size_t m = 0; m < m_mm; m++) {
178  m_elementmolefracs[m] = 0.0;
179  for (size_t k = 0; k < m_kk; k++) {
180  m_elementmolefracs[m] += nAtoms(k,m) * m_molefractions[k];
181  if (m_molefractions[k] < 0.0) {
182  throw CanteraError("update",
183  "negative mole fraction for "+s.speciesName(k)+
184  ": "+fp2str(m_molefractions[k]));
185  }
186  }
187  sum += m_elementmolefracs[m];
188  }
189  // Store the sum for later use
190  m_elementTotalSum = sum;
191  // normalize the element mole fractions
192  for (size_t m = 0; m < m_mm; m++) {
193  m_elementmolefracs[m] /= sum;
194  }
195 }
196 
198  int loglevel)
199 {
200  int iok = 0;
201  if (loglevel > 0) {
202  beginLogGroup("ChemEquil::setInitialMoles");
203  }
204  try {
205  MultiPhase mp;
206  mp.addPhase(&s, 1.0);
207  mp.init();
208  MultiPhaseEquil e(&mp, true, loglevel-1);
209  e.setInitialMixMoles(loglevel-1);
210 
211  // store component indices
212  if (m_nComponents > m_kk) {
214  }
215  for (size_t m = 0; m < m_nComponents; m++) {
216  m_component[m] = e.componentIndex(m);
217  }
218  for (size_t k = 0; k < m_kk; k++) {
219  if (s.moleFraction(k) > 0.0) {
220  if (loglevel > 0)
222  s.moleFraction(k));
223  }
224  }
225  /*
226  * Update the current values of the temp, density, and
227  * mole fraction, and element abundance vectors kept
228  * within the ChemEquil object.
229  */
230  update(s);
231 
232 #ifdef DEBUG_MODE
233  if (ChemEquil_print_lvl > 0) {
234  writelog("setInitialMoles: Estimated Mole Fractions\n");
235  writelogf(" Temperature = %g\n", s.temperature());
236  writelogf(" Pressure = %g\n", s.pressure());
237  for (size_t k = 0; k < m_kk; k++) {
238  string nnn = s.speciesName(k);
239  double mf = s.moleFraction(k);
240  writelogf(" %-12s % -10.5g\n", nnn.c_str(), mf);
241  }
242  writelog(" Element_Name ElementGoal ElementMF\n");
243  for (size_t m = 0; m < m_mm; m++) {
244  string nnn = s.elementName(m);
245  writelogf(" %-12s % -10.5g% -10.5g\n",
246  nnn.c_str(), elMoleGoal[m], m_elementmolefracs[m]);
247  }
248  }
249 #endif
250 
251  iok = 0;
252  } catch (CanteraError& err) {
253  err.save();
254  iok = -1;
255  }
256  if (loglevel > 0) {
257  endLogGroup();
258  }
259  return iok;
260 }
261 
263  vector_fp& elMolesGoal, int loglevel)
264 {
265  if (loglevel > 0) {
266  beginLogGroup("estimateElementPotentials");
267  }
268  //for (k = 0; k < m_kk; k++) {
269  // if (m_molefractions[k] > 0.0) {
270  // m_molefractions[k] = fmaxx(m_molefractions[k], 0.05);
271  // }
272  //}
273  //s.setState_PX(s.pressure(), m_molefractions.begin());
274 
275 
276  vector_fp b(m_mm, -999.0);
277  vector_fp mu_RT(m_kk, 0.0);
278  vector_fp xMF_est(m_kk, 0.0);
279 
280  s.getMoleFractions(DATA_PTR(xMF_est));
281  for (size_t n = 0; n < s.nSpecies(); n++) {
282  if (xMF_est[n] < 1.0E-20) {
283  xMF_est[n] = 1.0E-20;
284  }
285  }
286  s.setMoleFractions(DATA_PTR(xMF_est));
287  s.getMoleFractions(DATA_PTR(xMF_est));
288 
289  MultiPhase mp;
290  mp.addPhase(&s, 1.0);
291  mp.init();
292  int usedZeroedSpecies = 0;
293  vector_fp formRxnMatrix;
294  m_nComponents = BasisOptimize(&usedZeroedSpecies, false,
295  &mp, m_orderVectorSpecies,
296  m_orderVectorElements, formRxnMatrix);
297 
298  for (size_t m = 0; m < m_nComponents; m++) {
299  size_t k = m_orderVectorSpecies[m];
300  m_component[m] = k;
301  if (xMF_est[k] < 1.0E-8) {
302  xMF_est[k] = 1.0E-8;
303  }
304  }
305  s.setMoleFractions(DATA_PTR(xMF_est));
306  s.getMoleFractions(DATA_PTR(xMF_est));
307 
308  size_t nct = Cantera::ElemRearrange(m_nComponents, elMolesGoal, &mp,
309  m_orderVectorSpecies, m_orderVectorElements);
310  if (nct != m_nComponents) {
311  throw CanteraError("ChemEquil::estimateElementPotentials",
312  "confused");
313  }
314 
315  s.getChemPotentials(DATA_PTR(mu_RT));
316  doublereal rrt = 1.0/(GasConstant* s.temperature());
317  scale(mu_RT.begin(), mu_RT.end(), mu_RT.begin(), rrt);
318 
319 #ifdef DEBUG_MODE
320  if (ChemEquil_print_lvl > 0) {
321  for (size_t m = 0; m < m_nComponents; m++) {
322  int isp = m_component[m];
323  string nnn = s.speciesName(isp);
324  writelogf("isp = %d, %s\n", isp, nnn.c_str());
325  }
326  double pres = s.pressure();
327  double temp = s.temperature();
328  writelogf("Pressure = %g\n", pres);
329  writelogf("Temperature = %g\n", temp);
330  writelog(" id Name MF mu/RT \n");
331  for (size_t n = 0; n < s.nSpecies(); n++) {
332  string nnn = s.speciesName(n);
333  writelogf("%10d %15s %10.5g %10.5g\n",
334  n, nnn.c_str(), xMF_est[n], mu_RT[n]);
335  }
336  }
337 #endif
338  DenseMatrix aa(m_nComponents, m_nComponents, 0.0);
339  for (size_t m = 0; m < m_nComponents; m++) {
340  for (size_t n = 0; n < m_nComponents; n++) {
341  aa(m,n) = nAtoms(m_component[m], m_orderVectorElements[n]);
342  }
343  b[m] = mu_RT[m_component[m]];
344  }
345 
346  int info = solve(aa, DATA_PTR(b));
347  if (info) {
348  if (loglevel > 0) {
349  addLogEntry("failed to estimate initial element potentials.");
350  }
351  info = -2;
352  }
353  for (size_t m = 0; m < m_nComponents; m++) {
354  lambda_RT[m_orderVectorElements[m]] = b[m];
355  }
356  for (size_t m = m_nComponents; m < m_mm; m++) {
357  lambda_RT[m_orderVectorElements[m]] = 0.0;
358  }
359  if (info == 0) {
360  if (loglevel > 0) {
361  for (size_t m = 0; m < m_mm; m++) {
362  addLogEntry(s.elementName(m),lambda_RT[m]);
363  }
364  }
365  }
366 
367 #ifdef DEBUG_MODE
368  if (ChemEquil_print_lvl > 0) {
369  writelog(" id CompSpecies ChemPot EstChemPot Diff\n");
370  for (size_t m = 0; m < m_nComponents; m++) {
371  int isp = m_component[m];
372  double tmp = 0.0;
373  string sname = s.speciesName(isp);
374  for (size_t n = 0; n < m_mm; n++) {
375  tmp += nAtoms(isp, n) * lambda_RT[n];
376  }
377  writelogf("%3d %16s %10.5g %10.5g %10.5g\n",
378  m, sname.c_str(), mu_RT[isp], tmp, tmp - mu_RT[isp]);
379  }
380 
381  writelog(" id ElName Lambda_RT\n");
382  for (size_t m = 0; m < m_mm; m++) {
383  string ename = s.elementName(m);
384  writelogf(" %3d %6s %10.5g\n", m, ename.c_str(), lambda_RT[m]);
385  }
386  }
387 #endif
388  if (loglevel > 0) {
389  endLogGroup();
390  }
391  return info;
392 }
393 
394 int ChemEquil::equilibrate(thermo_t& s, const char* XY,
395  bool useThermoPhaseElementPotentials, int loglevel)
396 {
397  vector_fp elMolesGoal(s.nElements());
398  initialize(s);
399  update(s);
400  copy(m_elementmolefracs.begin(), m_elementmolefracs.end(),
401  elMolesGoal.begin());
402  return equilibrate(s, XY, elMolesGoal, useThermoPhaseElementPotentials,
403  loglevel-1);
404 }
405 
406 int ChemEquil::equilibrate(thermo_t& s, const char* XYstr,
407  vector_fp& elMolesGoal,
408  bool useThermoPhaseElementPotentials,
409  int loglevel)
410 {
411  doublereal xval, yval, tmp;
412  int fail = 0;
413 
414  bool tempFixed = true;
415  int XY = _equilflag(XYstr);
416 
417  vector_fp state;
418  s.saveState(state);
419 
420  /*
421  * Check Compatibility
422  */
423  if (m_mm != s.nElements() || m_kk != s.nSpecies()) {
424  throw CanteraError("ChemEquil::equilibrate ERROR",
425  "Input ThermoPhase is incompatible with initialization");
426  }
427 
428 #ifdef DEBUG_MODE
429  int n;
430  const vector<string>& eNames = s.elementNames();
431 #endif
432  if (loglevel > 0) {
433  beginLogGroup("ChemEquil::equilibrate");
434  }
435  initialize(s);
436  update(s);
437  switch (XY) {
438  case TP:
439  case PT:
440  m_p1.reset(new TemperatureCalculator<thermo_t>);
441  m_p2.reset(new PressureCalculator<thermo_t>);
442  break;
443  case HP:
444  case PH:
445  tempFixed = false;
446  m_p1.reset(new EnthalpyCalculator<thermo_t>);
447  m_p2.reset(new PressureCalculator<thermo_t>);
448  break;
449  case SP:
450  case PS:
451  tempFixed = false;
452  m_p1.reset(new EntropyCalculator<thermo_t>);
453  m_p2.reset(new PressureCalculator<thermo_t>);
454  break;
455  case SV:
456  case VS:
457  tempFixed = false;
458  m_p1.reset(new EntropyCalculator<thermo_t>);
459  m_p2.reset(new DensityCalculator<thermo_t>);
460  break;
461  case TV:
462  case VT:
463  m_p1.reset(new TemperatureCalculator<thermo_t>);
464  m_p2.reset(new DensityCalculator<thermo_t>);
465  break;
466  case UV:
467  case VU:
468  tempFixed = false;
469  m_p1.reset(new IntEnergyCalculator<thermo_t>);
470  m_p2.reset(new DensityCalculator<thermo_t>);
471  break;
472  default:
473  if (loglevel > 0) {
474  endLogGroup("ChemEquil::equilibrate");
475  }
476  throw CanteraError("equilibrate","illegal property pair.");
477  }
478  if (loglevel > 0) {
479  addLogEntry("Problem type","fixed "+m_p1->symbol()+", "+m_p2->symbol());
480  addLogEntry(m_p1->symbol(), m_p1->value(s));
481  addLogEntry(m_p2->symbol(), m_p2->value(s));
482  }
483  // If the temperature is one of the specified variables, and
484  // it is outside the valid range, throw an exception.
485  if (tempFixed) {
486  double tfixed = s.temperature();
487  if (tfixed > s.maxTemp() + 1.0 || tfixed < s.minTemp() - 1.0) {
488  if (loglevel > 0) {
489  endLogGroup("ChemEquil::equilibrate");
490  }
491  throw CanteraError("ChemEquil","Specified temperature ("
492  +fp2str(s.temperature())+" K) outside "
493  "valid range of "+fp2str(s.minTemp())+" K to "
494  +fp2str(s.maxTemp())+" K\n");
495  }
496  }
497 
498  /*
499  * Before we do anything to change the ThermoPhase object,
500  * we calculate and store the two specified thermodynamic
501  * properties that we are after.
502  */
503  xval = m_p1->value(s);
504  yval = m_p2->value(s);
505 
506  size_t mm = m_mm;
507  size_t nvar = mm + 1;
508  DenseMatrix jac(nvar, nvar); // jacobian
509  vector_fp x(nvar, -102.0); // solution vector
510  vector_fp res_trial(nvar, 0.0); // residual
511 
512  /*
513  * Replace one of the element abundance fraction equations
514  * with the specified property calculation.
515  *
516  * We choose the equation of the element with the highest element
517  * abundance.
518  */
519  size_t m;
520  tmp = -1.0;
521  for (size_t im = 0; im < m_nComponents; im++) {
522  m = m_orderVectorElements[im];
523  if (elMolesGoal[m] > tmp) {
524  m_skip = m;
525  tmp = elMolesGoal[m];
526  }
527  }
528  if (tmp <= 0.0) {
529  throw CanteraError("ChemEquil",
530  "Element Abundance Vector is zeroed");
531  }
532 
533  // start with a composition with everything non-zero. Note
534  // that since we have already save the target element moles,
535  // changing the composition at this point only affects the
536  // starting point, not the final solution.
537  vector_fp xmm(m_kk, 0.0);
538  for (size_t k = 0; k < m_kk; k++) {
539  xmm[k] = s.moleFraction(k) + 1.0E-32;
540  }
541  s.setMoleFractions(DATA_PTR(xmm));
542 
543  /*
544  * Update the internally stored values of m_temp,
545  * m_dens, and the element mole fractions.
546  */
547  update(s);
548 
549  doublereal tmaxPhase = s.maxTemp();
550  doublereal tminPhase = s.minTemp();
551  // loop to estimate T
552  if (!tempFixed) {
553  if (loglevel > 0) {
554  beginLogGroup("Initial T Estimate");
555  }
556 
557  doublereal tmin;
558  doublereal tmax;
559 
560  tmin = s.temperature();
561  if (tmin < tminPhase) {
562  tmin = tminPhase;
563  }
564  if (tmin > tmaxPhase) {
565  tmin = tmaxPhase - 20;
566  }
567  tmax = tmin + 10.;
568  if (tmax > tmaxPhase) {
569  tmax = tmaxPhase;
570  }
571  if (tmax < tminPhase) {
572  tmax = tminPhase + 20;
573  }
574 
575  doublereal slope, phigh, plow, pval, dt;
576 
577  // first get the property values at the upper and lower
578  // temperature limits. Since p1 (h, s, or u) is monotonic
579  // in T, these values determine the upper and lower
580  // bounnds (phigh, plow) for p1.
581 
582  s.setTemperature(tmax);
583  setInitialMoles(s, elMolesGoal, loglevel - 1);
584  phigh = m_p1->value(s);
585 
586  s.setTemperature(tmin);
587  setInitialMoles(s, elMolesGoal, loglevel - 1);
588  plow = m_p1->value(s);
589 
590  // start with T at the midpoint of the range
591  doublereal t0 = 0.5*(tmin + tmax);
592  s.setTemperature(t0);
593 
594  // loop up to 5 times
595  for (int it = 0; it < 10; it++) {
596 
597  // set the composition and get p1
598  setInitialMoles(s, elMolesGoal, loglevel - 1);
599  pval = m_p1->value(s);
600 
601  // If this value of p1 is greater than the specified
602  // property value, then the current temperature is too
603  // high. Use it as the new upper bound. Otherwise, it
604  // is too low, so use it as the new lower bound.
605  if (pval > xval) {
606  tmax = t0;
607  phigh = pval;
608  } else {
609  tmin = t0;
610  plow = pval;
611  }
612 
613  // Determine the new T estimate by linearly interpolating
614  // between the upper and lower bounds
615  slope = (phigh - plow)/(tmax - tmin);
616  dt = (xval - pval)/slope;
617 
618  // If within 50 K, terminate the search
619  if (fabs(dt) < 50.0) {
620  break;
621  }
622  if (dt > 200.) {
623  dt = 200.;
624  }
625  if (dt < -200.) {
626  dt = -200.;
627  }
628  if ((t0 + dt) < tminPhase) {
629  dt = 0.5*((t0) + tminPhase) - t0;
630  }
631  if ((t0 + dt) > tmaxPhase) {
632  dt = 0.5*((t0) + tmaxPhase) - t0;
633  }
634  // update the T estimate
635  t0 = t0 + dt;
636  if (t0 <= tminPhase || t0 >= tmaxPhase) {
637  printf("We shouldn't be here\n");
638  exit(EXIT_FAILURE);
639  }
640  if (loglevel > 0) {
641  addLogEntry("new T estimate", t0);
642  }
643  if (t0 < 100.) {
644  printf("t0 - we are here %g\n", t0);
645  exit(EXIT_FAILURE);
646  }
647  s.setTemperature(t0);
648  }
649  if (loglevel > 0) {
650  endLogGroup("Initial T Estimate"); // initial T estimate
651  }
652  }
653 
654 
655  setInitialMoles(s, elMolesGoal,loglevel);
656 
657  /*
658  * If requested, get the initial estimate for the
659  * chemical potentials from the ThermoPhase object
660  * itself. Or else, create our own estimate.
661  */
662  if (useThermoPhaseElementPotentials) {
663  bool haveEm = s.getElementPotentials(DATA_PTR(x));
664  if (haveEm) {
665  doublereal rt = GasConstant * s.temperature();
666  if (s.temperature() < 100.) {
667  printf("we are here %g\n", s.temperature());
668  }
669  for (m = 0; m < m_mm; m++) {
670  x[m] /= rt;
671  }
672  } else {
673  estimateElementPotentials(s, x, elMolesGoal);
674  }
675  } else {
676  /*
677  * Calculate initial estimates of the element potentials.
678  * This algorithm uese the MultiPhaseEquil object's
679  * initialization capabilities to calculate an initial
680  * estimate of the mole fractions for a set of linearly
681  * independent component species. Then, the element
682  * potentials are solved for based on the chemical
683  * potentials of the component species.
684  */
685  estimateElementPotentials(s, x, elMolesGoal);
686  }
687 
688 
689  /*
690  * Do a better estimate of the element potentials.
691  * We have found that the current estimate may not be good
692  * enough to avoid drastic numerical issues associated with
693  * the use of a numerically generated jacobian.
694  *
695  * The Brinkley algorithm assumes a constant T, P system
696  * and uses a linearized analytical Jacobian that turns out
697  * to be very stable.
698  */
699  int info = estimateEP_Brinkley(s, x, elMolesGoal);
700  if (info != 0) {
701  if (info == 1) {
702  addLogEntry("estimateEP_Brinkley didn't converge in given max iterations");
703  } else if (info == -3) {
704  addLogEntry("estimateEP_Brinkley had a singular Jacobian. Continuing anyway");
705  }
706  } else {
707  setToEquilState(s, x, s.temperature());
708  // Tempting -> However, nonideal is a problem. Turn on if not worried
709  // about nonideality and you are having problems with the main
710  // algorithm.
711  //if (XY == TP) {
712  // endLogGroup("ChemEquil::equilibrate");
713  // return 0;
714  //}
715  }
716 
717  /*
718  * Install the log(temp) into the last solution unknown
719  * slot.
720  */
721  x[m_mm] = log(s.temperature());
722 
723  /*
724  * Setting the max and min values for x[]. Also, if element
725  * abundance vector is zero, setting x[] to -1000. This
726  * effectively zeroes out all species containing that element.
727  */
728  vector_fp above(nvar);
729  vector_fp below(nvar);
730  for (m = 0; m < mm; m++) {
731  above[m] = 200.0;
732  below[m] = -2000.0;
733  if (elMolesGoal[m] < m_elemFracCutoff && m != m_eloc) {
734  x[m] = -1000.0;
735  }
736  }
737  /*
738  * Set the temperature bounds to be 25 degrees different than the max and min
739  * temperatures.
740  */
741  above[mm] = log(s.maxTemp() + 25.0);
742  below[mm] = log(s.minTemp() - 25.0);
743 
744  vector_fp grad(nvar, 0.0); // gradient of f = F*F/2
745  vector_fp oldx(nvar, 0.0); // old solution
746  vector_fp oldresid(nvar, 0.0);
747  doublereal f, oldf;
748 
749  int iter = 0;
750  doublereal fctr = 1.0, newval;
751 
752  goto converge;
753 next:
754 
755  iter++;
756  if (iter > 1 && loglevel > 0) {
757  endLogGroup("Iteration "+int2str(iter-1)); // iteration
758  }
759  if (loglevel > 0) {
760  beginLogGroup("Iteration "+int2str(iter));
761  }
762 
763  // compute the residual and the jacobian using the current
764  // solution vector
765  equilResidual(s, x, elMolesGoal, res_trial, xval, yval);
766  f = 0.5*dot(res_trial.begin(), res_trial.end(), res_trial.begin());
767  addLogEntry("Residual norm", f);
768 
769  // Compute the Jacobian matrix
770  equilJacobian(s, x, elMolesGoal, jac, xval, yval);
771 
772 #ifdef DEBUG_MODE
773  if (ChemEquil_print_lvl > 0) {
774  writelogf("Jacobian matrix %d:\n", iter);
775  for (m = 0; m <= m_mm; m++) {
776  writelog(" [ ");
777  for (n = 0; n <= m_mm; n++) {
778  writelogf("%10.5g ", jac(m,n));
779  }
780  writelog(" ]");
781  char xName[32];
782  if (m < m_mm) {
783  string nnn = eNames[m];
784  sprintf(xName, "x_%-10s", nnn.c_str());
785  } else {
786  sprintf(xName, "x_XX");
787  }
788  if (m_eloc == m) {
789  sprintf(xName, "x_ELOC");
790  }
791  if (m == m_skip) {
792  sprintf(xName, "x_YY");
793  }
794  writelogf("%-12s", xName);
795  writelogf(" = - (%10.5g)\n", res_trial[m]);
796  }
797  }
798 #endif
799 
800  copy(x.begin(), x.end(), oldx.begin());
801  oldf = f;
802  scale(res_trial.begin(), res_trial.end(), res_trial.begin(), -1.0);
803 
804  /*
805  * Solve the system
806  */
807  try {
808  info = solve(jac, DATA_PTR(res_trial));
809  } catch (CanteraError& err) {
810  err.save();
811  addLogEntry("Jacobian is singular.");
812  endLogGroup(); // iteration
813  endLogGroup(); // equilibrate
814  s.restoreState(state);
815 
816  throw CanteraError("equilibrate",
817  "Jacobian is singular. \nTry adding more species, "
818  "changing the elemental composition slightly, \nor removing "
819  "unused elements.");
820  //return -3;
821  }
822 
823  // find the factor by which the Newton step can be multiplied
824  // to keep the solution within bounds.
825  fctr = 1.0;
826  for (m = 0; m < nvar; m++) {
827  newval = x[m] + res_trial[m];
828  if (newval > above[m]) {
829  fctr = std::max(0.0,
830  std::min(fctr,0.8*(above[m] - x[m])/(newval - x[m])));
831  } else if (newval < below[m]) {
832  if (m < m_mm && (m != m_skip)) {
833  res_trial[m] = -50;
834  if (x[m] < below[m] + 50.) {
835  res_trial[m] = below[m] - x[m];
836  }
837  } else {
838  fctr = std::min(fctr, 0.8*(x[m] - below[m])/(x[m] - newval));
839  }
840  }
841  // Delta Damping
842  if (m == mm) {
843  if (fabs(res_trial[mm]) > 0.2) {
844  fctr = std::min(fctr, 0.2/fabs(res_trial[mm]));
845  }
846  }
847  }
848  if (fctr != 1.0) {
849  addLogEntry("WARNING: factor to keep solution in bounds", fctr);
850 #ifdef DEBUG_MODE
851  if (ChemEquil_print_lvl > 0) {
852  writelogf("WARNING Soln Damping because of bounds: %g\n", fctr);
853  }
854 #endif
855  }
856 
857  // multiply the step by the scaling factor
858  scale(res_trial.begin(), res_trial.end(), res_trial.begin(), fctr);
859 
860  if (!dampStep(s, oldx, oldf, grad, res_trial,
861  x, f, elMolesGoal , xval, yval)) {
862  fail++;
863  if (fail > 3) {
864  addLogEntry("dampStep","Failed 3 times. Giving up.");
865  endLogGroup(); // iteration
866  endLogGroup(); // equilibrate
867  s.restoreState(state);
868  throw CanteraError("equilibrate",
869  "Cannot find an acceptable Newton damping coefficient.");
870  //return -4;
871  }
872  } else {
873  fail = 0;
874  }
875 
876 converge:
877 
878  // check for convergence.
879  equilResidual(s, x, elMolesGoal, res_trial, xval, yval);
880  f = 0.5*dot(res_trial.begin(), res_trial.end(), res_trial.begin());
881  doublereal xx, yy, deltax, deltay;
882  xx = m_p1->value(s);
883  yy = m_p2->value(s);
884  deltax = (xx - xval)/xval;
885  deltay = (yy - yval)/yval;
886  doublereal rmax = 0.0;
887  bool passThis = true;
888  for (m = 0; m < nvar; m++) {
889  double tval = options.relTolerance;
890  if (m < mm) {
891  /*
892  * Special case convergence requirements for electron element.
893  * This is a special case because the element coefficients may
894  * be both positive and negative. And, typically they sum to 0.0.
895  * Therefore, there is no natural absolute value for this quantity.
896  * We supply the absolute value tolerance here. Note, this is
897  * made easier since the element abundances are normalized to one
898  * within this routine.
899  *
900  * Note, the 1.0E-13 value was recently relaxed from 1.0E-15, because
901  * convergence failures were found to occur for the lower value
902  * at small pressure (0.01 pascal).
903  */
904  if (m == m_eloc) {
905  tval = elMolesGoal[m] * options.relTolerance + options.absElemTol
906  + 1.0E-13;
907  } else {
908  tval = elMolesGoal[m] * options.relTolerance + options.absElemTol;
909  }
910  }
911  if (fabs(res_trial[m]) > tval) {
912  passThis = false;
913  }
914  }
915  if (iter > 0 && passThis
916  && fabs(deltax) < options.relTolerance
917  && fabs(deltay) < options.relTolerance) {
918  options.iterations = iter;
919  if (loglevel > 0) {
920  endLogGroup("Iteration "+int2str(iter)); // iteration
921  beginLogGroup("Converged solution");
922  addLogEntry("Iterations",iter);
923  addLogEntry("Relative error in "+m_p1->symbol(),deltax);
924  addLogEntry("Relative error in "+m_p2->symbol(),deltay);
925  addLogEntry("Max residual",rmax);
926  beginLogGroup("Element potentials");
927  }
928  doublereal rt = GasConstant* s.temperature();
929  for (m = 0; m < m_mm; m++) {
930  m_lambda[m] = x[m]*rt;
931  if (loglevel > 0) {
932  addLogEntry("element "+ s.elementName(m), fp2str(x[m]));
933  }
934  }
935 
936  if (m_eloc != npos) {
937  adjustEloc(s, elMolesGoal);
938  }
939  /*
940  * Save the calculated and converged element potentials
941  * to the original ThermoPhase object.
942  */
944  if (loglevel > 0) {
945  addLogEntry("Saving Element Potentials to ThermoPhase Object");
946  endLogGroup("Element potentials");
947  }
948  if (s.temperature() > s.maxTemp() + 1.0 ||
949  s.temperature() < s.minTemp() - 1.0) {
950  writelog("Warning: Temperature ("
951  +fp2str(s.temperature())+" K) outside "
952  "valid range of "+fp2str(s.minTemp())+" K to "
953  +fp2str(s.maxTemp())+" K\n");
954  }
955  if (loglevel > 0) {
956  endLogGroup("Converged solution");
957  endLogGroup("ChemEquil::equilibrate");
958  }
959  return 0;
960  }
961 
962  // no convergence
963 
964  if (iter > options.maxIterations) {
965  if (loglevel > 0) {
966  addLogEntry("equilibrate","no convergence");
967  endLogGroup("Iteration "+int2str(iter));
968  endLogGroup("ChemEquil::equilibrate");
969  }
970  s.restoreState(state);
971  throw CanteraError("equilibrate",
972  "no convergence in "+int2str(options.maxIterations)
973  +" iterations.");
974  //return -1;
975  }
976  goto next;
977 }
978 
980  double oldf, vector_fp& grad, vector_fp& step, vector_fp& x,
981  double& f, vector_fp& elmols, double xval, double yval)
982 {
983  double damp;
984 
985  /*
986  * Carry out a delta damping approach on the dimensionless element potentials.
987  */
988  damp = 1.0;
989  for (size_t m = 0; m < m_mm; m++) {
990  if (m == m_eloc) {
991  if (step[m] > 1.25) {
992  damp = std::min(damp, 1.25 /step[m]);
993  }
994  if (step[m] < -1.25) {
995  damp = std::min(damp, -1.25 / step[m]);
996  }
997  } else {
998  if (step[m] > 0.75) {
999  damp = std::min(damp, 0.75 /step[m]);
1000  }
1001  if (step[m] < -0.75) {
1002  damp = std::min(damp, -0.75 / step[m]);
1003  }
1004  }
1005  }
1006 
1007  /*
1008  * Update the solution unknown
1009  */
1010  for (size_t m = 0; m < x.size(); m++) {
1011  x[m] = oldx[m] + damp * step[m];
1012  }
1013 #ifdef DEBUG_MODE
1014  if (ChemEquil_print_lvl > 0) {
1015  writelogf("Solution Unknowns: damp = %g\n", damp);
1016  writelog(" X_new X_old Step\n");
1017  for (size_t m = 0; m < m_mm; m++) {
1018  writelogf(" % -10.5g % -10.5g % -10.5g\n", x[m], oldx[m], step[m]);
1019  }
1020  }
1021 #endif
1022  return 1;
1023 }
1024 
1026  const vector_fp& elmFracGoal, vector_fp& resid,
1027  doublereal xval, doublereal yval, int loglevel)
1028 {
1029  if (loglevel > 0) {
1030  beginLogGroup("ChemEquil::equilResidual");
1031  }
1032  doublereal xx, yy;
1033  doublereal temp = exp(x[m_mm]);
1034  setToEquilState(s, x, temp);
1035 
1036  // residuals are the total element moles
1037  vector_fp& elmFrac = m_elementmolefracs;
1038  for (size_t n = 0; n < m_mm; n++) {
1039  size_t m = m_orderVectorElements[n];
1040  // drive element potential for absent elements to -1000
1041  if (elmFracGoal[m] < m_elemFracCutoff && m != m_eloc) {
1042  resid[m] = x[m] + 1000.0;
1043  } else if (n >= m_nComponents) {
1044  resid[m] = x[m];
1045  } else {
1046  /*
1047  * Change the calculation for small element number, using
1048  * L'Hopital's rule.
1049  * The log formulation is unstable.
1050  */
1051  if (elmFracGoal[m] < 1.0E-10 || elmFrac[m] < 1.0E-10 || m == m_eloc) {
1052  resid[m] = elmFracGoal[m] - elmFrac[m];
1053  } else {
1054  resid[m] = log((1.0 + elmFracGoal[m]) / (1.0 + elmFrac[m]));
1055  }
1056  }
1057  if (loglevel > 0)
1058  addLogEntry(s.elementName(m),fp2str(elmFrac[m])+" ("
1059  +fp2str(elmFracGoal[m])+")");
1060  }
1061 
1062 #ifdef DEBUG_MODE
1063  if (ChemEquil_print_lvl > 0 && !m_doResPerturb) {
1064  writelog("Residual: ElFracGoal ElFracCurrent Resid\n");
1065  for (int n = 0; n < m_mm; n++) {
1066  writelogf(" % -14.7E % -14.7E % -10.5E\n",
1067  elmFracGoal[n], elmFrac[n], resid[n]);
1068  }
1069  }
1070 #endif
1071 
1072  xx = m_p1->value(s);
1073  yy = m_p2->value(s);
1074  resid[m_mm] = xx/xval - 1.0;
1075  resid[m_skip] = yy/yval - 1.0;
1076  if (loglevel > 0) {
1077  string xstr = fp2str(xx)+" ("+fp2str(xval)+")";
1078  addLogEntry(m_p1->symbol(), xstr);
1079  string ystr = fp2str(yy)+" ("+fp2str(yval)+")";
1080  addLogEntry(m_p2->symbol(), ystr);
1081  endLogGroup("ChemEquil::equilResidual");
1082  }
1083 
1084 #ifdef DEBUG_MODE
1085  if (ChemEquil_print_lvl > 0 && !m_doResPerturb) {
1086  writelog(" Goal Xvalue Resid\n");
1087  writelogf(" XX : % -14.7E % -14.7E % -10.5E\n", xval, xx, resid[m_mm]);
1088  writelogf(" YY(%1d): % -14.7E % -14.7E % -10.5E\n", m_skip, yval, yy, resid[m_skip]);
1089  }
1090 #endif
1091 }
1092 
1093 void ChemEquil::equilJacobian(thermo_t& s, vector_fp& x,
1094  const vector_fp& elmols, DenseMatrix& jac,
1095  doublereal xval, doublereal yval, int loglevel)
1096 {
1097  if (loglevel > 0) {
1098  beginLogGroup("equilJacobian");
1099  }
1100  vector_fp& r0 = m_jwork1;
1101  vector_fp& r1 = m_jwork2;
1102  size_t len = x.size();
1103  r0.resize(len);
1104  r1.resize(len);
1105  size_t n, m;
1106  doublereal rdx, dx, xsave, dx2;
1107  doublereal atol = 1.e-10;
1108 
1109  equilResidual(s, x, elmols, r0, xval, yval, loglevel-1);
1110 
1111  m_doResPerturb = false;
1112  for (n = 0; n < len; n++) {
1113  xsave = x[n];
1114  dx = atol;
1115  dx2 = fabs(xsave) * 1.0E-7;
1116  if (dx2 > dx) {
1117  dx = dx2;
1118  }
1119  x[n] = xsave + dx;
1120  dx = x[n] - xsave;
1121  rdx = 1.0/dx;
1122 
1123  // calculate perturbed residual
1124 
1125  equilResidual(s, x, elmols, r1, xval, yval, loglevel-1);
1126 
1127  // compute nth column of Jacobian
1128 
1129  for (m = 0; m < x.size(); m++) {
1130  jac(m, n) = (r1[m] - r0[m])*rdx;
1131  }
1132  x[n] = xsave;
1133  }
1134  m_doResPerturb = false;
1135  if (loglevel > 0) {
1136  endLogGroup("equilJacobian");
1137  }
1138 }
1139 
1140 double ChemEquil::calcEmoles(thermo_t& s, vector_fp& x, const double& n_t,
1141  const vector_fp& Xmol_i_calc,
1142  vector_fp& eMolesCalc, vector_fp& n_i_calc,
1143  double pressureConst)
1144 {
1145  double n_t_calc = 0.0;
1146  double tmp;
1147  /*
1148  * Calculate the activity coefficients of the solution, at the
1149  * previous solution state.
1150  */
1151  vector_fp actCoeff(m_kk, 1.0);
1152  s.setMoleFractions(DATA_PTR(Xmol_i_calc));
1153  s.setPressure(pressureConst);
1154  s.getActivityCoefficients(DATA_PTR(actCoeff));
1155 
1156  for (size_t k = 0; k < m_kk; k++) {
1157  tmp = - (m_muSS_RT[k] + log(actCoeff[k]));
1158  for (size_t m = 0; m < m_mm; m++) {
1159  tmp += nAtoms(k,m) * x[m];
1160  }
1161  if (tmp > 100.) {
1162  tmp = 100.;
1163  }
1164  if (tmp < -300.) {
1165  n_i_calc[k] = 0.0;
1166  } else {
1167  n_i_calc[k] = n_t * exp(tmp);
1168  }
1169  n_t_calc += n_i_calc[k];
1170  }
1171  for (size_t m = 0; m < m_mm; m++) {
1172  eMolesCalc[m] = 0.0;
1173  for (size_t k = 0; k < m_kk; k++) {
1174  eMolesCalc[m] += nAtoms(k,m) * n_i_calc[k];
1175  }
1176  }
1177  return n_t_calc;
1178 }
1179 
1181  vector_fp& elMoles)
1182 {
1183  /*
1184  * Before we do anything, we will save the state of the solution.
1185  * Then, if things go drastically wrong, we will restore the
1186  * saved state.
1187  */
1188  vector_fp state;
1189  s.saveState(state);
1190  double tmp, sum;
1191  bool modifiedMatrix = false;
1192  size_t neq = m_mm+1;
1193  int retn = 1;
1194  size_t m, n, k, im;
1195  DenseMatrix a1(neq, neq, 0.0);
1196  vector_fp b(neq, 0.0);
1197  vector_fp n_i(m_kk,0.0);
1198  vector_fp n_i_calc(m_kk,0.0);
1199  vector_fp actCoeff(m_kk, 1.0);
1200 
1201  vector_fp Xmol_i_calc(m_kk,0.0);
1202  double beta = 1.0;
1203 
1204  s.getMoleFractions(DATA_PTR(n_i));
1205  double pressureConst = s.pressure();
1206  copy(n_i.begin(), n_i.end(), Xmol_i_calc.begin());
1207 
1208  vector_fp x_old(m_mm+1, 0.0);
1209  vector_fp resid(m_mm+1, 0.0);
1210  vector_int lumpSum(m_mm+1, 0);
1211 
1212  /*
1213  * Get the nondimensional Gibbs functions for the species
1214  * at their standard states of solution at the current T and P
1215  * of the solution.
1216  */
1218 
1219 
1220  vector_fp eMolesCalc(m_mm, 0.0);
1221  vector_fp eMolesFix(m_mm, 0.0);
1222  double elMolesTotal = 0.0;
1223  for (m = 0; m < m_mm; m++) {
1224  elMolesTotal += elMoles[m];
1225  for (k = 0; k < m_kk; k++) {
1226  eMolesFix[m] += nAtoms(k,m) * n_i[k];
1227  }
1228  }
1229 
1230  for (m = 0; m < m_mm; m++) {
1231  if (x[m] > 50.0) {
1232  x[m] = 50.;
1233  }
1234  if (elMoles[m] > 1.0E-70) {
1235  if (x[m] < -100) {
1236  x[m] = -100.;
1237  }
1238  } else {
1239  if (x[m] < -1000.) {
1240  x[m] = -1000.;
1241  }
1242  }
1243  }
1244 
1245 
1246  double n_t = 0.0;
1247  double sum2 = 0.0;
1248  double nAtomsMax = 1.0;
1249  s.setMoleFractions(DATA_PTR(Xmol_i_calc));
1250  s.setPressure(pressureConst);
1251  s.getActivityCoefficients(DATA_PTR(actCoeff));
1252  for (k = 0; k < m_kk; k++) {
1253  tmp = - (m_muSS_RT[k] + log(actCoeff[k]));
1254  sum2 = 0.0;
1255  for (m = 0; m < m_mm; m++) {
1256  sum = nAtoms(k,m);
1257  tmp += sum * x[m];
1258  sum2 += sum;
1259  if (sum2 > nAtomsMax) {
1260  nAtomsMax = sum2;
1261  }
1262  }
1263  if (tmp > 100.) {
1264  n_t += 2.8E43;
1265  } else {
1266  n_t += exp(tmp);
1267  }
1268  }
1269 
1270 
1271 #ifdef DEBUG_MODE
1272  const vector<string>& eNames = s.elementNames();
1273  if (ChemEquil_print_lvl > 0) {
1274  writelog("estimateEP_Brinkley::\n\n");
1275  double temp = s.temperature();
1276  double pres = s.pressure();
1277  writelogf("temp = %g\n", temp);
1278  writelogf("pres = %g\n", pres);
1279  writelog("Initial mole numbers and mu_SS:\n");
1280  writelog(" Name MoleNum mu_SS actCoeff\n");
1281  for (k = 0; k < m_kk; k++) {
1282  string nnn = s.speciesName(k);
1283  writelogf("%15s %13.5g %13.5g %13.5g\n",
1284  nnn.c_str(), n_i[k], m_muSS_RT[k], actCoeff[k]);
1285  }
1286  writelogf("Initial n_t = %10.5g\n", n_t);
1287  writelog("Comparison of Goal Element Abundance with Initial Guess:\n");
1288  writelog(" eName eCurrent eGoal\n");
1289  for (m = 0; m < m_mm; m++) {
1290  string nnn = s.elementName(m);
1291  writelogf("%5s %13.5g %13.5g\n",nnn.c_str(), eMolesFix[m], elMoles[m]);
1292  }
1293  }
1294 #endif
1295  for (m = 0; m < m_mm; m++) {
1296  if (m != m_eloc) {
1297  if (elMoles[m] <= options.absElemTol) {
1298  x[m] = -200.;
1299  }
1300  }
1301  }
1302 
1303  /*
1304  * -------------------------------------------------------------------
1305  * Main Loop.
1306  */
1307  for (int iter = 0; iter < 20* options.maxIterations; iter++) {
1308  /*
1309  * Save the old solution
1310  */
1311  for (m = 0; m < m_mm; m++) {
1312  x_old[m] = x[m];
1313  }
1314  x_old[m_mm] = n_t;
1315  /*
1316  * Calculate the mole numbers of species
1317  */
1318 #ifdef DEBUG_MODE
1319  if (ChemEquil_print_lvl > 0) {
1320  writelogf("START ITERATION %d:\n", iter);
1321  }
1322 #endif
1323  /*
1324  * Calculate the mole numbers of species and elements.
1325  */
1326  double n_t_calc = calcEmoles(s, x, n_t, Xmol_i_calc, eMolesCalc, n_i_calc,
1327  pressureConst);
1328 
1329  for (k = 0; k < m_kk; k++) {
1330  Xmol_i_calc[k] = n_i_calc[k]/n_t_calc;
1331  }
1332 
1333 
1334 #ifdef DEBUG_MODE
1335  if (ChemEquil_print_lvl > 0) {
1336  writelog(" Species: Calculated_Moles Calculated_Mole_Fraction\n");
1337  for (k = 0; k < m_kk; k++) {
1338  string nnn = s.speciesName(k);
1339  writelogf("%15s: %10.5g %10.5g\n", nnn.c_str(), n_i_calc[k], Xmol_i_calc[k]);
1340  }
1341  writelogf("%15s: %10.5g\n", "Total Molar Sum", n_t_calc);
1342  writelogf("(iter %d) element moles bal: Goal Calculated\n", iter);
1343  for (m = 0; m < m_mm; m++) {
1344  string nnn = eNames[m];
1345  writelogf(" %8s: %10.5g %10.5g \n", nnn.c_str(), elMoles[m], eMolesCalc[m]);
1346  }
1347  }
1348 #endif
1349 
1350  double nCutoff;
1351 
1352  bool normalStep = true;
1353  /*
1354  * Decide if we are to do a normal step or a modified step
1355  */
1356  size_t iM = npos;
1357  for (m = 0; m < m_mm; m++) {
1358  if (elMoles[m] > 0.001 * elMolesTotal) {
1359  if (eMolesCalc[m] > 1000. * elMoles[m]) {
1360  normalStep = false;
1361  iM = m;
1362  }
1363  if (1000 * eMolesCalc[m] < elMoles[m]) {
1364  normalStep = false;
1365  iM = m;
1366  }
1367  }
1368  }
1369  if (DEBUG_MODE_ENABLED && ChemEquil_print_lvl > 0) {
1370  if (!normalStep) {
1371  writelogf(" NOTE: iter(%d) Doing an abnormal step due to row %d\n", iter, iM);
1372  }
1373  }
1374  if (!normalStep) {
1375  beta = 1.0;
1376  resid[m_mm] = 0.0;
1377  for (im = 0; im < m_mm; im++) {
1378  m = m_orderVectorElements[im];
1379  resid[m] = 0.0;
1380  if (im < m_nComponents) {
1381  if (elMoles[m] > 0.001 * elMolesTotal) {
1382  if (eMolesCalc[m] > 1000. * elMoles[m]) {
1383  resid[m] = -0.5;
1384  resid[m_mm] -= 0.5;
1385  }
1386  if (1000 * eMolesCalc[m] < elMoles[m]) {
1387  resid[m] = 0.5;
1388  resid[m_mm] += 0.5;
1389  }
1390  }
1391  }
1392  }
1393  if (n_t < (elMolesTotal / nAtomsMax)) {
1394  if (resid[m_mm] < 0.0) {
1395  resid[m_mm] = 0.1;
1396  }
1397  } else if (n_t > elMolesTotal) {
1398  if (resid[m_mm] > 0.0) {
1399  resid[m_mm] = 0.0;
1400  }
1401  }
1402  goto updateSolnVector;
1403  }
1404 
1405 
1406  /*
1407  * Determine whether the matrix should be dumbed down because
1408  * the coefficient matrix of species (with significant concentrations)
1409  * is rank deficient.
1410  *
1411  * The basic idea is that at any time during the calculation only a
1412  * small subset of species with sufficient concentration matters.
1413  * If the rank of the element coefficient matrix for that subset of species
1414  * is less than the number of elements, then the matrix created by
1415  * the Brinkley method below may become singular.
1416  *
1417  * The logic below looks for obvious cases where the current element
1418  * coefficient matrix is rank deficient.
1419  *
1420  * The way around rank-deficiency is to lump-sum the corresponding row
1421  * of the matrix. Note, lump-summing seems to work very well in terms of
1422  * its stability properties, i.e., it heads in the right direction,
1423  * albeit with lousy convergence rates.
1424  *
1425  * NOTE: This probably should be extended to a full blown Gauss-Jordan
1426  * factorization scheme in the future. For Example
1427  * the scheme below would fail for the set: HCl NH4Cl, NH3.
1428  * Hopefully, it's caught by the equal rows logic below.
1429  */
1430  for (m = 0; m < m_mm; m++) {
1431  lumpSum[m] = 1;
1432  }
1433 
1434  nCutoff = 1.0E-9 * n_t_calc;
1435 #ifdef DEBUG_MODE
1436  writelog(" Lump Sum Elements Calculation: \n", ChemEquil_print_lvl);
1437 #endif
1438  for (m = 0; m < m_mm; m++) {
1439  size_t kMSp = npos;
1440  size_t kMSp2 = npos;
1441  int nSpeciesWithElem = 0;
1442  for (k = 0; k < m_kk; k++) {
1443  if (n_i_calc[k] > nCutoff) {
1444  if (fabs(nAtoms(k,m)) > 0.001) {
1445  nSpeciesWithElem++;
1446  if (kMSp != npos) {
1447  kMSp2 = k;
1448  double factor = fabs(nAtoms(kMSp,m) / nAtoms(kMSp2,m));
1449  for (n = 0; n < m_mm; n++) {
1450  if (fabs(factor * nAtoms(kMSp2,n) - nAtoms(kMSp,n)) > 1.0E-8) {
1451  lumpSum[m] = 0;
1452  break;
1453  }
1454  }
1455  } else {
1456  kMSp = k;
1457  }
1458  }
1459  }
1460  }
1461 #ifdef DEBUG_MODE
1462  if (ChemEquil_print_lvl > 0) {
1463  string nnn = eNames[m];
1464  writelogf(" %5s %3d : %5d %5d\n",nnn.c_str(), lumpSum[m], kMSp, kMSp2);
1465  }
1466 #endif
1467  }
1468 
1469  /*
1470  * Formulate the matrix.
1471  */
1472  for (im = 0; im < m_mm; im++) {
1473  m = m_orderVectorElements[im];
1474  if (im < m_nComponents) {
1475  for (n = 0; n < m_mm; n++) {
1476  a1(m,n) = 0.0;
1477  for (k = 0; k < m_kk; k++) {
1478  a1(m,n) += nAtoms(k,m) * nAtoms(k,n) * n_i_calc[k];
1479  }
1480  }
1481  a1(m,m_mm) = eMolesCalc[m];
1482  a1(m_mm, m) = eMolesCalc[m];
1483  } else {
1484  for (n = 0; n <= m_mm; n++) {
1485  a1(m,n) = 0.0;
1486  }
1487  a1(m,m) = 1.0;
1488  }
1489  }
1490  a1(m_mm, m_mm) = 0.0;
1491 
1492  /*
1493  * Formulate the residual, resid, and the estimate for the convergence criteria, sum
1494  */
1495  sum = 0.0;
1496  for (im = 0; im < m_mm; im++) {
1497  m = m_orderVectorElements[im];
1498  if (im < m_nComponents) {
1499  resid[m] = elMoles[m] - eMolesCalc[m];
1500  } else {
1501  resid[m] = 0.0;
1502  }
1503  /*
1504  * For equations with positive and negative coefficients, (electronic charge),
1505  * we must mitigate the convergence criteria by a condition limited by
1506  * finite precision of inverting a matrix.
1507  * Other equations with just positive coefficients aren't limited by this.
1508  */
1509  if (m == m_eloc) {
1510  tmp = resid[m] / (elMoles[m] + elMolesTotal*1.0E-6 + options.absElemTol);
1511  } else {
1512  tmp = resid[m] / (elMoles[m] + options.absElemTol);
1513  }
1514  sum += tmp * tmp;
1515  }
1516 
1517  for (m = 0; m < m_mm; m++) {
1518  if (a1(m,m) < 1.0E-50) {
1519 #ifdef DEBUG_MODE
1520  if (ChemEquil_print_lvl > 0) {
1521  writelogf(" NOTE: Diagonalizing the analytical Jac row %d\n", m);
1522  }
1523 #endif
1524  for (n = 0; n < m_mm; n++) {
1525  a1(m,n) = 0.0;
1526  }
1527  a1(m,m) = 1.0;
1528  if (resid[m] > 0.0) {
1529  resid[m] = 1.0;
1530  } else if (resid[m] < 0.0) {
1531  resid[m] = -1.0;
1532  } else {
1533  resid[m] = 0.0;
1534  }
1535  }
1536  }
1537 
1538 
1539  resid[m_mm] = n_t - n_t_calc;
1540 
1541 #ifdef DEBUG_MODE
1542  if (ChemEquil_print_lvl > 0) {
1543  writelog("Matrix:\n");
1544  for (m = 0; m <= m_mm; m++) {
1545  writelog(" [");
1546  for (n = 0; n <= m_mm; n++) {
1547  writelogf(" %10.5g", a1(m,n));
1548  }
1549  writelogf("] = %10.5g\n", resid[m]);
1550  }
1551  }
1552 #endif
1553 
1554  tmp = resid[m_mm] /(n_t + 1.0E-15);
1555  sum += tmp * tmp;
1556 #ifdef DEBUG_MODE
1557  if (ChemEquil_print_lvl > 0) {
1558  writelogf("(it %d) Convergence = %g\n", iter, sum);
1559  }
1560 #endif
1561  /*
1562  * Insist on 20x accuracy compared to the top routine.
1563  * There are instances, for ill-conditioned or
1564  * singular matrices where this is needed to move
1565  * the system to a point where the matrices aren't
1566  * singular.
1567  */
1568  if (sum < 0.05 * options.relTolerance) {
1569  retn = 0;
1570  goto exit;
1571  }
1572 
1573  /*
1574  * Row Sum scaling
1575  */
1576  for (m = 0; m <= m_mm; m++) {
1577  tmp = 0.0;
1578  for (n = 0; n <= m_mm; n++) {
1579  tmp += fabs(a1(m,n));
1580  }
1581  if (m < m_mm && tmp < 1.0E-30) {
1582 #ifdef DEBUG_MODE
1583  if (ChemEquil_print_lvl > 0) {
1584  writelogf(" NOTE: Diagonalizing row %d\n", m);
1585  }
1586 #endif
1587  for (n = 0; n <= m_mm; n++) {
1588  if (n != m) {
1589  a1(m,n) = 0.0;
1590  a1(n,m) = 0.0;
1591  }
1592  }
1593  }
1594  tmp = 1.0/tmp;
1595  for (n = 0; n <= m_mm; n++) {
1596  a1(m,n) *= tmp;
1597  }
1598  resid[m] *= tmp;
1599  }
1600 
1601 #ifdef DEBUG_MODE
1602  if (ChemEquil_print_lvl > 0) {
1603  writelog("Row Summed Matrix:\n");
1604  for (m = 0; m <= m_mm; m++) {
1605  writelog(" [");
1606  for (n = 0; n <= m_mm; n++) {
1607  writelogf(" %10.5g", a1(m,n));
1608  }
1609  writelogf("] = %10.5g\n", resid[m]);
1610  }
1611  }
1612 #endif
1613  /*
1614  * Next Step: We have row-summed the equations.
1615  * However, there are some degenerate cases where two
1616  * rows will be multiplies of each other in terms of
1617  * 0 < m, 0 < m part of the matrix. This occurs on a case
1618  * by case basis, and depends upon the current state of the
1619  * element potential values, which affect the concentrations
1620  * of species.
1621  * So, the way we have found to eliminate this problem is to
1622  * lump-sum one of the rows of the matrix, except for the
1623  * last column, and stick it all on the diagonal.
1624  * Then, we at least have a non-singular matrix, and the
1625  * modified equation moves the corresponding unknown in the
1626  * correct direction.
1627  * The previous row-sum operation has made the identification
1628  * of identical rows much simpler.
1629  *
1630  * Note at least 6E-4 is necessary for the comparison.
1631  * I'm guessing 1.0E-3. If two rows are anywhere close to being
1632  * equivalent, the algorithm can get stuck in an oscillatory mode.
1633  */
1634  modifiedMatrix = false;
1635  for (m = 0; m < m_mm; m++) {
1636  size_t sameAsRow = npos;
1637  for (size_t im = 0; im < m; im++) {
1638  bool theSame = true;
1639  for (n = 0; n < m_mm; n++) {
1640  if (fabs(a1(m,n) - a1(im,n)) > 1.0E-7) {
1641  theSame = false;
1642  break;
1643  }
1644  }
1645  if (theSame) {
1646  sameAsRow = im;
1647  }
1648  }
1649  if (sameAsRow != npos || lumpSum[m]) {
1650 #ifdef DEBUG_MODE
1651  if (ChemEquil_print_lvl > 0) {
1652  if (lumpSum[m]) {
1653  writelogf("Lump summing row %d, due to rank deficiency analysis\n", m);
1654  } else if (sameAsRow != npos) {
1655  writelogf("Identified that rows %d and %d are the same\n", m, sameAsRow);
1656  }
1657  }
1658 #endif
1659  modifiedMatrix = true;
1660  for (n = 0; n < m_mm; n++) {
1661  if (n != m) {
1662  a1(m,m) += fabs(a1(m,n));
1663  a1(m,n) = 0.0;
1664  }
1665  }
1666  }
1667  }
1668 
1669  if (DEBUG_MODE_ENABLED && ChemEquil_print_lvl > 0 && modifiedMatrix) {
1670  writelog("Row Summed, MODIFIED Matrix:\n");
1671  for (m = 0; m <= m_mm; m++) {
1672  writelog(" [");
1673  for (n = 0; n <= m_mm; n++) {
1674  writelogf(" %10.5g", a1(m,n));
1675  }
1676  writelogf("] = %10.5g\n", resid[m]);
1677  }
1678  }
1679 
1680  try {
1681  solve(a1, DATA_PTR(resid));
1682  } catch (CanteraError& err) {
1683  err.save();
1684  addLogEntry("estimateEP_Brinkley:Jacobian is singular.");
1685 #ifdef DEBUG_MODE
1686  writelog("Matrix is SINGULAR.ERROR\n", ChemEquil_print_lvl);
1687 #endif
1688  s.restoreState(state);
1689  throw CanteraError("equilibrate:estimateEP_Brinkley()",
1690  "Jacobian is singular. \nTry adding more species, "
1691  "changing the elemental composition slightly, \nor removing "
1692  "unused elements.");
1693  //return -3;
1694  }
1695 
1696  /*
1697  * Figure out the damping coefficient: Use a delta damping
1698  * coefficient formulation: magnitude of change is capped
1699  * to exp(1).
1700  */
1701  beta = 1.0;
1702  for (m = 0; m < m_mm; m++) {
1703  if (resid[m] > 1.0) {
1704  double betat = 1.0 / resid[m];
1705  if (betat < beta) {
1706  beta = betat;
1707  }
1708  }
1709  if (resid[m] < -1.0) {
1710  double betat = -1.0 / resid[m];
1711  if (betat < beta) {
1712  beta = betat;
1713  }
1714  }
1715  }
1716 #ifdef DEBUG_MODE
1717  if (ChemEquil_print_lvl > 0) {
1718  if (beta != 1.0) {
1719  writelogf("(it %d) Beta = %g\n", iter, beta);
1720  }
1721  }
1722 #endif
1723 
1724  /*
1725  * Update the solution vector
1726  */
1727 updateSolnVector:
1728  for (m = 0; m < m_mm; m++) {
1729  x[m] += beta * resid[m];
1730  }
1731  n_t *= exp(beta * resid[m_mm]);
1732 
1733 
1734 #ifdef DEBUG_MODE
1735  if (ChemEquil_print_lvl > 0) {
1736  writelogf("(it %d) OLD_SOLUTION NEW SOLUTION (undamped updated)\n", iter);
1737  for (m = 0; m < m_mm; m++) {
1738  string eee = eNames[m];
1739  writelogf(" %5s %10.5g %10.5g %10.5g\n", eee.c_str(), x_old[m], x[m], resid[m]);
1740  }
1741  writelogf(" n_t %10.5g %10.5g %10.5g \n", x_old[m_mm], n_t, exp(resid[m_mm]));
1742  }
1743 #endif
1744  }
1745 exit:
1746 #ifdef DEBUG_MODE
1747  if (ChemEquil_print_lvl > 0) {
1748  double temp = s.temperature();
1749  double pres = s.pressure();
1750 
1751  if (retn == 0) {
1752  writelogf(" ChemEquil::estimateEP_Brinkley() SUCCESS: equilibrium found at T = %g, Pres = %g\n",
1753  temp, pres);
1754  } else {
1755  writelogf(" ChemEquil::estimateEP_Brinkley() FAILURE: equilibrium not found at T = %g, Pres = %g\n",
1756  temp, pres);
1757  }
1758  }
1759 #endif
1760  return retn;
1761 }
1762 
1763 
1764 void ChemEquil::adjustEloc(thermo_t& s, vector_fp& elMolesGoal)
1765 {
1766  if (m_eloc == npos) {
1767  return;
1768  }
1769  if (fabs(elMolesGoal[m_eloc]) > 1.0E-20) {
1770  return;
1771  }
1773  size_t k;
1774 
1775 #ifdef DEBUG_MODE
1776  int maxPosEloc = -1;
1777  int maxNegEloc = -1;
1778  double maxPosVal = -1.0;
1779  double maxNegVal = -1.0;
1780  if (ChemEquil_print_lvl > 0) {
1781  for (k = 0; k < m_kk; k++) {
1782  if (nAtoms(k,m_eloc) > 0.0) {
1783  if (m_molefractions[k] > maxPosVal && m_molefractions[k] > 0.0) {
1784  maxPosVal = m_molefractions[k];
1785  maxPosEloc = k;
1786  }
1787  }
1788  if (nAtoms(k,m_eloc) < 0.0) {
1789  if (m_molefractions[k] > maxNegVal && m_molefractions[k] > 0.0) {
1790  maxNegVal = m_molefractions[k];
1791  maxNegEloc = k;
1792  }
1793  }
1794  }
1795  }
1796 #endif
1797 
1798  double sumPos = 0.0;
1799  double sumNeg = 0.0;
1800  for (k = 0; k < m_kk; k++) {
1801  if (nAtoms(k,m_eloc) > 0.0) {
1802  sumPos += nAtoms(k,m_eloc) * m_molefractions[k];
1803  }
1804  if (nAtoms(k,m_eloc) < 0.0) {
1805  sumNeg += nAtoms(k,m_eloc) * m_molefractions[k];
1806  }
1807  }
1808  sumNeg = - sumNeg;
1809 
1810  if (sumPos >= sumNeg) {
1811  if (sumPos <= 0.0) {
1812  return;
1813  }
1814  double factor = (elMolesGoal[m_eloc] + sumNeg) / sumPos;
1815 #ifdef DEBUG_MODE
1816  if (ChemEquil_print_lvl > 0) {
1817  if (factor < 0.9999999999) {
1818  string nnn = s.speciesName(maxPosEloc);
1819  writelogf("adjustEloc: adjusted %s and friends from %g to %g to ensure neutrality condition\n",
1820  nnn.c_str(),
1821  m_molefractions[maxPosEloc], m_molefractions[maxPosEloc]*factor);
1822  }
1823  }
1824 #endif
1825  for (k = 0; k < m_kk; k++) {
1826  if (nAtoms(k,m_eloc) > 0.0) {
1827  m_molefractions[k] *= factor;
1828  }
1829  }
1830  } else {
1831  double factor = (-elMolesGoal[m_eloc] + sumPos) / sumNeg;
1832 #ifdef DEBUG_MODE
1833  if (ChemEquil_print_lvl > 0) {
1834  if (factor < 0.9999999999) {
1835  string nnn = s.speciesName(maxNegEloc);
1836  writelogf("adjustEloc: adjusted %s and friends from %g to %g to ensure neutrality condition\n",
1837  nnn.c_str(),
1838  m_molefractions[maxNegEloc], m_molefractions[maxNegEloc]*factor);
1839  }
1840  }
1841 #endif
1842  for (k = 0; k < m_kk; k++) {
1843  if (nAtoms(k,m_eloc) < 0.0) {
1844  m_molefractions[k] *= factor;
1845  }
1846  }
1847  }
1848 
1851 
1852 }
1853 
1854 } // namespace
doublereal nAtoms(size_t k, size_t m) const
Number of atoms of element m in species k.
Definition: Phase.cpp:215
EquilOpt options
Options controlling how the calculation is carried out.
Definition: ChemEquil.h:150
virtual doublereal density() const
Density (kg/m^3).
Definition: Phase.h:534
std::string int2str(const int n, const std::string &fmt)
Convert an int to a string using a format converter.
Definition: stringUtils.cpp:40
void restoreState(const vector_fp &state)
Restore a state saved on a previous call to saveState.
Definition: Phase.cpp:289
vector_fp m_muSS_RT
Dimensionless values of the gibbs free energy for the standard state of each species, at the temperature and pressure of the solution (the star standard state).
Definition: ChemEquil.h:336
size_t BasisOptimize(int *usedZeroedSpecies, bool doFormRxn, MultiPhase *mphase, std::vector< size_t > &orderVectorSpecies, std::vector< size_t > &orderVectorElements, vector_fp &formRxnMatrix)
Choose the optimum basis of species for the equilibrium calculations.
thermo_t * m_phase
Pointer to the ThermoPhase object used to initialize this object.
Definition: ChemEquil.h:163
doublereal relTolerance
Relative tolerance.
Definition: ChemEquil.h:38
size_t nElements() const
Number of elements.
Definition: Phase.cpp:139
const doublereal OneAtm
One atmosphere [Pa].
Definition: ct_defs.h:71
void beginLogGroup(const std::string &title, int loglevel)
Create a new group for log messages.
Definition: global.cpp:90
int setInitialMoles(thermo_t &s, vector_fp &elMoleGoal, int loglevel=0)
Estimate the initial mole numbers.
Definition: ChemEquil.cpp:197
const size_t npos
index returned by functions to indicate "no position"
Definition: ct_defs.h:173
size_t m_kk
number of species in the phase
Definition: ChemEquil.h:276
void addPhase(ThermoPhase *p, doublereal moles)
Add a phase to the mixture.
Definition: MultiPhase.cpp:98
int equilibrate(thermo_t &s, const char *XY, bool useThermoPhaseElementPotentials=false, int loglevel=0)
Definition: ChemEquil.cpp:394
doublereal atomicWeight(size_t m) const
Atomic weight of element m.
Definition: Phase.cpp:179
This file contains definitions of terms that are used in internal routines and are unlikely to need m...
void initialize(thermo_t &s)
Definition: ChemEquil.cpp:72
size_t m_nComponents
This is equal to the rank of the stoichiometric coefficient matrix when it is computed.
Definition: ChemEquil.h:283
int solve(DenseMatrix &A, double *b)
Solve Ax = b. Array b is overwritten on exit with x.
This file contains definitions for utility functions and text for modules, inputfiles, logs, textlogs, HTML_logs (see Input File Handling, Diagnostic Output, Writing messages to the screen and Writing HTML Logfiles).
void getMoleFractions(doublereal *const x) const
Get the species mole fraction vector.
Definition: Phase.cpp:519
Headers for the MultiPhase object that is used to set up multiphase equilibrium problems (see Equilib...
void equilResidual(thermo_t &s, const vector_fp &x, const vector_fp &elmtotal, vector_fp &resid, double xval, double yval, int loglevel=0)
Evaluates the residual vector F, of length m_mm.
Definition: ChemEquil.cpp:1025
void update(const thermo_t &s)
Update internally stored state information.
Definition: ChemEquil.cpp:167
int iterations
Iteration counter.
Definition: ChemEquil.h:41
Base class for a phase with thermodynamic properties.
Definition: ThermoPhase.h:101
std::vector< int > vector_int
Vector of ints.
Definition: ct_defs.h:167
bool getElementPotentials(doublereal *lambda) const
Returns the element potentials stored in the ThermoPhase object.
void setToEquilState(thermo_t &s, const vector_fp &x, doublereal t)
Definition: ChemEquil.cpp:147
A class for multiphase mixtures.
Definition: MultiPhase.h:54
void writelogf(const char *fmt,...)
Write a formatted message to the screen.
Definition: global.cpp:48
size_t m_eloc
Index of the element id corresponding to the electric charge of each species.
Definition: ChemEquil.h:324
virtual doublereal maxTemp(size_t k=npos) const
Maximum temperature for which the thermodynamic data for the species are valid.
Definition: ThermoPhase.h:244
std::string fp2str(const double x, const std::string &fmt)
Convert a double into a c++ string.
Definition: stringUtils.cpp:29
virtual void getChemPotentials(doublereal *mu) const
Get the species chemical potentials. Units: J/kmol.
Definition: ThermoPhase.h:595
doublereal absElemTol
Abs Tol in element number.
Definition: ChemEquil.h:39
double m_elemFracCutoff
element fractional cutoff, below which the element will be zeroed.
Definition: ChemEquil.h:340
int estimateEP_Brinkley(thermo_t &s, vector_fp &lambda, vector_fp &elMoles)
Definition: ChemEquil.cpp:1180
Base class for exceptions thrown by Cantera classes.
Definition: ctexceptions.h:68
void setElementPotentials(const vector_fp &lambda)
Stores the element potentials in the ThermoPhase object.
virtual void setMoleFractions(const doublereal *const x)
Set the mole fractions to the specified values There is no restriction on the sum of the mole fractio...
Definition: Phase.cpp:306
const std::vector< std::string > & elementNames() const
Return a read-only reference to the vector of element names.
Definition: Phase.cpp:174
int estimateElementPotentials(thermo_t &s, vector_fp &lambda, vector_fp &elMolesGoal, int loglevel=0)
Generate a starting estimate for the element potentials.
Definition: ChemEquil.cpp:262
size_t m_mm
number of elements in the phase
Definition: ChemEquil.h:275
doublereal dot(InputIter x_begin, InputIter x_end, InputIter2 y_begin)
Function that calculates a templated inner product.
Definition: utilities.h:131
virtual doublereal refPressure() const
Returns the reference pressure in Pa.
Definition: ThermoPhase.h:159
void endLogGroup(const std::string &title)
Close the current group of log messages.
Definition: global.cpp:115
size_t nSpecies() const
Returns the number of species in the phase.
Definition: Phase.h:252
void addLogEntry(const std::string &tag, const std::string &value)
Add an entry to an HTML log file.
Definition: global.cpp:95
virtual doublereal pressure() const
Return the thermodynamic pressure (Pa).
Definition: ThermoPhase.h:314
doublereal moleFraction(size_t k) const
Return the mole fraction of a single species.
Definition: Phase.cpp:524
doublereal temperature() const
Temperature (K).
Definition: Phase.h:528
virtual void setPressure(doublereal p)
Set the internally stored pressure (Pa) at constant temperature and composition.
Definition: ThermoPhase.h:331
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:165
Templates for operations on vector-like objects.
virtual void setTemperature(const doublereal temp)
Set the internally stored temperature of the phase (K).
Definition: Phase.h:562
void scale(InputIter begin, InputIter end, OutputIter out, S scale_factor)
Multiply elements of an array by a scale factor.
Definition: utilities.h:154
vector_fp m_lambda
Current value of the dimensional element potentials -> length = m_mm.
Definition: ChemEquil.h:296
Headers for the DenseMatrix object, which deals with dense rectangular matrices and description of th...
const doublereal GasConstant
Universal Gas Constant. [J/kmol/K].
Definition: ct_defs.h:66
void init()
Process phases and build atomic composition array.
Definition: MultiPhase.cpp:173
virtual void getActivityCoefficients(doublereal *ac) const
Get the array of non-dimensional molar-based activity coefficients at the current solution temperatur...
Definition: ThermoPhase.h:553
doublereal nAtoms(size_t k, size_t m) const
number of atoms of element m in species k.
Definition: ChemEquil.h:166
Contains declarations for string manipulation functions within Cantera.
#define DATA_PTR(vec)
Creates a pointer to the start of the raw data for a vector.
Definition: ct_defs.h:36
std::string elementName(size_t m) const
Name of the element with index m.
Definition: Phase.cpp:158
void saveState(vector_fp &state) const
Save the current internal state of the phase Write to vector 'state' the current internal state...
Definition: Phase.cpp:277
vector_fp m_molefractions
Current value of the mole fractions in the single phase.
Definition: ChemEquil.h:291
double calcEmoles(thermo_t &s, vector_fp &x, const double &n_t, const vector_fp &Xmol_i_calc, vector_fp &eMolesCalc, vector_fp &n_i_calc, double pressureConst)
Given a vector of dimensionless element abundances, this routine calculates the moles of the elements...
Definition: ChemEquil.cpp:1140
int maxIterations
Maximum number of iterations.
Definition: ChemEquil.h:40
void writelog(const std::string &msg)
Write a message to the screen.
Definition: global.cpp:43
int _equilflag(const char *xy)
map property strings to integers
Definition: ChemEquil.cpp:29
size_t ElemRearrange(size_t nComponents, const vector_fp &elementAbundances, MultiPhase *mphase, std::vector< size_t > &orderVectorSpecies, std::vector< size_t > &orderVectorElements)
Handles the potential rearrangement of the constraint equations represented by the Formula Matrix...
virtual void setToEquilState(const doublereal *lambda_RT)
This method is used by the ChemEquil equilibrium solver.
Definition: ThermoPhase.h:1198
std::string speciesName(size_t k) const
Name of the species with index k.
Definition: Phase.cpp:246
Chemical equilibrium.
Definitions for the classes that are thrown when Cantera experiences an error condition (also contain...
virtual void getGibbs_RT(doublereal *grt) const
Get the nondimensional Gibbs functions for the species in their standard states at the current T and ...
Definition: ThermoPhase.h:741
A class for full (non-sparse) matrices with Fortran-compatible data storage, which adds matrix operat...
Definition: DenseMatrix.h:70
int dampStep(thermo_t &s, vector_fp &oldx, double oldf, vector_fp &grad, vector_fp &step, vector_fp &x, double &f, vector_fp &elmols, double xval, double yval)
Find an acceptable step size and take it.
Definition: ChemEquil.cpp:979
virtual doublereal minTemp(size_t k=npos) const
Minimum temperature for which the thermodynamic data for the species or phase are valid...
Definition: ThermoPhase.h:175
void save()
Function to put this error onto Cantera's error stack.