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