Cantera  2.0
MultiPhaseEquil.cpp
Go to the documentation of this file.
1 /**
2  * @file MultiPhaseEquil.cpp
3  */
4 #include "cantera/equil/MultiPhaseEquil.h"
7 #include "cantera/base/global.h"
9 
10 #include <math.h>
11 #include <iostream>
12 #include <cstdio>
13 
14 using namespace std;
15 
16 namespace Cantera
17 {
18 
19 const doublereal TINY = 1.0e-20;
20 
21 #if defined(WITH_HTML_LOGS)
22 /// Used to print reaction equations. Given a stoichiometric
23 /// coefficient 'nu' and a chemical symbol 'sym', return a string
24 /// for this species in the reaction.
25 /// @param first if this is false, then a " + " string will be
26 /// added to the beginning of the string.
27 /// @param nu Stoichiometric coefficient. May be positive or negative. The
28 /// absolute value will be used in the string.
29 /// @param sym Species chemical symbol.
30 ///
31 static string coeffString(bool first, doublereal nu, string sym)
32 {
33  if (nu == 0.0) {
34  return "";
35  }
36  string strt = " + ";
37  if (first) {
38  strt = "";
39  }
40  if (nu == 1.0 || nu == -1.0) {
41  return strt + sym;
42  }
43  string s = fp2str(fabs(nu));
44  return strt + s + " " + sym;
45 }
46 #endif
47 
48 /// Constructor. Construct a multiphase equilibrium manager for a
49 /// multiphase mixture.
50 /// @param mix Pointer to a multiphase mixture object.
51 /// @param start If true, the initial composition will be
52 /// determined by a linear Gibbs minimization, otherwise the
53 /// initial mixture composition will be used.
54 MultiPhaseEquil::MultiPhaseEquil(MultiPhase* mix, bool start, int loglevel) : m_mix(mix)
55 {
56  // the multi-phase mixture
57  // m_mix = mix;
58 
59  // store some mixture parameters locally
60  m_nel_mix = mix->nElements();
61  m_nsp_mix = mix->nSpecies();
62  m_np = mix->nPhases();
63  m_press = mix->pressure();
64  m_temp = mix->temperature();
65 
66  index_t m, k;
67  m_force = true;
68  m_nel = 0;
69  m_nsp = 0;
70  m_eloc = 1000;
71  m_incl_species.resize(m_nsp_mix,1);
72  m_incl_element.resize(m_nel_mix,1);
73  for (m = 0; m < m_nel_mix; m++) {
74  string enm = mix->elementName(m);
75  // element 'E' or 'e' represents an electron; this
76  // requires special handling, so save its index
77  // for later use
78  if (enm == "E" || enm == "e") {
79  m_eloc = m;
80  }
81  // if an element other than electrons is not present in
82  // the mixture, then exclude it and all species containing
83  // it from the calculation. Electrons are a special case,
84  // since a species can have a negative number of 'atoms'
85  // of electrons (positive ions).
86  if (m_mix->elementMoles(m) <= 0.0) {
87  if (m != m_eloc) {
88  m_incl_element[m] = 0;
89  for (k = 0; k < m_nsp_mix; k++) {
90  if (m_mix->nAtoms(k,m) != 0.0) {
91  m_incl_species[k] = 0;
92  }
93  }
94  }
95  }
96  }
97 
98  // Now build the list of elements to be included, starting with
99  // electrons, if they are present.
100  if (m_eloc < m_nel_mix) {
101  m_element.push_back(m_eloc);
102  m_nel++;
103  }
104  // add the included elements other than electrons
105  for (m = 0; m < m_nel_mix; m++) {
106  if (m_incl_element[m] == 1 && m != m_eloc) {
107  m_nel++;
108  m_element.push_back(m);
109  }
110  }
111 
112  // include pure single-constituent phases only if their thermo
113  // data are valid for this temperature. This is necessary,
114  // since some thermo polynomial fits are done only for a
115  // limited temperature range. For example, using the NASA
116  // polynomial fits for solid ice and liquid water, if this
117  // were not done the calculation would predict solid ice to be
118  // present far above its melting point, since the thermo
119  // polynomial fits only extend to 273.15 K, and give
120  // unphysical results above this temperature, leading
121  // (incorrectly) to Gibbs free energies at high temperature
122  // lower than for liquid water.
123  index_t ip;
124  for (k = 0; k < m_nsp_mix; k++) {
125  ip = m_mix->speciesPhaseIndex(k);
126  if (!m_mix->solutionSpecies(k) &&
127  !m_mix->tempOK(ip)) {
128  m_incl_species[k] = 0;
129  if (m_mix->speciesMoles(k) > 0.0) {
130  throw CanteraError("MultiPhaseEquil",
131  "condensed-phase species"+ m_mix->speciesName(k)
132  + " is excluded since its thermo properties are \n"
133  "not valid at this temperature, but it has "
134  "non-zero moles in the initial state.");
135  }
136  }
137  }
138 
139  // Now build the list of all species to be included in the
140  // calculation.
141  for (k = 0; k < m_nsp_mix; k++) {
142  if (m_incl_species[k] ==1) {
143  m_nsp++;
144  m_species.push_back(k);
145  }
146  }
147 
148  // some work arrays for internal use
149  m_work.resize(m_nsp);
150  m_work2.resize(m_nsp);
151  m_work3.resize(m_nsp_mix);
152  m_mu.resize(m_nsp_mix);
153 
154  // number of moles of each species
155  m_moles.resize(m_nsp);
156  m_lastmoles.resize(m_nsp);
157  m_dxi.resize(nFree());
158 
159  // initialize the mole numbers to the mixture composition
160  index_t ik;
161  for (ik = 0; ik < m_nsp; ik++) {
162  m_moles[ik] = m_mix->speciesMoles(m_species[ik]);
163  }
164 
165  // Delta G / RT for each reaction
166  m_deltaG_RT.resize(nFree(), 0.0);
167 
168  m_majorsp.resize(m_nsp);
169  m_sortindex.resize(m_nsp,0);
170  m_lastsort.resize(m_nel);
171  m_solnrxn.resize(nFree());
172  m_A.resize(m_nel, m_nsp, 0.0);
173  m_N.resize(m_nsp, nFree());
174  m_order.resize(std::max(m_nsp, m_nel), 0);
175  for (k = 0; k < m_nsp; k++) {
176  m_order[k] = k;
177  }
178 
179  // if the 'start' flag is set, estimate the initial mole
180  // numbers by doing a linear Gibbs minimization. In this case,
181  // only the elemental composition of the initial mixture state
182  // matters.
183  if (start) {
184  setInitialMoles(loglevel-1);
185  }
186  computeN();
187 
188  // Take a very small step in composition space, so that no
189  // species has precisely zero moles.
190  vector_fp dxi(nFree(), 1.0e-20);
191  if (!dxi.empty()) {
192  multiply(m_N, DATA_PTR(dxi), DATA_PTR(m_work));
193  unsort(m_work);
194  }
195 
196  for (k = 0; k < m_nsp; k++) {
197  m_moles[k] += m_work[k];
198  m_lastmoles[k] = m_moles[k];
199  if (m_mix->solutionSpecies(m_species[k])) {
200  m_dsoln.push_back(1);
201  } else {
202  m_dsoln.push_back(0);
203  }
204  }
205  m_force = false;
206  updateMixMoles();
207 
208  // At this point, the instance has been created, the species
209  // to be included have been determined, and an initial
210  // composition has been selected that has all non-zero mole
211  // numbers for the included species.
212 }
213 
214 
215 doublereal MultiPhaseEquil::equilibrate(int XY, doublereal err,
216  int maxsteps, int loglevel)
217 {
218  int i;
219  m_iter = 0;
220  string iterstr;
221  if (loglevel > 0) {
222  beginLogGroup("MultiPhaseEquil::equilibrate", loglevel);
223  }
224 
225  for (i = 0; i < maxsteps; i++) {
226  if (loglevel > 0) {
227  iterstr = "iteration "+int2str(i);
228  beginLogGroup(iterstr);
229  }
230  stepComposition(loglevel-1);
231  if (loglevel > 0) {
232  addLogEntry("error",fp2str(error()));
233  endLogGroup(iterstr);
234  }
235  if (error() < err) {
236  break;
237  }
238  }
239  if (i >= maxsteps) {
240  if (loglevel > 0) {
241  addLogEntry("Error","no convergence in "+int2str(maxsteps)
242  +" iterations");
243  endLogGroup("MultiPhaseEquil::equilibrate");
244  }
245  throw CanteraError("MultiPhaseEquil::equilibrate",
246  "no convergence in " + int2str(maxsteps) +
247  " iterations. Error = " + fp2str(error()));
248  }
249  if (loglevel > 0) {
250  addLogEntry("iterations",int2str(iterations()));
251  addLogEntry("error tolerance",fp2str(err));
252  addLogEntry("error",fp2str(error()));
253  endLogGroup("MultiPhaseEquil::equilibrate");
254  }
255  finish();
256  return error();
257 }
258 
259 void MultiPhaseEquil::updateMixMoles()
260 {
261  fill(m_work3.begin(), m_work3.end(), 0.0);
262  index_t k;
263  for (k = 0; k < m_nsp; k++) {
264  m_work3[m_species[k]] = m_moles[k];
265  }
266  m_mix->setMoles(DATA_PTR(m_work3));
267 }
268 
269 /// Clean up the composition. The solution algorithm can leave
270 /// some species in stoichiometric condensed phases with very
271 /// small negative mole numbers. This method simply sets these to
272 /// zero.
274 {
275  fill(m_work3.begin(), m_work3.end(), 0.0);
276  index_t k;
277  for (k = 0; k < m_nsp; k++) {
278  m_work3[m_species[k]] = (m_moles[k] > 0.0 ? m_moles[k] : 0.0);
279  }
280  m_mix->setMoles(DATA_PTR(m_work3));
281 }
282 
283 
284 /// Estimate the initial mole numbers. This is done by running
285 /// each reaction as far forward or backward as possible, subject
286 /// to the constraint that all mole numbers remain
287 /// non-negative. Reactions for which \f$ \Delta \mu^0 \f$ are
288 /// positive are run in reverse, and ones for which it is negative
289 /// are run in the forward direction. The end result is equivalent
290 /// to solving the linear programming problem of minimizing the
291 /// linear Gibbs function subject to the element and
292 /// non-negativity constraints.
294 {
295  index_t ik, j;
296 
297  double not_mu = 1.0e12;
298  if (loglevel > 0) {
299  beginLogGroup("MultiPhaseEquil::setInitialMoles");
300  }
301 
302  m_mix->getValidChemPotentials(not_mu, DATA_PTR(m_mu), true);
303  doublereal dg_rt;
304 
305  int idir;
306  double nu;
307  double delta_xi, dxi_min = 1.0e10;
308  bool redo = true;
309  int iter = 0;
310 
311  while (redo) {
312 
313  // choose a set of components based on the current
314  // composition
315  computeN();
316  if (loglevel > 0) {
317  addLogEntry("iteration",iter);
318  }
319  redo = false;
320  iter++;
321  if (iter > 4) {
322  break;
323  }
324 
325  // loop over all reactions
326  for (j = 0; j < nFree(); j++) {
327  dg_rt = 0.0;
328  dxi_min = 1.0e10;
329  for (ik = 0; ik < m_nsp; ik++) {
330  dg_rt += mu(ik) * m_N(ik,j);
331  }
332 
333  // fwd or rev direction
334  idir = (dg_rt < 0.0 ? 1 : -1);
335 
336  for (ik = 0; ik < m_nsp; ik++) {
337  nu = m_N(ik, j);
338 
339  // set max change in progress variable by
340  // non-negativity requirement
341  // -> Note, 0.99 factor is so that difference of 2 numbers
342  // isn't zero. This causes differences between
343  // optimized and debug versions of the code
344  if (nu*idir < 0) {
345  delta_xi = fabs(0.99*moles(ik)/nu);
346  // if a component has nearly zero moles, redo
347  // with a new set of components
348  if (!redo && delta_xi < 1.0e-10 && ik < m_nel) {
349  if (loglevel > 0) {
350  addLogEntry("component too small",speciesName(ik));
351  }
352  redo = true;
353  }
354  if (delta_xi < dxi_min) {
355  dxi_min = delta_xi;
356  }
357  }
358  }
359  // step the composition by dxi_min
360  for (ik = 0; ik < m_nsp; ik++) {
361  moles(ik) += m_N(ik, j) * idir*dxi_min;
362  }
363  }
364  // set the moles of the phase objects to match
365  updateMixMoles();
366  }
367  for (ik = 0; ik < m_nsp; ik++)
368  if (moles(ik) != 0.0) {
369  addLogEntry(speciesName(ik), moles(ik));
370  }
371  if (loglevel > 0) {
372  endLogGroup("MultiPhaseEquil::setInitialMoles");
373  }
374  return 0;
375 }
376 
377 
378 /// This method finds a set of component species and a complete
379 /// set of formation reactions for the non-components in terms of
380 /// the components. Note that in most cases, many different
381 /// component sets are possible, and therefore neither the
382 /// components returned by this method nor the formation
383 /// reactions are unique. The algorithm used here is described in
384 /// Smith and Missen, Chemical Reaction Equilibrium Analysis.
385 ///
386 /// The component species are taken to be the first M species
387 /// in array 'species' that have linearly-independent compositions.
388 ///
389 /// @param order On entry, vector \a order should contain species
390 /// index numbers in the order of decreasing desirability as a
391 /// component. For example, if it is desired to choose the
392 /// components from among the major species, this array might
393 /// list species index numbers in decreasing order of mole
394 /// fraction. If array 'species' does not have length =
395 /// nSpecies(), then the species will be considered as candidates
396 /// to be components in declaration order, beginning with the
397 /// first phase added.
398 ///
399 void MultiPhaseEquil::getComponents(const std::vector<size_t>& order)
400 {
401  index_t m, k, j;
402 
403  // if the input species array has the wrong size, ignore it
404  // and consider the species for components in declaration order.
405  if (order.size() != m_nsp) {
406  for (k = 0; k < m_nsp; k++) {
407  m_order[k] = k;
408  }
409  } else {
410  for (k = 0; k < m_nsp; k++) {
411  m_order[k] = order[k];
412  }
413  }
414 
415  index_t nRows = m_nel;
416  index_t nColumns = m_nsp;
417  doublereal fctr;
418 
419  // set up the atomic composition matrix
420  for (m = 0; m < nRows; m++) {
421  for (k = 0; k < nColumns; k++) {
422  m_A(m, k) = m_mix->nAtoms(m_species[m_order[k]], m_element[m]);
423  }
424  }
425 
426  // Do Gaussian elimination
427  for (m = 0; m < nRows; m++) {
428  // Check for rows that are zero
429  bool isZeroRow = true;
430  for (k = m; k < nColumns; k++) {
431  if (fabs(m_A(m,k)) > sqrt(TINY)) {
432  isZeroRow = false;
433  break;
434  }
435  }
436  if (isZeroRow) {
437  // Find the last non-zero row
438  index_t n = nRows - 1;
439  bool foundSwapCandidate = false;
440  for (; n > m; n--) {
441  for (k = m; k < nColumns; k++) {
442  if (fabs(m_A(n,k)) > sqrt(TINY)) {
443  foundSwapCandidate = true;
444  break;
445  }
446  }
447  if (foundSwapCandidate) {
448  break;
449  }
450  }
451  if (m != n) {
452  // Swap this row with the last non-zero row
453  for (k = 0; k < nColumns; k++) {
454  std::swap(m_A(n,k), m_A(m,k));
455  }
456  } else {
457  // All remaining rows are zero. Elimination is complete.
458  break;
459  }
460  }
461 
462  // If a pivot is zero, exchange columns. This occurs when
463  // a species has an elemental composition that is not
464  // linearly independent of the component species that have
465  // already been assigned
466  if (m < nColumns && m_A(m,m) == 0.0) {
467  // First, we need to find a good candidate for a
468  // component species to swap in for the one that has
469  // zero pivot. It must contain element m, be linearly
470  // independent of the components processed so far
471  // (m_A(m,k) != 0), and should be a major species if
472  // possible. We'll choose the species with greatest
473  // mole fraction that satisfies these criteria.
474  doublereal maxmoles = -999.0;
475  index_t kmax = 0;
476  for (k = m+1; k < nColumns; k++) {
477  if (m_A(m,k) != 0.0) {
478  if (fabs(m_moles[m_order[k]]) > maxmoles) {
479  kmax = k;
480  maxmoles = fabs(m_moles[m_order[k]]);
481  }
482  }
483  }
484 
485  // Now exchange the column with zero pivot with the
486  // column for this major species
487  for (size_t n = 0; n < nRows; n++) {
488  std::swap(m_A(n, m), m_A(n, kmax));
489  }
490 
491  // exchange the species labels on the columns
492  std::swap(m_order[m], m_order[kmax]);
493  }
494 
495  // scale row m so that the diagonal element is unity
496  fctr = 1.0/m_A(m,m);
497  for (k = 0; k < nColumns; k++) {
498  m_A(m,k) *= fctr;
499  }
500 
501  // For all rows below the diagonal, subtract A(n,m)/A(m,m)
502  // * (row m) from row n, so that A(n,m) = 0.
503  for (size_t n = m+1; n < m_nel; n++) {
504  fctr = m_A(n,m)/m_A(m,m);
505  for (k = 0; k < m_nsp; k++) {
506  m_A(n,k) -= m_A(m,k)*fctr;
507  }
508  }
509  }
510 
511  // The left m_nel columns of A are now upper-diagonal. Now
512  // reduce the m_nel columns to diagonal form by back-solving
513  for (m = std::min(nRows,nColumns)-1; m > 0; m--) {
514  for (size_t n = m-1; n != npos; n--) {
515  if (m_A(n,m) != 0.0) {
516  fctr = m_A(n,m);
517  for (k = m; k < m_nsp; k++) {
518  m_A(n,k) -= fctr*m_A(m,k);
519  }
520  }
521  }
522  }
523 
524  // create stoichiometric coefficient matrix.
525  for (size_t n = 0; n < m_nsp; n++) {
526  if (n < m_nel)
527  for (k = 0; k < nFree(); k++) {
528  m_N(n, k) = -m_A(n, k + m_nel);
529  }
530  else {
531  for (k = 0; k < nFree(); k++) {
532  m_N(n, k) = 0.0;
533  }
534  m_N(n, n - m_nel) = 1.0;
535  }
536  }
537 
538  // find reactions involving solution phase species
539  for (j = 0; j < nFree(); j++) {
540  m_solnrxn[j] = false;
541  for (k = 0; k < m_nsp; k++) {
542  if (m_N(k, j) != 0)
543  if (m_mix->solutionSpecies(m_species[m_order[k]])) {
544  m_solnrxn[j] = true;
545  }
546  }
547  }
548 }
549 
550 
551 
552 
553 /// Re-arrange a vector of species properties in sorted form
554 /// (components first) into unsorted, sequential form.
556 {
557  copy(x.begin(), x.end(), m_work2.begin());
558  index_t k;
559  for (k = 0; k < m_nsp; k++) {
560  x[m_order[k]] = m_work2[k];
561  }
562 }
563 
564 #if defined(WITH_HTML_LOGS)
565 void MultiPhaseEquil::printInfo(int loglevel)
566 {
567  index_t m, ik, k;
568  if (loglevel > 0) {
569  beginLogGroup("info");
570  beginLogGroup("components");
571  }
572  for (m = 0; m < m_nel; m++) {
573  ik = m_order[m];
574  k = m_species[ik];
575  if (loglevel > 0) {
576  addLogEntry(m_mix->speciesName(k), fp2str(m_moles[ik]));
577  }
578  }
579  if (loglevel > 0) {
580  endLogGroup("components");
581  beginLogGroup("non-components");
582  }
583  for (m = m_nel; m < m_nsp; m++) {
584  ik = m_order[m];
585  k = m_species[ik];
586  if (loglevel > 0) {
587  addLogEntry(m_mix->speciesName(k), fp2str(m_moles[ik]));
588  }
589  }
590  if (loglevel > 0) {
591  endLogGroup("non-components");
592  addLogEntry("Error",fp2str(error()));
593  beginLogGroup("Delta G / RT");
594  }
595  for (k = 0; k < nFree(); k++) {
596  if (loglevel > 0) {
597  addLogEntry(reactionString(k), fp2str(m_deltaG_RT[k]));
598  }
599  }
600  if (loglevel > 0) {
601  endLogGroup("Delta G / RT");
602  endLogGroup("info");
603  }
604 }
605 
606 /// Return a string specifying the jth reaction.
608 {
609  string sr = "", sp = "";
610  index_t i, k;
611  bool rstrt = true;
612  bool pstrt = true;
613  doublereal nu;
614  for (i = 0; i < m_nsp; i++) {
615  nu = m_N(i, j);
616  k = m_species[m_order[i]];
617  if (nu < 0.0) {
618  sr += coeffString(rstrt, nu, m_mix->speciesName(k));
619  rstrt = false;
620  }
621  if (nu > 0.0) {
622  sp += coeffString(pstrt, nu, m_mix->speciesName(k));
623  pstrt = false;
624  }
625  }
626  return sr + " <=> " + sp;
627 }
628 #endif
629 
630 void MultiPhaseEquil::step(doublereal omega, vector_fp& deltaN,
631  int loglevel)
632 {
633  index_t k, ik;
634  if (loglevel > 0) {
635  beginLogGroup("MultiPhaseEquil::step");
636  }
637  if (omega < 0.0) {
638  throw CanteraError("step","negative omega");
639  }
640 
641  for (ik = 0; ik < m_nel; ik++) {
642  k = m_order[ik];
643  m_lastmoles[k] = m_moles[k];
644  if (loglevel > 0) {
645  addLogEntry("component "+m_mix->speciesName(m_species[k])+" moles",
646  m_moles[k]);
647  addLogEntry("component "+m_mix->speciesName(m_species[k])+" step",
648  omega*deltaN[k]);
649  }
650  m_moles[k] += omega * deltaN[k];
651  }
652 
653  for (ik = m_nel; ik < m_nsp; ik++) {
654  k = m_order[ik];
655  m_lastmoles[k] = m_moles[k];
656  if (m_majorsp[k]) {
657  m_moles[k] += omega * deltaN[k];
658  } else {
659  m_moles[k] = fabs(m_moles[k])*std::min(10.0,
660  exp(-m_deltaG_RT[ik - m_nel]));
661  }
662  }
663  updateMixMoles();
664  if (loglevel > 0) {
665  endLogGroup("MultiPhaseEquil::step");
666  }
667 }
668 
669 
670 /// Take one step in composition, given the gradient of G at the
671 /// starting point, and a vector of reaction steps dxi.
672 doublereal MultiPhaseEquil::
673 stepComposition(int loglevel)
674 {
675  if (loglevel > 0) {
676  beginLogGroup("MultiPhaseEquil::stepComposition");
677  }
678 
679  m_iter++;
680  index_t ik, k = 0;
681  doublereal grad0 = computeReactionSteps(m_dxi);
682 
683  // compute the mole fraction changes.
684  if (nFree()) {
685  multiply(m_N, DATA_PTR(m_dxi), DATA_PTR(m_work));
686  }
687 
688  // change to sequential form
689  unsort(m_work);
690 
691  // scale omega to keep the major species non-negative
692  doublereal FCTR = 0.99;
693  const doublereal MAJOR_THRESHOLD = 1.0e-12;
694 
695  doublereal omega = 1.0, omax, omegamax = 1.0;
696  for (ik = 0; ik < m_nsp; ik++) {
697  k = m_order[ik];
698  if (ik < m_nel) {
699  FCTR = 0.99;
700  if (m_moles[k] < MAJOR_THRESHOLD) {
701  m_force = true;
702  }
703  } else {
704  FCTR = 0.9;
705  }
706  // if species k is in a multi-species solution phase, then its
707  // mole number must remain positive, unless the entire phase
708  // goes away. First we'll determine an upper bound on omega,
709  // such that all
710  if (m_dsoln[k] == 1) {
711 
712  if ((m_moles[k] > MAJOR_THRESHOLD) || (ik < m_nel)) {
713  if (m_moles[k] < MAJOR_THRESHOLD) {
714  m_force = true;
715  }
716  omax = m_moles[k]*FCTR/(fabs(m_work[k]) + TINY);
717  if (m_work[k] < 0.0 && omax < omegamax) {
718  omegamax = omax;
719  if (omegamax < 1.0e-5) {
720  m_force = true;
721  }
722  }
723  m_majorsp[k] = true;
724  } else {
725  m_majorsp[k] = false;
726  }
727  } else {
728  if (m_work[k] < 0.0 && m_moles[k] > 0.0) {
729  omax = -m_moles[k]/m_work[k];
730  if (omax < omegamax) {
731  omegamax = omax; //*1.000001;
732  if (omegamax < 1.0e-5) {
733  m_force = true;
734  }
735  }
736  }
737  if (m_moles[k] < -Tiny) {
738  if (loglevel > 0)
739  addLogEntry("Negative moles for "
740  +m_mix->speciesName(m_species[k]), fp2str(m_moles[k]));
741  }
742  m_majorsp[k] = true;
743  }
744  }
745 
746  // now take a step with this scaled omega
747  if (loglevel > 0) {
748  addLogEntry("Stepping by ", fp2str(omegamax));
749  }
750  step(omegamax, m_work);
751  // compute the gradient of G at this new position in the
752  // current direction. If it is positive, then we have overshot
753  // the minimum. In this case, interpolate back.
754  doublereal not_mu = 1.0e12;
755  m_mix->getValidChemPotentials(not_mu, DATA_PTR(m_mu));
756  doublereal grad1 = 0.0;
757  for (k = 0; k < m_nsp; k++) {
758  grad1 += m_work[k] * m_mu[m_species[k]];
759  }
760 
761  omega = omegamax;
762  if (grad1 > 0.0) {
763  omega *= fabs(grad0) / (grad1 + fabs(grad0));
764  for (k = 0; k < m_nsp; k++) {
765  m_moles[k] = m_lastmoles[k];
766  }
767  if (loglevel > 0) {
768  addLogEntry("Stepped over minimum. Take smaller step ", fp2str(omega));
769  }
770  step(omega, m_work);
771  }
772  printInfo(loglevel);
773  if (loglevel > 0) {
774  endLogGroup("MultiPhaseEquil::stepComposition");
775  }
776  return omega;
777 }
778 
779 
780 /// Compute the change in extent of reaction for each reaction.
781 
783 {
784 
785  index_t j, k, ik, kc, ip;
786  doublereal stoich, nmoles, csum, term1, fctr, rfctr;
787  vector_fp nu;
788  const doublereal TINY = 1.0e-20;
789  doublereal grad = 0.0;
790 
791  dxi.resize(nFree());
792  computeN();
793  doublereal not_mu = 1.0e12;
794  m_mix->getValidChemPotentials(not_mu, DATA_PTR(m_mu));
795 
796  for (j = 0; j < nFree(); j++) {
797 
798  // get stoichiometric vector
799  getStoichVector(j, nu);
800 
801  // compute Delta G
802  doublereal dg_rt = 0.0;
803  for (k = 0; k < m_nsp; k++) {
804  dg_rt += m_mu[m_species[k]] * nu[k];
805  }
806  dg_rt /= (m_temp * GasConstant);
807 
808  m_deltaG_RT[j] = dg_rt;
809  fctr = 1.0;
810 
811  // if this is a formation reaction for a single-component phase,
812  // check whether reaction should be included
813  ik = j + m_nel;
814  k = m_order[ik];
815  if (!m_dsoln[k]) {
816  if (m_moles[k] <= 0.0 && dg_rt > 0.0) {
817  fctr = 0.0;
818  } else {
819  fctr = 0.5;
820  }
821  } else if (!m_solnrxn[j]) {
822  fctr = 1.0;
823  } else {
824 
825  // component sum
826  csum = 0.0;
827  for (k = 0; k < m_nel; k++) {
828  kc = m_order[k];
829  stoich = nu[kc];
830  nmoles = fabs(m_mix->speciesMoles(m_species[kc])) + TINY;
831  csum += stoich*stoich*m_dsoln[kc]/nmoles;
832  }
833 
834  // noncomponent term
835  kc = m_order[j + m_nel];
836  nmoles = fabs(m_mix->speciesMoles(m_species[kc])) + TINY;
837  term1 = m_dsoln[kc]/nmoles;
838 
839  // sum over solution phases
840  doublereal sum = 0.0, psum;
841  for (ip = 0; ip < m_np; ip++) {
842  ThermoPhase& p = m_mix->phase(ip);
843  if (p.nSpecies() > 1) {
844  psum = 0.0;
845  for (k = 0; k < m_nsp; k++) {
846  kc = m_species[k];
847  if (m_mix->speciesPhaseIndex(kc) == ip) {
848  // bug fixed 7/12/06 DGG
849  stoich = nu[k]; // nu[kc];
850  psum += stoich * stoich;
851  }
852  }
853  sum -= psum / (fabs(m_mix->phaseMoles(ip)) + TINY);
854  }
855  }
856  rfctr = term1 + csum + sum;
857  if (fabs(rfctr) < TINY) {
858  fctr = 1.0;
859  } else {
860  fctr = 1.0/(term1 + csum + sum);
861  }
862  }
863  dxi[j] = -fctr*dg_rt;
864 
865  index_t m;
866  for (m = 0; m < m_nel; m++) {
867  if (m_moles[m_order[m]] <= 0.0 && (m_N(m, j)*dxi[j] < 0.0)) {
868  dxi[j] = 0.0;
869  }
870  }
871  grad += dxi[j]*dg_rt;
872 
873  }
874  return grad*GasConstant*m_temp;
875 }
876 
877 void MultiPhaseEquil::computeN()
878 {
879  // Sort the list of species by mole fraction (decreasing order)
880  std::vector<std::pair<double, size_t> > moleFractions(m_nsp);
881  for (size_t k = 0; k < m_nsp; k++) {
882  // use -Xk to generate reversed sort order
883  moleFractions[k].first = - m_mix->speciesMoles(m_species[k]);
884  moleFractions[k].second = k;
885  }
886  std::sort(moleFractions.begin(), moleFractions.end());
887  for (size_t k = 0; k < m_nsp; k++) {
888  m_sortindex[k] = moleFractions[k].second;
889  }
890 
891  bool ok;
892  for (size_t m = 0; m < m_nel; m++) {
893  size_t k = 0;
894  for (size_t ik = 0; ik < m_nsp; ik++) {
895  k = m_sortindex[ik];
896  if (m_mix->nAtoms(m_species[k],m_element[m]) != 0) {
897  break;
898  }
899  }
900  ok = false;
901  for (size_t ij = 0; ij < m_nel; ij++) {
902  if (k == m_order[ij]) {
903  ok = true;
904  }
905  }
906  if (!ok || m_force) {
907  getComponents(m_sortindex);
908  m_force = true;
909  break;
910  }
911  }
912 }
913 
914 doublereal MultiPhaseEquil::error()
915 {
916  doublereal err, maxerr = 0.0;
917 
918  // examine every reaction
919  for (size_t j = 0; j < nFree(); j++) {
920  size_t ik = j + m_nel;
921 
922  // don't require formation reactions for solution species
923  // present in trace amounts to be equilibrated
924  if (!isStoichPhase(ik) && fabs(moles(ik)) <= SmallNumber) {
925  err = 0.0;
926  }
927 
928  // for stoichiometric phase species, no error if not present and
929  // delta G for the formation reaction is positive
930  if (isStoichPhase(ik) && moles(ik) <= 0.0 &&
931  m_deltaG_RT[j] >= 0.0) {
932  err = 0.0;
933  } else {
934  err = fabs(m_deltaG_RT[j]);
935  }
936  if (err > maxerr) {
937  maxerr = err;
938  }
939  }
940  return maxerr;
941 }
942 
943 double MultiPhaseEquil::phaseMoles(index_t iph) const
944 {
945  return m_mix->phaseMoles(iph);
946 }
947 
948 void MultiPhaseEquil::reportCSV(const std::string& reportFile)
949 {
950  size_t k;
951  size_t istart;
952  size_t nSpecies;
953 
954  double vol = 0.0;
955  string sName;
956  size_t nphase = m_np;
957 
958  FILE* FP = fopen(reportFile.c_str(), "w");
959  if (!FP) {
960  printf("Failure to open file\n");
961  exit(EXIT_FAILURE);
962  }
963  double Temp = m_mix->temperature();
964  double pres = m_mix->pressure();
965  vector<double> mf(m_nsp_mix, 1.0);
966  vector<double> fe(m_nsp_mix, 0.0);
967 
968  std::vector<double> VolPM;
969  std::vector<double> activity;
970  std::vector<double> ac;
971  std::vector<double> mu;
972  std::vector<double> mu0;
973  std::vector<double> molalities;
974 
975 
976  vol = 0.0;
977  for (size_t iphase = 0; iphase < nphase; iphase++) {
978  istart = m_mix->speciesIndex(0, iphase);
979  ThermoPhase& tref = m_mix->phase(iphase);
980  nSpecies = tref.nSpecies();
981  VolPM.resize(nSpecies, 0.0);
982  tref.getMoleFractions(&mf[istart]);
983  tref.getPartialMolarVolumes(DATA_PTR(VolPM));
984  //vcs_VolPhase *volP = m_vprob->VPhaseList[iphase];
985 
986  double TMolesPhase = phaseMoles(iphase);
987  double VolPhaseVolumes = 0.0;
988  for (k = 0; k < nSpecies; k++) {
989  VolPhaseVolumes += VolPM[k] * mf[istart + k];
990  }
991  VolPhaseVolumes *= TMolesPhase;
992  vol += VolPhaseVolumes;
993  }
994  fprintf(FP,"--------------------- VCS_MULTIPHASE_EQUIL FINAL REPORT"
995  " -----------------------------\n");
996  fprintf(FP,"Temperature = %11.5g kelvin\n", Temp);
997  fprintf(FP,"Pressure = %11.5g Pascal\n", pres);
998  fprintf(FP,"Total Volume = %11.5g m**3\n", vol);
999  // fprintf(FP,"Number Basis optimizations = %d\n", m_vprob->m_NumBasisOptimizations);
1000  // fprintf(FP,"Number VCS iterations = %d\n", m_vprob->m_Iterations);
1001 
1002  for (size_t iphase = 0; iphase < nphase; iphase++) {
1003  istart = m_mix->speciesIndex(0, iphase);
1004 
1005  ThermoPhase& tref = m_mix->phase(iphase);
1006  ThermoPhase* tp = &tref;
1007  tp->getMoleFractions(&mf[istart]);
1008  string phaseName = tref.name();
1009  // vcs_VolPhase *volP = m_vprob->VPhaseList[iphase];
1010  double TMolesPhase = phaseMoles(iphase);
1011  //AssertTrace(TMolesPhase == m_mix->phaseMoles(iphase));
1012  nSpecies = tref.nSpecies();
1013  activity.resize(nSpecies, 0.0);
1014  ac.resize(nSpecies, 0.0);
1015 
1016  mu0.resize(nSpecies, 0.0);
1017  mu.resize(nSpecies, 0.0);
1018  VolPM.resize(nSpecies, 0.0);
1019  molalities.resize(nSpecies, 0.0);
1020 
1021  int actConvention = tp->activityConvention();
1022  tp->getActivities(DATA_PTR(activity));
1023  tp->getActivityCoefficients(DATA_PTR(ac));
1024  tp->getStandardChemPotentials(DATA_PTR(mu0));
1025 
1026  tp->getPartialMolarVolumes(DATA_PTR(VolPM));
1027  tp->getChemPotentials(DATA_PTR(mu));
1028  double VolPhaseVolumes = 0.0;
1029  for (k = 0; k < nSpecies; k++) {
1030  VolPhaseVolumes += VolPM[k] * mf[istart + k];
1031  }
1032  VolPhaseVolumes *= TMolesPhase;
1033  vol += VolPhaseVolumes;
1034  if (actConvention == 1) {
1035  MolalityVPSSTP* mTP = static_cast<MolalityVPSSTP*>(tp);
1036  mTP->getMolalities(DATA_PTR(molalities));
1037  tp->getChemPotentials(DATA_PTR(mu));
1038 
1039  if (iphase == 0) {
1040  fprintf(FP," Name, Phase, PhaseMoles, Mole_Fract, "
1041  "Molalities, ActCoeff, Activity,"
1042  "ChemPot_SS0, ChemPot, mole_num, PMVol, Phase_Volume\n");
1043 
1044  fprintf(FP," , , (kmol), , "
1045  ", , ,"
1046  " (kJ/gmol), (kJ/gmol), (kmol), (m**3/kmol), (m**3)\n");
1047  }
1048  for (k = 0; k < nSpecies; k++) {
1049  sName = tp->speciesName(k);
1050  fprintf(FP,"%12s, %11s, %11.3e, %11.3e, %11.3e, %11.3e, %11.3e,"
1051  "%11.3e, %11.3e, %11.3e, %11.3e, %11.3e\n",
1052  sName.c_str(),
1053  phaseName.c_str(), TMolesPhase,
1054  mf[istart + k], molalities[k], ac[k], activity[k],
1055  mu0[k]*1.0E-6, mu[k]*1.0E-6,
1056  mf[istart + k] * TMolesPhase,
1057  VolPM[k], VolPhaseVolumes);
1058  }
1059 
1060  } else {
1061  if (iphase == 0) {
1062  fprintf(FP," Name, Phase, PhaseMoles, Mole_Fract, "
1063  "Molalities, ActCoeff, Activity,"
1064  " ChemPotSS0, ChemPot, mole_num, PMVol, Phase_Volume\n");
1065 
1066  fprintf(FP," , , (kmol), , "
1067  ", , ,"
1068  " (kJ/gmol), (kJ/gmol), (kmol), (m**3/kmol), (m**3)\n");
1069  }
1070  for (k = 0; k < nSpecies; k++) {
1071  molalities[k] = 0.0;
1072  }
1073  for (k = 0; k < nSpecies; k++) {
1074  sName = tp->speciesName(k);
1075  fprintf(FP,"%12s, %11s, %11.3e, %11.3e, %11.3e, %11.3e, %11.3e, "
1076  "%11.3e, %11.3e,% 11.3e, %11.3e, %11.3e\n",
1077  sName.c_str(),
1078  phaseName.c_str(), TMolesPhase,
1079  mf[istart + k], molalities[k], ac[k],
1080  activity[k], mu0[k]*1.0E-6, mu[k]*1.0E-6,
1081  mf[istart + k] * TMolesPhase,
1082  VolPM[k], VolPhaseVolumes);
1083  }
1084  }
1085 #ifdef DEBUG_MODE
1086  /*
1087  * Check consistency: These should be equal
1088  */
1089  tp->getChemPotentials(&(fe[istart]));
1090  for (k = 0; k < nSpecies; k++) {
1091  //if (!vcs_doubleEqual(fe[istart+k], mu[k])) {
1092  // fprintf(FP,"ERROR: incompatibility!\n");
1093  // fclose(FP);
1094  // printf("ERROR: incompatibility!\n");
1095  // exit(EXIT_FAILURE);
1096  // }
1097  }
1098 #endif
1099 
1100  }
1101  fclose(FP);
1102 }
1103 
1104 
1105 
1106 
1107 }