Cantera  2.2.1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
vcs_solve.h
Go to the documentation of this file.
1 /**
2  * @file vcs_solve.h Header file for the internal object that holds the vcs
3  * equilibrium problem (see Class \link Cantera::VCS_SOLVE
4  * VCS_SOLVE\endlink and \ref equilfunctions ).
5  */
6 /*
7  * Copyright (2005) Sandia Corporation. Under the terms of
8  * Contract DE-AC04-94AL85000 with Sandia Corporation, the
9  * U.S. Government retains certain rights in this software.
10  */
11 
12 #ifndef _VCS_SOLVE_H
13 #define _VCS_SOLVE_H
14 
15 /*
16 * Index of Symbols
17 * -------------------
18 * irxn -> refers to the species or rxn between the species and
19 * the components in the problem
20 * k -> refers to the species
21 * j -> refers to the element or component
22 *
23 * ### -> to be eliminated
24 */
25 
26 #include "cantera/base/ct_defs.h"
27 #include "cantera/equil/vcs_defs.h"
29 #include "cantera/base/Array.h"
30 
31 namespace Cantera
32 {
33 /*
34  * Forward references
35  */
36 class vcs_VolPhase;
37 class VCS_SPECIES_THERMO;
38 class VCS_PROB;
39 class VCS_COUNTERS;
40 
41 
42 //! This is the main structure used to hold the internal data
43 //! used in vcs_solve_TP(), and to solve TP systems.
44 /*!
45  * The indices of information in this structure may change when the species
46  * basis changes or when phases pop in and out of existence. Both of these
47  * operations change the species ordering.
48  */
49 class VCS_SOLVE
50 {
51 public:
52  VCS_SOLVE();
53 
54  ~VCS_SOLVE();
55 
56  //! Initialize the sizes within the VCS_SOLVE object
57  /*!
58  * This resizes all of the internal arrays within the object. This
59  * routine operates in two modes. If all of the parameters are the same
60  * as it currently exists in the object, nothing is done by this routine;
61  * a quick exit is carried out and all of the data in the object
62  * persists.
63  *
64  * If any of the parameters are different than currently exists in the
65  * object, then all of the data in the object must be redone. It may not
66  * be zeroed, but it must be redone.
67  *
68  * @param nspecies0 Number of species within the object
69  * @param nelements Number of element constraints within the problem
70  * @param nphase0 Number of phases defined within the problem.
71  */
72  void vcs_initSizes(const size_t nspecies0, const size_t nelements, const size_t nphase0);
73 
74  //! Solve an equilibrium problem
75  /*!
76  * This is the main interface routine to the equilibrium solver
77  *
78  * Input:
79  * @param vprob Object containing the equilibrium Problem statement
80  *
81  * @param ifunc Determines the operation to be done: Valid values:
82  * 0 -> Solve a new problem by initializing structures
83  * first. An initial estimate may or may not have
84  * been already determined. This is indicated in the
85  * VCS_PROB structure.
86  * 1 -> The problem has already been initialized and
87  * set up. We call this routine to resolve it
88  * using the problem statement and
89  * solution estimate contained in
90  * the VCS_PROB structure.
91  * 2 -> Don't solve a problem. Destroy all the private
92  * structures.
93  *
94  * @param ipr Printing of results
95  * ipr = 1 -> Print problem statement and final results to
96  * standard output
97  * 0 -> don't report on anything
98  * @param ip1 Printing of intermediate results
99  * IP1 = 1 -> Print intermediate results.
100  *
101  * @param maxit Maximum number of iterations for the algorithm
102  *
103  * Output:
104  *
105  * @return
106  * nonzero value: failure to solve the problem at hand.
107  * zero : success
108  */
109  int vcs(VCS_PROB* vprob, int ifunc, int ipr, int ip1, int maxit);
110 
111  //! Main routine that solves for equilibrium at constant T and P
112  //! using a variant of the VCS method
113  /*!
114  * This is the main routine that solves for equilibrium at constant T and P
115  * using a variant of the VCS method. Nonideal phases can be accommodated
116  * as well.
117  *
118  * Any number of single-species phases and multi-species phases
119  * can be handled by the present version.
120  *
121  * @param print_lvl 1 -> Print results to standard output;
122  * 0 -> don't report on anything
123  * @param printDetails 1 -> Print intermediate results.
124  * @param maxit Maximum number of iterations for the algorithm
125  *
126  * @return
127  * * 0 = Equilibrium Achieved
128  * * 1 = Range space error encountered. The element abundance criteria
129  * are only partially satisfied. Specifically, the first NC= (number
130  * of components) conditions are satisfied. However, the full NE
131  * (number of elements) conditions are not satisfied. The
132  * equilibrium condition is returned.
133  * * -1 = Maximum number of iterations is exceeded. Convergence was
134  * not found.
135  */
136  int vcs_solve_TP(int print_lvl, int printDetails, int maxit);
137 
138  int vcs_PS(VCS_PROB* vprob, int iph, int printLvl, double& feStable);
139 
140  /*!
141  * We make decisions on the initial mole number, and major-minor status
142  * here. We also fix up the total moles in a phase.
143  *
144  * irxn = id of the noncomponent species formation reaction for the
145  * species to be added in.
146  *
147  * The algorithm proceeds to implement these decisions in the previous
148  * position of the species. Then, vcs_switch_pos is called to move the
149  * species into the last active species slot, incrementing the number
150  * of active species at the same time.
151  *
152  * This routine is responsible for the global data manipulation only.
153  */
154  void vcs_reinsert_deleted(size_t kspec);
155 
156  //! Choose the optimum species basis for the calculations
157  /*!
158  * This is done by choosing the species with the largest mole fraction not
159  * currently a linear combination of the previous components. Then,
160  * calculate the stoichiometric coefficient matrix for that basis.
161  *
162  * Rearranges the solution data to put the component data at the
163  * front of the species list.
164  *
165  * Then, calculates m_stoichCoeffRxnMatrix(jcomp,irxn) the formation
166  * reactions for all noncomponent species in the mechanism. Also
167  * calculates DNG(I) and DNL(I), the net mole change for each formation
168  * reaction. Also, initializes IR(I) to the default state.
169  *
170  * @param[in] doJustComponents If true, the m_stoichCoeffRxnMatrix and
171  * m_deltaMolNumPhase are not calculated.
172  *
173  * @param[in] aw Vector of mole fractions which will be used to construct an
174  * optimal basis from.
175  *
176  * @param[in] sa Gram-Schmidt orthog work space (nc in length) sa[j]
177  * @param[in] ss Gram-Schmidt orthog work space (nc in length) ss[j]
178  * @param[in] sm QR matrix work space (nc*ne in length) sm[i+j*ne]
179  * @param[in] test This is a small negative number dependent upon whether
180  * an estimate is supplied or not.
181  * @param[out] usedZeroedSpecies If true, then a species with a zero
182  * concentration was used as a component.
183  * The problem may be converged. Or, the
184  * problem may have a range space error and
185  * may not have a proper solution.
186  * @return Returns VCS_SUCCESS if everything went ok. Returns
187  * VCS_FAILED_CONVERGENCE if there is a problem.
188  *
189  * ### Internal Variables calculated by this routine:
190  *
191  * - #m_numComponents: Number of component species. This routine
192  * calculates the #m_numComponents species. It switches their positions
193  * in the species vector so that they occupy the first #m_numComponents
194  * spots in the species vector.
195  * - #m_stoichCoeffRxnMatrix(jcomp,irxn) Stoichiometric coefficient
196  * matrix for the reaction mechanism expressed in Reduced Canonical
197  * Form. jcomp refers to the component number, and irxn refers to the
198  * irxn_th non-component species.
199  * - #m_deltaMolNumPhase(iphase,irxn): Change in the number of moles in
200  * phase, iphase, due to the noncomponent formation reaction, irxn.
201  * - #m_phaseParticipation(iphase,irxn): This is 1 if the phase, iphase,
202  * participates in the formation reaction, irxn, and zero otherwise.
203  */
204  int vcs_basopt(const bool doJustComponents, double aw[], double sa[], double sm[],
205  double ss[], double test, bool* const usedZeroedSpecies);
206 
207  //! Choose a species to test for the next component
208  /*!
209  * We make the choice based on testing (molNum[i] * spSize[i]) for its
210  * maximum value. Preference for single species phases is also made.
211  *
212  * @param molNum Mole number vector
213  * @param j index into molNum[] that indicates where the search
214  * will start from Previous successful components are
215  * swapped into the front of molNum[].
216  * @param n Length of molNum[]
217  */
218  size_t vcs_basisOptMax(const double* const molNum, const size_t j, const size_t n);
219 
220  //! Evaluate the species category for the indicated species
221  /*!
222  * All evaluations are done using the "old" version of the solution.
223  *
224  * @param kspec Species to be evaluated
225  *
226  * @return Returns the calculated species type
227  */
228  int vcs_species_type(const size_t kspec) const;
229 
230  //! This routine evaluates the species type for all species
231  /*!
232  * - #VCS_SPECIES_MAJOR: Major species
233  * - #VCS_SPECIES_MINOR: Minor species
234  * - #VCS_SPECIES_SMALLMS: The species lies in a multicomponent phase
235  * that exists. Its concentration is currently very low, necessitating
236  * a different method of calculation.
237  * - #VCS_SPECIES_ZEROEDMS: The species lies in a multicomponent phase
238  * which currently doesn't exist. Its concentration is currently zero.
239  * - #VCS_SPECIES_ZEROEDSS: Species lies in a single-species phase which
240  * is currently zeroed out.
241  * - #VCS_SPECIES_DELETED: Species has such a small mole fraction it is
242  * deleted even though its phase may possibly exist. The species is
243  * believed to have such a small mole fraction that it best to throw
244  * the calculation of it out. It will be added back in at the end of
245  * the calculation.
246  * - #VCS_SPECIES_INTERFACIALVOLTAGE: Species refers to an electron in
247  * the metal The unknown is equal to the interfacial voltage drop
248  * across the interface on the SHE (standard hydrogen electrode) scale
249  * (volts).
250  * - #VCS_SPECIES_ZEROEDPHASE: Species lies in a multicomponent phase
251  * that is zeroed atm and will stay deleted due to a choice from a
252  * higher level. These species will formally always have zero mole
253  * numbers in the solution vector.
254  * - #VCS_SPECIES_ACTIVEBUTZERO: The species lies in a multicomponent
255  * phase which currently does exist. Its concentration is currently
256  * identically zero, though the phase exists. Note, this is a temporary
257  * condition that exists at the start of an equilibrium problem. The
258  * species is soon "birthed" or "deleted".
259  * - #VCS_SPECIES_STOICHZERO: The species lies in a multicomponent phase
260  * which currently does exist. Its concentration is currently
261  * identically zero, though the phase exists. This is a permanent
262  * condition due to stoich constraints
263  */
265 
266  //! We calculate the dimensionless chemical potentials of all species
267  //! in a single phase.
268  /*!
269  * Note, for multispecies phases which are currently zeroed out, the
270  * chemical potential is filled out with the standard chemical potential.
271  *
272  * For species in multispecies phases whose concentration is zero, we need
273  * to set the mole fraction to a very low value. Its chemical potential is
274  * then calculated using the #VCS_DELETE_MINORSPECIES_CUTOFF concentration
275  * to keep numbers positive.
276  *
277  * # Formula:
278  *
279  * ## Ideal Mixtures:
280  *
281  * m_feSpecies(I) = m_SSfeSpecies(I) + ln(z(I)) - ln(m_tPhaseMoles[iph])
282  * + m_chargeSpecies[I] * Faraday_dim * m_phasePhi[iphase];
283  *
284  * (This is equivalent to the adding the log of the mole fraction onto
285  * the standard chemical potential. )
286  *
287  * ## Non-Ideal Mixtures:
288  *
289  * ### ActivityConvention = 0: molarity activity formulation
290  *
291  * m_feSpecies(I) = m_SSfeSpecies(I)
292  * + ln(ActCoeff[I] * z(I)) - ln(m_tPhaseMoles[iph])
293  * + m_chargeSpecies[I] * Faraday_dim * m_phasePhi[iphase];
294  *
295  * ( This is equivalent to the adding the log of the mole fraction
296  * multiplied by the activity coefficient onto the standard chemical
297  * potential. )
298  *
299  * ### ActivityConvention = 1: molality activity formulation
300  *
301  * m_feSpecies(I) = m_SSfeSpecies(I)
302  * + ln(ActCoeff[I] * z(I)) - ln(m_tPhaseMoles[iph])
303  * - ln(Mnaught * m_units)
304  * + m_chargeSpecies[I] * Faraday_dim * m_phasePhi[iphase];
305  *
306  * Note: `m_SSfeSpecies(I)` is the molality based standard state. However,
307  * `ActCoeff[I]` is the molar based activity coefficient We have used the
308  * formulas:
309  *
310  * ActCoeff_M[I] = ActCoeff[I] / Xmol[N]
311  *
312  * where `Xmol[N]` is the mole fraction of the solvent and `ActCoeff_M[I]`
313  * is the molality based act coeff.
314  *
315  * Note: This is equivalent to the "normal" molality formulation:
316  *
317  * m_feSpecies(I) = m_SSfeSpecies(I)
318  * + ln(ActCoeff_M[I] * m(I))
319  * + m_chargeSpecies[I] * Faraday_dim * m_phasePhi[iphase]
320  *
321  * where `m[I]` is the molality of the ith solute
322  *
323  * m[I] = Xmol[I] / ( Xmol[N] * Mnaught * m_units)
324  *
325  * `z(I)/tPhMoles_ptr[iph] = Xmol[i]` is the mole fraction of i in the phase.
326  *
327  *
328  * NOTE: As per the discussion in vcs_dfe(), for small species where the
329  * mole fraction is small:
330  *
331  * z(i) < VCS_DELETE_MINORSPECIES_CUTOFF
332  *
333  * The chemical potential is calculated as:
334  *
335  * m_feSpecies(I) = m_SSfeSpecies(I)
336  * + ln(ActCoeff[i](VCS_DELETE_MINORSPECIES_CUTOFF))
337  *
338  * @param[in] iph Phase to be calculated
339  * @param[in] molNum Number of moles of species i (VCS species order)
340  * @param[out] ac Activity coefficients for species in phase (VCS
341  * species order)
342  * @param[out] mu_i Dimensionless chemical potentials for phase
343  * species (VCS species order)
344  */
345  void vcs_chemPotPhase(const int stateCalc, const size_t iph, const double* const molNum,
346  double* const ac, double* const mu_i,
347  const bool do_deleted = false);
348 
349  //! Calculate the dimensionless chemical potentials of all species or
350  //! of certain groups of species, at a fixed temperature and pressure.
351  /*!
352  * We calculate the dimensionless chemical potentials of all species
353  * or certain groups of species here, at a fixed temperature and pressure,
354  * for the input mole vector z[] in the parameter list.
355  * Nondimensionalization is achieved by division by RT.
356  *
357  * Note, for multispecies phases which are currently zeroed out,
358  * the chemical potential is filled out with the standard chemical
359  * potential.
360  *
361  * For species in multispecies phases whose concentration is zero, we need
362  * to set the mole fraction to a very low value. Its chemical potential is
363  * then calculated using the VCS_DELETE_MINORSPECIES_CUTOFF concentration
364  * to keep numbers positive.
365  *
366  * For formulas, see vcs_chemPotPhase().
367  *
368  * Handling of Small Species:
369  * ------------------------------
370  * As per the discussion above, for small species where the mole fraction
371  *
372  * z(i) < VCS_DELETE_MINORSPECIES_CUTOFF
373  *
374  * The chemical potential is calculated as:
375  *
376  * m_feSpecies(I)(I) = m_SSfeSpecies(I) + ln(ActCoeff[i](VCS_DELETE_MINORSPECIES_CUTOFF))
377  *
378  * Species in the following categories are treated as "small species"
379  * - #VCS_SPECIES_DELETED
380  * - #VCS_SPECIES_ACTIVEBUTZERO
381  *
382  * For species in multispecies phases which are currently not active, the
383  * treatment is different. These species are in the following species
384  * categories:
385  * - #VCS_SPECIES_ZEROEDMS
386  * - #VCS_SPECIES_ZEROEDPHASE
387  *
388  * For these species, the `ln( ActCoeff[I] X[I])` term is dropped
389  * altogether. The following equation is used:
390  *
391  * m_feSpecies(I) = m_SSfeSpecies(I)
392  * + Charge[I] * Faraday_dim * phasePhi[iphase];
393  *
394  * Handling of "Species" Representing Interfacial Voltages
395  * ---------------------------------------------------------
396  *
397  * These species have species types of
398  * #VCS_SPECIES_TYPE_INTERFACIALVOLTAGE The chemical potentials for these
399  * "species" refer to electrons in metal electrodes. They have the
400  * following formula
401  *
402  * m_feSpecies(I) = m_SSfeSpecies(I) - F z[I] / RT
403  *
404  * - `F` is Faraday's constant.
405  * - `R` = gas constant
406  * - `T` = temperature
407  * - `V` = potential of the interface = phi_electrode - phi_solution
408  *
409  * For these species, the solution vector unknown, z[I], is V, the phase
410  * voltage, in volts.
411  *
412  * @param ll Determine which group of species gets updated
413  * - `ll = 0`: Calculate for all species
414  * - `ll < 0`: calculate for components and for major non-components
415  * - `ll = 1`: calculate for components and for minor non-components
416  *
417  * @param lbot Restricts the calculation of the chemical potential
418  * to the species between LBOT <= i < LTOP. Usually
419  * LBOT and LTOP will be equal to 0 and MR, respectively.
420  * @param ltop Top value of the loops
421  *
422  * @param stateCalc Determines whether z is old or new or tentative:
423  * - 1: Use the tentative values for the total number of
424  * moles in the phases, i.e., use TG1 instead of TG etc.
425  * - 0: Use the base values of the total number of
426  * moles in each system.
427  *
428  * Also needed:
429  * ff : standard state chemical potentials. These are the
430  * chemical potentials of the standard states at
431  * the same T and P as the solution.
432  * tg : Total Number of moles in the phase.
433  */
434  void vcs_dfe(const int stateCalc, const int ll, const size_t lbot, const size_t ltop);
435 
436  //! Print out a table of chemical potentials
437  /*!
438  * @param stateCalc Determines where to get the mole numbers from.
439  * - VCS_STATECALC_OLD -> from m_molNumSpecies_old
440  * - VCS_STATECALC_NEW -> from m_molNumSpecies_new
441  */
442  void vcs_printSpeciesChemPot(const int stateCalc) const;
443 
444  //! This routine uploads the state of the system into all of the
445  //! vcs_VolumePhase objects in the current problem.
446  /*!
447  * @param stateCalc Determines where to get the mole numbers from.
448  * - VCS_STATECALC_OLD -> from m_molNumSpecies_old
449  * - VCS_STATECALC_NEW -> from m_molNumSpecies_new
450  */
451  void vcs_updateVP(const int stateCalc);
452 
453  //! Utility function that evaluates whether a phase can be popped
454  //! into existence
455  /*!
456  * A phase can be popped iff the stoichiometric coefficients for the
457  * component species, whose concentrations will be lowered during the
458  * process, are positive by at least a small degree.
459  *
460  * If one of the phase species is a zeroed component, then the phase can
461  * be popped if the component increases in mole number as the phase moles
462  * are increased.
463  *
464  * @param iphasePop id of the phase, which is currently zeroed,
465  *
466  * @return Returns true if the phase can come into existence
467  * and false otherwise.
468  */
469  bool vcs_popPhasePossible(const size_t iphasePop) const;
470 
471  //! Determine the list of problems that need to be checked to see if there are any phases pops
472  /*!
473  * This routine evaluates and fills in #phasePopProblemLists_. Need to
474  * work in species that are zeroed by element constraints.
475  *
476  * @return Returns the number of problems that must be checked.
477  */
479 
480  //! Decision as to whether a phase pops back into existence
481  /*!
482  * @param phasePopPhaseIDs Vector containing the phase ids of the phases
483  * that will be popped this step.
484  *
485  * @return returns the phase id of the phase that pops back into
486  * existence. Returns -1 if there are no phases
487  */
488  size_t vcs_popPhaseID(std::vector<size_t> &phasePopPhaseIDs);
489 
490  //! Calculates the deltas of the reactions due to phases popping
491  //! into existence
492  /*!
493  * Updates #m_deltaMolNumSpecies[irxn] : reaction adjustments, where irxn
494  * refers to the irxn'th species formation reaction. This adjustment is
495  * for species irxn + M, where M is the number of components.
496  *
497  * @param iphasePop Phase id of the phase that will come into existence
498  *
499  * @return Returns an int representing the status of the step
500  * - 0 : normal return
501  * - 1 : A single species phase species has been zeroed out
502  * in this routine. The species is a noncomponent
503  * - 2 : Same as one but, the zeroed species is a component.
504  * - 3 : Nothing was done because the phase couldn't be birthed
505  * because a needed component is zero.
506  */
507  int vcs_popPhaseRxnStepSizes(const size_t iphasePop);
508 
509  //! Calculates formation reaction step sizes.
510  /*!
511  * This is equation 6.4-16, p. 143 in Smith and Missen.
512  *
513  * Output
514  * -------
515  * m_deltaMolNumSpecies(irxn) : reaction adjustments, where irxn refers to
516  * the irxn'th species formation reaction.
517  * This adjustment is for species irxn + M,
518  * where M is the number of components.
519  *
520  * Special branching occurs sometimes. This causes the component basis
521  * to be reevaluated
522  *
523  * @param forceComponentCalc integer flagging whether a component
524  * recalculation needs to be carried out.
525  * @param kSpecial species number of phase being zeroed.
526  *
527  * @return Returns an int representing which phase may need to be zeroed
528  */
529  size_t vcs_RxnStepSizes(int& forceComponentCalc, size_t& kSpecial);
530 
531  //! Calculates the total number of moles of species in all phases.
532  /*!
533  * Also updates the variable m_totalMolNum and Reconciles Phase existence
534  * flags with total moles in each phase.
535  */
536  double vcs_tmoles();
537 #ifdef DEBUG_MODE
538  void check_tmoles() const;
539 #endif
540 
541  //! This subroutine calculates reaction free energy changes for
542  //! all noncomponent formation reactions.
543  /*!
544  * Formation reactions are
545  * reactions which create each noncomponent species from the component
546  * species. m_stoichCoeffRxnMatrix(jcomp,irxn) are the stoichiometric
547  * coefficients for these reactions. A stoichiometric coefficient of
548  * one is assumed for species irxn in this reaction.
549  *
550  * @param l
551  * - `L < 0`: Calculate reactions corresponding to major noncomponent
552  * and zeroed species only
553  * - `L = 0`: Do all noncomponent reactions, i, between
554  * 0 <= i < irxnl
555  * - `L > 0`: Calculate reactions corresponding to minor noncomponent
556  * and zeroed species only
557  *
558  * @param doDeleted Do deleted species
559  * @param vcsState Calculate deltaG corresponding to either old or new
560  * free energies
561  * @param alterZeroedPhases boolean indicating whether we should
562  * add in a special section for zeroed phases.
563  *
564  * Note we special case one important issue. If the component has zero
565  * moles, then we do not allow deltaG < 0.0 for formation reactions which
566  * would lead to the loss of more of that same component. This dG < 0.0
567  * condition feeds back into the algorithm in several places, and leads
568  * to a infinite loop in at least one case.
569  */
570  void vcs_deltag(const int l, const bool doDeleted, const int vcsState,
571  const bool alterZeroedPhases = true);
572 
573  void vcs_printDeltaG(const int stateCalc);
574 
575  //! Calculate deltag of formation for all species in a single phase.
576  /*!
577  * Calculate deltag of formation for all species in a single phase. It is
578  * assumed that the fe[] is up to date for all species. However, if the
579  * phase is currently zeroed out, a subproblem is calculated to solve for
580  * AC[i] and pseudo-X[i] for that phase.
581  *
582  * @param iphase phase index of the phase to be calculated
583  * @param doDeleted boolean indicating whether to do deleted
584  * species or not
585  * @param stateCalc integer describing which set of free energies
586  * to use and where to stick the results.
587  * @param alterZeroedPhases boolean indicating whether we should
588  * add in a special section for zeroed phases.
589  *
590  * NOTE: this is currently not used used anywhere.
591  * It may be in the future?
592  */
593  void vcs_deltag_Phase(const size_t iphase, const bool doDeleted,
594  const int stateCalc, const bool alterZeroedPhases = true);
595 
596  //! Swaps the indices for all of the global data for two species, k1
597  //! and k2.
598  /*!
599  * @param ifunc: If true, switch the species data and the noncomponent
600  * reaction data. This must be called for a non-component
601  * species only. If false, switch the species data only.
602  * Typically, we use this option when determining the
603  * component species and at the end of the calculation,
604  * when we want to return unscrambled results. All rxn
605  * data will be out-of-date.
606  *
607  * @param k1 First species index
608  * @param k2 Second species index
609  */
610  void vcs_switch_pos(const bool ifunc, const size_t k1, const size_t k2);
611 
612  //! Birth guess returns the number of moles of a species
613  //! that is coming back to life.
614  /*!
615  * Birth guess returns the number of moles of a species that is coming
616  * back to life. Note, this routine is not applicable if the whole phase
617  * is coming back to life, not just one species in that phase.
618  *
619  * Do a minor alt calculation. But, cap the mole numbers at 1.0E-15. For
620  * SS phases use VCS_DELETE_SPECIES_CUTOFF * 100.
621  *
622  * The routine makes sure the guess doesn't reduce the concentration of a
623  * component by more than 1/3. Note this may mean that the vlaue coming
624  * back from this routine is zero or a very small number.
625  *
626  * @param kspec Species number that is coming back to life
627  * @return Returns the number of kmol that the species should have.
628  */
629  double vcs_birthGuess(const int kspec);
630 
631  //! Routine that independently determines whether a phase should be popped
632  //! under the current conditions.
633  /*
634  * This is the main routine that solves for equilibrium at constant T and
635  * P using a variant of the VCS method. Nonideal phases can be
636  * accommodated as well. Any number of single-species phases and multi-
637  * species phases can be handled by the present version.
638  *
639  * @param print_lvl 1 -> Print results to standard output; -> don't
640  * report on anything
641  * @param printDetails 1 -> Print intermediate results.
642  * @param maxit Maximum number of iterations for the algorithm
643  *
644  * @return
645  * - 0 = Equilibrium Achieved
646  * - 1 = Range space error encountered. The element abundance criteria are
647  * only partially satisfied. Specifically, the first NC= (number of
648  * components) conditions are satisfied. However, the full NE (number of
649  * elements) conditions are not satisfied. The equilibrium condition is
650  * returned.
651  * - -1 = Maximum number of iterations is exceeded. Convergence was not
652  * found.
653  */
654  int vcs_solve_phaseStability(const int iphase, int ifunc, double& funcval, int print_lvl);
655 
656  //! Main program to test whether a deleted phase should be brought
657  //! back into existence
658  /*!
659  * @param iph Phase id of the deleted phase
660  */
661  double vcs_phaseStabilityTest(const size_t iph);
662 
663  //! Solve an equilibrium problem at a particular fixed temperature
664  //! and pressure
665  /*!
666  * The actual problem statement is assumed to be in the structure
667  * already. This is a wrapper around the solve_TP() function. In this
668  * wrapper, we nondimensionalize the system we calculate the standard
669  * state Gibbs free energies of the species, and we decide whether to we
670  * need to use the initial guess algorithm.
671  *
672  * @param ipr = 1 -> Print results to standard output;
673  * 0 -> don't report on anything
674  * @param ip1 = 1 -> Print intermediate results;
675  * 0 -> Dont print any intermediate results
676  * @param maxit Maximum number of iterations for the algorithm
677  * @param T Value of the Temperature (Kelvin)
678  * @param pres Value of the Pressure (units given by m_VCS_UnitsFormat variable
679  *
680  * @return Returns an integer representing the success of the algorithm
681  * * 0 = Equilibrium Achieved
682  * * 1 = Range space error encountered. The element abundance criteria are
683  * only partially satisfied. Specifically, the first NC= (number of
684  * components) conditions are satisfied. However, the full NE (number of
685  * elements) conditions are not satisfied. The equilibrium condition is
686  * returned.
687  * * -1 = Maximum number of iterations is exceeded. Convergence was not
688  * found.
689  */
690  int vcs_TP(int ipr, int ip1, int maxit, double T, double pres);
691 
692  /*!
693  * Evaluate the standard state free energies at the current temperature
694  * and pressure. Ideal gas pressure contribution is added in here.
695  *
696  * @param ipr 1 -> Print results to standard output; 0 -> don't report
697  * on anything
698  * @param ip1 1 -> Print intermediate results; 0 -> don't.
699  * @param Temp Temperature (Kelvin)
700  * @param pres Pressure (Pascal)
701  */
702  int vcs_evalSS_TP(int ipr, int ip1, double Temp, double pres);
703 
704  //! Initialize the chemical potential of single species phases
705  /*!
706  * For single species phases, initialize the chemical potential with the
707  * value of the standard state chemical potential. This value doesn't
708  * change during the calculation
709  */
710  void vcs_fePrep_TP();
711 
712  //! Calculation of the total volume and the partial molar volumes
713  /*!
714  * This function calculates the partial molar volume for all species,
715  * kspec, in the thermo problem at the temperature TKelvin and pressure,
716  * Pres, pres is in atm. And, it calculates the total volume of the
717  * combined system.
718  *
719  * @param[in] tkelvin Temperature in kelvin()
720  * @param[in] pres Pressure in Pascal
721  * @param[in] w w[] is the vector containing the current mole
722  * numbers in units of kmol.
723  * @param[out] volPM[] For species in all phase, the entries are the
724  * partial molar volumes units of M**3 / kmol.
725  * @return The return value is the total volume of
726  * the entire system in units of m**3.
727  */
728  double vcs_VolTotal(const double tkelvin, const double pres,
729  const double w[], double volPM[]);
730 
731  //! This routine is mostly concerned with changing the private data
732  //! to be consistent with what's needed for solution. It is called one
733  //! time for each new problem structure definition.
734  /*!
735  * This routine is always followed by vcs_prep(). Therefore, tasks
736  * that need to be done for every call to vcsc() should be placed in
737  * vcs_prep() and not in this routine.
738  *
739  * The problem structure refers to:
740  *
741  * - the number and identity of the species.
742  * - the formula matrix and thus the number of components.
743  * - the number and identity of the phases.
744  * - the equation of state
745  * - the method and parameters for determining the standard state
746  * - The method and parameters for determining the activity coefficients.
747  *
748  * Tasks:
749  * 1. Fill in the SSPhase[] array.
750  * 2. Check to see if any multispecies phases actually have only one
751  * species in that phase. If true, reassign that phase and species
752  * to be a single-species phase.
753  * 3. Determine the number of components in the problem if not already
754  * done so. During this process the order of the species is changed
755  * in the private data structure. All references to the species
756  * properties must employ the ind[] index vector.
757  *
758  * @param printLvl Print level of the routine
759  * @return VCS_SUCCESS = everything went OK
760  */
761  int vcs_prep_oneTime(int printLvl);
762 
763  //! Prepare the object for solution
764  /*!
765  * This routine is mostly concerned with changing the private data
766  * to be consistent with that needed for solution. It is called for
767  * every invocation of the vcs_solve() except for the cleanup invocation.
768  *
769  * Tasks:
770  * 1. Initialization of arrays to zero.
771  *
772  * @return
773  * VCS_SUCCESS = everything went OK;
774  * VCS_PUB_BAD = There is an irreconcilable difference in the
775  * public data structure from when the problem was
776  * initially set up.
777  */
778  int vcs_prep();
779 
780  //! In this routine, we check for things that will cause the algorithm
781  //! to fail.
782  /*!
783  * We check to see if the problem is well posed. If it is not, we return
784  * false and print out error conditions.
785  *
786  * Current there is one condition. If all the element abundances are
787  * zero, the algorithm will fail.
788  *
789  * @param vprob VCS_PROB pointer to the definition of the equilibrium
790  * problem
791  *
792  * @return If true, the problem is well-posed. If false, the problem
793  * is not well posed.
794  */
795  bool vcs_wellPosed(VCS_PROB* vprob);
796 
797  //! Rearrange the constraint equations represented by the Formula
798  //! Matrix so that the operational ones are in the front
799  /*!
800  * This subroutine handles the rearrangement of the constraint equations
801  * represented by the Formula Matrix. Rearrangement is only necessary when
802  * the number of components is less than the number of elements. For this
803  * case, some constraints can never be satisfied exactly, because the
804  * range space represented by the Formula Matrix of the components can't
805  * span the extra space. These constraints, which are out of the range
806  * space of the component Formula matrix entries, are migrated to the back
807  * of the Formula matrix.
808  *
809  * A prototypical example is an extra element column in FormulaMatrix[],
810  * which is identically zero. For example, let's say that argon is has an
811  * element column in FormulaMatrix[], but no species in the mechanism
812  * actually contains argon. Then, nc < ne. Also, without perturbation of
813  * FormulaMatrix[] vcs_basopt[] would produce a zero pivot because the
814  * matrix would be singular (unless the argon element column was already
815  * the last column of FormulaMatrix[].
816  *
817  * This routine borrows heavily from vcs_basopt's algorithm. It finds nc
818  * constraints which span the range space of the Component Formula matrix,
819  * and assigns them as the first nc components in the formula matrix. This
820  * guarantees that vcs_basopt[] has a nonsingular matrix to invert.
821  *
822  * Other Variables
823  * @param aw Mole fraction work space (ne in length)
824  * @param sa Gram-Schmidt orthog work space (ne in length)
825  * @param sm QR matrix work space (ne*ne in length)
826  * @param ss Gram-Schmidt orthog work space (ne in length)
827  */
828  int vcs_elem_rearrange(double* const aw, double* const sa,
829  double* const sm, double* const ss);
830 
831  //! Swaps the indices for all of the global data for two elements, ipos
832  //! and jpos.
833  /*!
834  * This function knows all of the element information with VCS_SOLVE, and
835  * can therefore switch element positions
836  *
837  * @param ipos first global element index
838  * @param jpos second global element index
839  */
840  void vcs_switch_elem_pos(size_t ipos, size_t jpos);
841 
842  //! Calculates reaction adjustments using a full Hessian approximation
843  /*!
844  * This does what equation 6.4-16, p. 143 in Smith and Missen is supposed
845  * to do. However, a full matrix is formed and then solved via a conjugate
846  * gradient algorithm. No preconditioning is done.
847  *
848  * If special branching is warranted, then the program bails out.
849  *
850  * Output
851  * -------
852  * DS(I) : reaction adjustment, where I refers to the Ith species
853  * Special branching occurs sometimes. This causes the component basis
854  * to be reevaluated
855  * return = 0 : normal return
856  * 1 : A single species phase species has been zeroed out
857  * in this routine. The species is a noncomponent
858  * 2 : Same as one but, the zeroed species is a component.
859  *
860  * Special attention is taken to flag cases where the direction of the
861  * update is contrary to the steepest descent rule. This is an important
862  * attribute of the regular vcs algorithm. We don't want to violate this.
863  *
864  * NOTE: currently this routine is not used.
865  */
866  int vcs_rxn_adj_cg(void);
867 
868  //! Calculates the diagonal contribution to the Hessian due to
869  //! the dependence of the activity coefficients on the mole numbers.
870  /*!
871  * (See framemaker notes, Eqn. 20 - VCS Equations document)
872  *
873  * We allow the diagonal to be increased positively to any degree.
874  * We allow the diagonal to be decreased to 1/3 of the ideal solution
875  * value, but no more -> it must remain positive.
876  *
877  * NOTE: currently this routine is not used
878  */
879  double vcs_Hessian_diag_adj(size_t irxn, double hessianDiag_Ideal);
880 
881  //! Calculates the diagonal contribution to the Hessian due to
882  //! the dependence of the activity coefficients on the mole numbers.
883  /*!
884  * (See framemaker notes, Eqn. 20 - VCS Equations document)
885  *
886  * NOTE: currently this routine is not used
887  */
888  double vcs_Hessian_actCoeff_diag(size_t irxn);
889 
890  //! Recalculate all of the activity coefficients in all of the phases
891  //! based on input mole numbers
892  /*
893  * @param moleSpeciesVCS kmol of species to be used in the update.
894  *
895  * NOTE: This routine needs to be regulated.
896  */
897  void vcs_CalcLnActCoeffJac(const double* const moleSpeciesVCS);
898 
899  //! A line search algorithm is carried out on one reaction
900  /*!
901  * In this routine we carry out a rough line search algorithm to make
902  * sure that the m_deltaGRxn_new doesn't switch signs prematurely.
903  *
904  * @param irxn Reaction number
905  * @param dx_orig Original step length
906  *
907  * @param ANOTE Output character string stating the conclusions of the
908  * line search
909  * @return Returns the optimized step length found by the search
910  */
911  double vcs_line_search(const size_t irxn, const double dx_orig,
912  char* const ANOTE=0);
913 
914  //! Print out a report on the state of the equilibrium problem to
915  //! standard output.
916  /*!
917  * @param iconv Indicator of convergence, to be printed out in the report:
918  * - 0 converged
919  * - 1 range space error
920  * - -1 not converged
921  */
922  int vcs_report(int iconv);
923 
924  //! Switch all species data back to the original order.
925  /*!
926  * This destroys the data based on reaction ordering.
927  */
928  int vcs_rearrange();
929 
930  //! Returns the multiplier for electric charge terms
931  /*
932  * This is basically equal to F/RT
933  *
934  * @param mu_units integer representing the dimensional units system
935  * @param TKelvin double Temperature in Kelvin
936  *
937  * @return Returns the value of F/RT
938  */
939  double vcs_nondim_Farad(int mu_units, double TKelvin) const;
940 
941  //! Returns the multiplier for the nondimensionalization of the equations
942  /*!
943  * This is basically equal to RT
944  *
945  * @param mu_units integer representing the dimensional units system
946  * @param TKelvin double Temperature in Kelvin
947  *
948  * @return Returns the value of RT
949  */
950  double vcs_nondimMult_TP(int mu_units, double TKelvin) const;
951 
952  //! Nondimensionalize the problem data
953  /*!
954  * Nondimensionalize the free energies using the divisor, R * T
955  *
956  * Essentially the internal data can either be in dimensional form
957  * or in nondimensional form. This routine switches the data from
958  * dimensional form into nondimensional form.
959  *
960  * @todo Add a scale factor based on the total mole numbers.
961  * The algorithm contains hard coded numbers based on the
962  * total mole number. If we ever were faced with a problem
963  * with significantly different total kmol numbers than one
964  * the algorithm would have problems.
965  */
966  void vcs_nondim_TP();
967 
968  //! Redimensionalize the problem data
969  /*!
970  * Reddimensionalize the free energies using the multiplier R * T
971  *
972  * Essentially the internal data can either be in dimensional form
973  * or in nondimensional form. This routine switches the data from
974  * nondimensional form into dimensional form.
975  */
976  void vcs_redim_TP();
977 
978  //! Print the string representing the Chemical potential units
979  /*!
980  * This gets printed using plogf()
981  *
982  * @param unitsFormat Integer representing the units system
983  */
984  void vcs_printChemPotUnits(int unitsFormat) const;
985 
986  //! Computes the current elemental abundances vector
987  /*!
988  * Computes the elemental abundances vector, m_elemAbundances[], and stores it
989  * back into the global structure
990  */
991  void vcs_elab();
992 
993  /*!
994  * Checks to see if the element abundances are in compliance. If they are,
995  * then TRUE is returned. If not, FALSE is returned. Note the number of
996  * constraints checked is usually equal to the number of components in the
997  * problem. This routine can check satisfaction of all of the constraints
998  * in the problem, which is equal to ne. However, the solver can't fix
999  * breakage of constraints above nc, because that nc is the range space by
1000  * definition. Satisfaction of extra constraints would have had to occur
1001  * in the problem specification.
1002  *
1003  * The constraints should be broken up into 2 sections. If a constraint
1004  * involves a formula matrix with positive and negative signs, and eaSet =
1005  * 0.0, then you can't expect that the sum will be zero. There may be
1006  * roundoff that inhibits this. However, if the formula matrix is all of
1007  * one sign, then this requires that all species with nonzero entries in
1008  * the formula matrix be identically zero. We put this into the logic
1009  * below.
1010  *
1011  * @param ibound 1: Checks constraints up to the number of elements;
1012  * 0: Checks constraints up to the number of components.
1013  */
1014  bool vcs_elabcheck(int ibound);
1015 
1016  /*!
1017  * Computes the elemental abundances vector for a single phase,
1018  * elemAbundPhase[], and returns it through the argument list. The mole
1019  * numbers of species are taken from the current value in
1020  * m_molNumSpecies_old[].
1021  */
1022  void vcs_elabPhase(size_t iphase, double* const elemAbundPhase);
1023 
1024  /*!
1025  * This subroutine corrects for element abundances. At the end of the
1026  * subroutine, the total moles in all phases are recalculated again,
1027  * because we have changed the number of moles in this routine.
1028  *
1029  * @param aa temporary work vector, length ne*ne
1030  * @param x temporary work vector, length ne
1031  *
1032  * @returns
1033  * - 0 = Nothing of significance happened, Element abundances were and
1034  * still are good.
1035  * - 1 = The solution changed significantly; The element abundances are
1036  * now good.
1037  * - 2 = The solution changed significantly, The element abundances are
1038  * still bad.
1039  * - 3 = The solution changed significantly, The element abundances are
1040  * still bad and a component species got zeroed out.
1041  *
1042  * Internal data to be worked on::
1043  * - ga Current element abundances
1044  * - m_elemAbundancesGoal Required elemental abundances
1045  * - m_molNumSpecies_old Current mole number of species.
1046  * - m_formulaMatrix[][] Formula matrix of the species
1047  * - ne Number of elements
1048  * - nc Number of components.
1049  *
1050  * NOTES:
1051  * This routine is turning out to be very problematic. There are
1052  * lots of special cases and problems with zeroing out species.
1053  *
1054  * Still need to check out when we do loops over nc vs. ne.
1055  *
1056  */
1057  int vcs_elcorr(double aa[], double x[]);
1058 
1059  //! Create an initial estimate of the solution to the thermodynamic
1060  //! equilibrium problem.
1061  /*!
1062  * @return Return value indicates success:
1063  * - 0: successful initial guess
1064  * - -1: Unsuccessful initial guess; the elemental abundances aren't
1065  * satisfied.
1066  */
1067  int vcs_inest_TP();
1068 
1069  //! Estimate the initial mole numbers by constrained linear programming
1070  /*!
1071  * This is done by running each reaction as far forward or backward as
1072  * possible, subject to the constraint that all mole numbers remain non-
1073  * negative. Reactions for which \f$ \Delta \mu^0 \f$ are positive are run
1074  * in reverse, and ones for which it is negative are run in the forward
1075  * direction. The end result is equivalent to solving the linear
1076  * programming problem of minimizing the linear Gibbs function subject to
1077  * the element and non-negativity constraints.
1078  */
1079  int vcs_setMolesLinProg();
1080 
1081  //! Calculate the total dimensionless Gibbs free energy
1082  /*!
1083  * Inert species are handled as if they had a standard free energy of
1084  * zero. Note, for this algorithm this function should be monotonically
1085  * decreasing.
1086  */
1087  double vcs_Total_Gibbs(double* w, double* fe, double* tPhMoles);
1088 
1089  //! Calculate the total dimensionless Gibbs free energy of a single phase
1090  /*!
1091  * Inert species are handled as if they had a standard free energy of
1092  * zero and if they obeyed ideal solution/gas theory.
1093  *
1094  * @param iphase ID of the phase
1095  * @param w Species mole number vector for all species
1096  * @param fe vector of partial molar free energies of all of the
1097  * species
1098  */
1099  double vcs_GibbsPhase(size_t iphase, const double* const w,
1100  const double* const fe);
1101 
1102  //! Transfer the results of the equilibrium calculation back to VCS_PROB
1103  /*!
1104  * The VCS_PROB structure is returned to the user.
1105  *
1106  * @param pub Pointer to VCS_PROB object that will get the results of the
1107  * equilibrium calculation transfered to it.
1108  */
1109  int vcs_prob_update(VCS_PROB* pub);
1110 
1111  //! Fully specify the problem to be solved using VCS_PROB
1112  /*!
1113  * Use the contents of the VCS_PROB to specify the contents of the
1114  * private data, VCS_SOLVE.
1115  *
1116  * @param pub Pointer to VCS_PROB that will be used to
1117  * initialize the current equilibrium problem
1118  */
1119  int vcs_prob_specifyFully(const VCS_PROB* pub);
1120 
1121  //! Specify the problem to be solved using VCS_PROB, incrementally
1122  /*!
1123  * Use the contents of the VCS_PROB to specify the contents of the
1124  * private data, VCS_SOLVE.
1125  *
1126  * It's assumed we are solving the same problem.
1127  *
1128  * @param pub Pointer to VCS_PROB that will be used to
1129  * initialize the current equilibrium problem
1130  */
1131  int vcs_prob_specify(const VCS_PROB* pub);
1132 
1133 private:
1134 
1135  //! Zero out the concentration of a species.
1136  /*!
1137  * Make sure to conserveelements and keep track of the total moles in all
1138  * phases.
1139  * - w[]
1140  * - m_tPhaseMoles_old[]
1141  *
1142  * @param kspec Species index
1143  *
1144  * @return:
1145  * 1: succeeded
1146  * 0: failed.
1147  */
1148  int vcs_zero_species(const size_t kspec);
1149 
1150  //! Change a single species from active to inactive status
1151  /*!
1152  * Rearrange data when species is added or removed. The Lth species is
1153  * moved to the back of the species vector. The back of the species
1154  * vector is indicated by the value of MR, the current number of
1155  * active species in the mechanism.
1156  *
1157  * @param kspec Species Index
1158  * @return
1159  * Returns 0 unless.
1160  * The return is 1 when the current number of
1161  * noncomponent species is equal to zero. A recheck of deleted species
1162  * is carried out in the main code.
1163  */
1164  int vcs_delete_species(const size_t kspec);
1165 
1166  //! This routine handles the bookkeeping involved with the
1167  //! deletion of multiphase phases from the problem.
1168  /*!
1169  * When they are deleted, all of their species become active
1170  * species, even though their mole numbers are set to zero.
1171  * The routine does not make the decision to eliminate multiphases.
1172  *
1173  * Note, species in phases with zero mole numbers are still
1174  * considered active. Whether the phase pops back into
1175  * existence or not is checked as part of the main iteration
1176  * loop.
1177  *
1178  * @param iph Phase to be deleted
1179  *
1180  * @return Returns whether the operation was successful or not
1181  */
1182  bool vcs_delete_multiphase(const size_t iph);
1183 
1184  //! Change the concentration of a species by delta moles.
1185  /*!
1186  * Make sure to conserve elements and keep track of the total kmoles in
1187  * all phases.
1188  *
1189  * @param kspec The species index
1190  * @param delta_ptr pointer to the delta for the species. This may
1191  * change during the calculation
1192  *
1193  * @return
1194  * 1: succeeded without change of dx
1195  * 0: Had to adjust dx, perhaps to zero, in order to do the delta.
1196  */
1197  int delta_species(const size_t kspec, double* const delta_ptr);
1198 
1199  //! Provide an estimate for the deleted species in phases that
1200  //! are not zeroed out
1201  /*!
1202  * Try to add back in all deleted species. An estimate of the kmol numbers
1203  * are obtained and the species is added back into the equation system,
1204  * into the old state vector.
1205  *
1206  * This routine is called at the end of the calculation, just before
1207  * returning to the user.
1208  */
1209  size_t vcs_add_all_deleted();
1210 
1211  //! Recheck deleted species in multispecies phases.
1212  /*!
1213  * We are checking the equation:
1214  *
1215  * sum_u = sum_j_comp [ sigma_i_j * u_j ]
1216  * = u_i_O + log((AC_i * W_i)/m_tPhaseMoles_old)
1217  *
1218  * by first evaluating:
1219  *
1220  * DG_i_O = u_i_O - sum_u.
1221  *
1222  * Then, if TL is zero, the phase pops into existence if DG_i_O < 0.
1223  * Also, if the phase exists, then we check to see if the species
1224  * can have a mole number larger than VCS_DELETE_SPECIES_CUTOFF
1225  * (default value = 1.0E-32).
1226  *
1227  */
1228  int vcs_recheck_deleted();
1229 
1230  //! Recheck deletion condition for multispecies phases.
1231  /*!
1232  * We assume here that DG_i_0 has been calculated for deleted species correctly
1233  *
1234  * m_feSpecies(I) = m_SSfeSpecies(I)
1235  * + ln(ActCoeff[I])
1236  * - ln(Mnaught * m_units)
1237  * + m_chargeSpecies[I] * Faraday_dim * m_phasePhi[iphase];
1238  *
1239  * sum_u = sum_j_comp [ sigma_i_j * u_j ]
1240  * = u_i_O + log((AC_i * W_i)/m_tPhaseMoles_old)
1241  *
1242  * DG_i_0 = m_feSpecies(I) - sum_m{ a_i_m DG_m }
1243  *
1244  * by first evaluating:
1245  *
1246  * DG_i_O = u_i_O - sum_u.
1247  *
1248  * Then, the phase pops into existence iff
1249  *
1250  * phaseDG = 1.0 - sum_i{exp(-DG_i_O)} < 0.0
1251  *
1252  * This formula works for both single species phases and for multispecies
1253  * phases. It's an overkill for single species phases.
1254  *
1255  * @param iphase Phase index number
1256  *
1257  * @return Returns true if the phase is currently deleted
1258  * but should be reinstated. Returns false otherwise.
1259  *
1260  * NOTE: this routine is currently not used in the code, and
1261  * contains some basic changes that are incompatible.
1262  *
1263  * assumptions:
1264  * 1. Vphase Existence is up to date
1265  * 2. Vphase->IndSpecies is up to date
1266  * 3. m_deltaGRxn_old[irxn] is up to date
1267  */
1268  bool recheck_deleted_phase(const int iphase);
1269 
1270  //! Minor species alternative calculation
1271  /*!
1272  * This is based upon the following approximation: The mole fraction
1273  * changes due to these reactions don't affect the mole numbers of the
1274  * component species. Therefore the following approximation is valid for
1275  * a small component of an ideal phase:
1276  *
1277  * 0 = m_deltaGRxn_old(I) + log(molNum_new(I)/molNum_old(I))
1278  *
1279  * `m_deltaGRxn_old` contains the contribution from
1280  *
1281  * m_feSpecies_old(I) =
1282  * m_SSfeSpecies(I) +
1283  * log(ActCoeff[i] * molNum_old(I) / m_tPhaseMoles_old(iph))
1284  * Thus,
1285  *
1286  * molNum_new(I)= molNum_old(I) * EXP(-m_deltaGRxn_old(I))
1287  *
1288  * Most of this section is mainly restricting the update to reasonable
1289  * values. We restrict the update a factor of 1.0E10 up and 1.0E-10 down
1290  * because we run into trouble with the addition operator due to roundoff
1291  * if we go larger than ~1.0E15. Roundoff will then sometimes produce
1292  * zero mole fractions.
1293  *
1294  * Note: This routine was generalized to incorporate nonideal phases and
1295  * phases on the molality basis
1296  *
1297  * @param[in] kspec The current species and corresponding formation
1298  * reaction number.
1299  * @param[in] irxn The current species and corresponding formation
1300  * reaction number.
1301  * @param[out] do_delete: BOOLEAN which if true on return, then we
1302  * branch to the section that deletes a species
1303  * from the current set of active species.
1304  */
1305  double vcs_minor_alt_calc(size_t kspec, size_t irxn, bool* do_delete,
1306  char* ANOTE=0) const;
1307 
1308  //! This routine optimizes the minimization of the total Gibbs free energy
1309  //! by making sure the slope of the following functional stays negative:
1310  /*!
1311  * The slope of the following functional is equivalent to the slope
1312  * of the total Gibbs free energy of the system:
1313  *
1314  * d_Gibbs/ds = sum_k( m_deltaGRxn * m_deltaMolNumSpecies[k] )
1315  *
1316  * along the current direction m_deltaMolNumSpecies[], by choosing a value, al: (0<al<1)
1317  * such that the a parabola approximation to Gibbs(al) fit to the
1318  * end points al = 0 and al = 1 is minimized.
1319  * - s1 = slope of Gibbs function at al = 0, which is the previous
1320  * solution = d(Gibbs)/d(al).
1321  * - s2 = slope of Gibbs function at al = 1, which is the current
1322  * solution = d(Gibbs)/d(al).
1323  *
1324  * Only if there has been an inflection point (i.e., s1 < 0 and s2 > 0),
1325  * does this code section kick in. It finds the point on the parabola
1326  * where the slope is equal to zero.
1327  */
1328  bool vcs_globStepDamp();
1329 
1330  //! Calculate the norm of a deltaGibbs free energy vector
1331  /*!
1332  * Positive DG for species which don't exist are ignored.
1333  *
1334  * @param dg Vector of local delta G's.
1335  */
1336  double l2normdg(double dg[]) const;
1337 
1338 #ifdef DEBUG_MODE
1339 
1340  //! Print out and check the elemental abundance vector
1341  void prneav() const;
1342 #endif
1343 
1344  void checkDelta1(double* const ds, double* const delTPhMoles, size_t kspec);
1345 
1346  //! Estimate equilibrium compositions
1347  /*!
1348  * Algorithm covered in a section of Smith and Missen's Book.
1349  *
1350  * Linear programming module is based on using dbolm.
1351  *
1352  * @param aw aw[i[ Mole fraction work space (ne in length)
1353  * @param sa sa[j] = Gram-Schmidt orthog work space (ne in length)
1354  * @param sm sm[i+j*ne] = QR matrix work space (ne*ne in length)
1355  * @param ss ss[j] = Gram-Schmidt orthog work space (ne in length)
1356  * @param test This is a small negative number.
1357  */
1358  void vcs_inest(double* const aw, double* const sa, double* const sm,
1359  double* const ss, double test);
1360 
1361  //! Calculate the status of single species phases.
1362  void vcs_SSPhase(void);
1363 
1364  //! This function recalculates the deltaG for reaction, irxn
1365  /*!
1366  * This function recalculates the deltaG for reaction irxn, given the mole
1367  * numbers in molNum. It uses the temporary space mu_i, to hold the
1368  * recalculated chemical potentials. It only recalculates the chemical
1369  * potentials for species in phases which participate in the irxn
1370  * reaction. This function is used by the vcs_line_search algorithm() and
1371  * should not be used widely due to the unknown state it leaves the
1372  * system.
1373  *
1374  * @param[in] irxn Reaction number
1375  * @param[in] molNum Current mole numbers of species to be used as input
1376  * to the calculation (units = kmol), (length = totalNuMSpecies)
1377  *
1378  * @param[out] ac Activity coefficients (length = totalNumSpecies) Note
1379  * this is only partially formed. Only species in phases that
1380  * participate in the reaction will be updated
1381  * @param[out] mu_i dimensionless chemical potentials (length
1382  * totalNumSpecies) Note this is only partially formed. Only
1383  * species in phases that participate in the reaction will be
1384  * updated
1385  *
1386  * @return Returns the dimensionless deltaG of the reaction
1387  */
1388  double deltaG_Recalc_Rxn(const int stateCalc,
1389  const size_t irxn, const double* const molNum,
1390  double* const ac, double* const mu_i);
1391 
1392  //! Delete memory that isn't just resizable STL containers
1393  /*!
1394  * This gets called by the destructor or by InitSizes().
1395  */
1396  void vcs_delete_memory();
1397 
1398  //! Initialize the internal counters
1399  /*!
1400  * Initialize the internal counters containing the subroutine call
1401  * values and times spent in the subroutines.
1402  *
1403  * ifunc = 0 Initialize only those counters appropriate for the top of
1404  * vcs_solve_TP().
1405  * = 1 Initialize all counters.
1406  */
1407  void vcs_counters_init(int ifunc);
1408 
1409  //! Create a report on the plog file containing timing and its information
1410  /*!
1411  * @param timing_print_lvl If 0, just report the iteration count.
1412  * If larger than zero, report the timing information
1413  */
1414  void vcs_TCounters_report(int timing_print_lvl = 1);
1415 
1416  void vcs_setFlagsVolPhases(const bool upToDate, const int stateCalc);
1417 
1418  void vcs_setFlagsVolPhase(const size_t iph, const bool upToDate, const int stateCalc);
1419 
1420  //! Update all underlying vcs_VolPhase objects
1421  /*!
1422  * Update the mole numbers and the phase voltages of all phases in the
1423  * vcs problem
1424  *
1425  * @param stateCalc Location of the update (either VCS_STATECALC_NEW or
1426  * VCS_STATECALC_OLD).
1427  */
1428  void vcs_updateMolNumVolPhases(const int stateCalc);
1429 
1430  // Helper functions used internally by vcs_solve_TP
1431  int solve_tp_component_calc(bool& allMinorZeroedSpecies);
1432  void solve_tp_inner(size_t& iti, size_t& it1, bool& uptodate_minors,
1433  bool& allMinorZeroedSpecies, int& forceComponentCalc,
1434  int& stage, bool printDetails, char* ANOTE);
1435  void solve_tp_equilib_check(bool& allMinorZeroedSpecies, bool& uptodate_minors,
1436  bool& giveUpOnElemAbund, int& solveFail,
1437  size_t& iti, size_t& it1, int maxit,
1438  int& stage, bool& lec);
1439  void solve_tp_elem_abund_check(size_t& iti, int& stage, bool& lec,
1440  bool& giveUpOnElemAbund,
1441  int& finalElemAbundAttempts,
1442  int& rangeErrorFound);
1443 
1444  // data used by vcs_solve_TP and it's helper functions
1445  std::vector<double> m_sm;
1446  std::vector<double> m_ss;
1447  std::vector<double> m_sa;
1448  std::vector<double> m_aw;
1449  std::vector<double> m_wx;
1450 
1451 public:
1452  //! Calculate the rank of a matrix and return the rows and columns that
1453  //! will generate an independent basis for that rank
1454  /*
1455  * Choose the optimum component species basis for the calculations,
1456  * finding the rank and set of linearly independent rows for that
1457  * calculation. Then find the set of linearly indepedent element columns
1458  * that can support that rank. This is done by taking the transpose of the
1459  * matrix and redoing the same calculation. (there may be a better way to
1460  * do this. I don't know.)
1461  *
1462  * @param[in] awtmp Vector of mole numbers which will be used to
1463  * construct a ranking for how to pick the basis species. This is
1464  * largely ignored here.
1465  * @param[in] numSpecies Number of species. This is the number of rows in
1466  * the matrix.
1467  * @param[in] matrix This is the formula matrix. Nominally, the rows are
1468  * species, while the columns are element compositions. However,
1469  * this routine is totally general, so that the rows and columns
1470  * can be anything.
1471  * @param[in] numElemConstraints Number of element constraints
1472  *
1473  * @param[out] usedZeroedSpecies If true, then a species with a zero
1474  * concentration was used as a component.
1475  * @param[out] compRes Vector of rows which are linearly independent.
1476  * (these are the components)
1477  * @param[out] elemComp Vector of columns which are linearly independent
1478  * (These are the actionable element constraints).
1479  *
1480  * @return Returns number of components. This is the rank of the matrix
1481  * @deprecated To be removed after Cantera 2.2.
1482  */
1483  int vcs_rank(const double* awtmp, size_t numSpecies, const double* matrix, size_t numElemConstraints,
1484  std::vector<size_t> &compRes, std::vector<size_t> &elemComp, int* const usedZeroedSpecies) const;
1485 
1486  //! value of the number of species used to malloc data structures
1487  size_t NSPECIES0;
1488 
1489  //! value of the number of phases used to malloc data structures
1490  size_t NPHASE0;
1491 
1492  //! Total number of species in the problems
1494 
1495  //! Number of element constraints in the problem
1496  /*!
1497  * This is typically equal to the number of elements in the problem
1498  */
1500 
1501  //! Number of components calculated for the problem
1503 
1504  //! Total number of non-component species in the problem
1505  size_t m_numRxnTot;
1506 
1507  //! Current number of species in the problems
1508  /*!
1509  * Species can be deleted if they aren't
1510  * stable under the current conditions
1511  */
1513 
1514  //! Current number of non-component species in the problem
1515  /*!
1516  * Species can be deleted if they aren't
1517  * stable under the current conditions
1518  */
1519  size_t m_numRxnRdc;
1520 
1521  //! Number of active species which are currently either treated as
1522  //! minor species
1524 
1525  //! Number of Phases in the problem
1526  size_t m_numPhases;
1527 
1528  //! Formula matrix for the problem
1529  /*!
1530  * FormulaMatrix(kspec,j) = Number of elements, j, in the kspec species
1531  *
1532  * Both element and species indices are swapped.
1533  */
1535 
1536  //! Stoichiometric coefficient matrix for the reaction mechanism expressed in Reduced Canonical Form.
1537  /*!
1538  * This is the stoichiometric coefficient matrix for the
1539  * reaction which forms species kspec from the component species. A
1540  * stoichiometric coefficient of one is assumed for the species kspec in this mechanism.
1541  *
1542  * NOTE: kspec = irxn + m_numComponents
1543  *
1544  * m_stoichCoeffRxnMatrix(j,irxn) :
1545  * j refers to the component number, and irxn refers to the irxn_th non-component species.
1546  * The stoichiometric coefficients multilplied by the Formula coefficients of the
1547  * component species add up to the negative value of the number of elements in
1548  * the species kspec.
1549  *
1550  * size = nelements0 x nspecies0
1551  */
1553 
1554  //! Absolute size of the stoichiometric coefficients
1555  /*!
1556  * scSize[irxn] = abs(Size) of the stoichiometric
1557  * coefficients. These are used to determine
1558  * whether a given species should be
1559  * handled by the alt_min treatment or
1560  * should be handled as a major species.
1561  */
1562  std::vector<double> m_scSize;
1563 
1564  //! total size of the species
1565  /*!
1566  * This is used as a multiplier to the mole number in figuring out which
1567  * species should be components.
1568  */
1569  std::vector<double> m_spSize;
1570 
1571  //! Standard state chemical potentials for species K at the current
1572  //! temperature and pressure.
1573  /*!
1574  * The first NC entries are for components. The following NR entries are
1575  * for the current non-component species in the mechanism.
1576  */
1577  std::vector<double> m_SSfeSpecies;
1578 
1579  //! Free energy vector from the start of the current iteration
1580  /*!
1581  * The free energies are saved at the start of the current iteration.
1582  * Length = number of species
1583  */
1584  std::vector<double> m_feSpecies_old;
1585 
1586  //! Dimensionless new free energy for all the species in the mechanism
1587  //! at the new tentatite T, P, and mole numbers.
1588  /*!
1589  * The first NC entries are for components. The following
1590  * NR entries are for the current non-component species in the mechanism.
1591  * Length = number of species
1592  */
1593  std::vector<double> m_feSpecies_new;
1594 
1595  //! Setting for whether to do an initial estimate
1596  /*!
1597  * Initial estimate:
1598  * * 0 Do not estimate the solution at all. Use the
1599  * * supplied mole numbers as is.
1600  * * 1 Only do an estimate if the element abundances
1601  * * aren't satisfied.
1602  * * -1 Force an estimate of the soln. Throw out the input
1603  * * mole numbers.
1604  */
1606 
1607  //! Total moles of the species
1608  /*!
1609  * Total number of moles of the kth species.
1610  * Length = Total number of species = m
1611  */
1612  std::vector<double> m_molNumSpecies_old;
1613 
1614  //! Specifies the species unknown type
1615  /*!
1616  * There are two types. One is the straightforward species, with the mole
1617  * number w[k], as the unknown. The second is the an interfacial voltage
1618  * where w[k] refers to the interfacial voltage in volts.
1619  *
1620  * These species types correspond to metallic electrons corresponding to
1621  * electrodes. The voltage and other interfacial conditions sets up an
1622  * interfacial current, which is set to zero in this initial treatment.
1623  * Later we may have non-zero interfacial currents.
1624  */
1625  std::vector<int> m_speciesUnknownType;
1626 
1627  //! Change in the number of moles of phase, iphase, due to the
1628  //! noncomponent formation reaction, irxn, for species, k:
1629  /*!
1630  * m_deltaMolNumPhase(iphase,irxn) = k = nc + irxn
1631  */
1633 
1634  //! This is 1 if the phase, iphase, participates in the formation reaction
1635  //! irxn, and zero otherwise. PhaseParticipation(iphase,irxn)
1637 
1638  //! electric potential of the iph phase
1639  std::vector<double> m_phasePhi;
1640 
1641  //! Tentative value of the mole number vector. It's also used to store the
1642  //! mole fraction vector.
1643  //std::vector<double> wt;
1644  std::vector<double> m_molNumSpecies_new;
1645 
1646  //! Delta G(irxn) for the noncomponent species in the mechanism.
1647  /*!
1648  * Computed by the subroutine deltaG. m_deltaGRxn is the free
1649  * energy change for the reaction which forms species K from the
1650  * component species. This vector has length equal to the number
1651  * of noncomponent species in the mechanism. It starts with
1652  * the first current noncomponent species in the mechanism.
1653  */
1654  std::vector<double> m_deltaGRxn_new;
1655 
1656  //! Last deltag[irxn] from the previous step
1657  std::vector<double> m_deltaGRxn_old;
1658 
1659  //! Last deltag[irxn] from the previous step with additions for
1660  //! possible births of zeroed phases.
1661  std::vector<double> m_deltaGRxn_Deficient;
1662 
1663  //! Temporary vector of Rxn DeltaG's
1664  /*!
1665  * This is used from time to time, for printing purposes
1666  */
1667  std::vector<double> m_deltaGRxn_tmp;
1668 
1669  //! Reaction Adjustments for each species during the current step
1670  /*!
1671  * delta Moles for each species during the current step.
1672  * Length = number of species
1673  */
1674  std::vector<double> m_deltaMolNumSpecies;
1675 
1676  //! Element abundances vector
1677  /*!
1678  * Vector of moles of each element actually in the solution
1679  * vector. Except for certain parts of the algorithm,
1680  * this is a constant.
1681  * Note other constraint conditions are added to this vector.
1682  * This is input from the input file and
1683  * is considered a constant from thereon.
1684  * units = kmoles
1685  */
1686  std::vector<double> m_elemAbundances;
1687 
1688  //! Element abundances vector Goals
1689  /*!
1690  * Vector of moles of each element that are the goals of the simulation.
1691  * This is a constant in the problem. Note other constraint conditions
1692  * are added to this vector. This is input from the input file and is
1693  * considered a constant from thereon. units = kmoles
1694  */
1695  std::vector<double> m_elemAbundancesGoal;
1696 
1697  //! Total number of kmoles in all phases
1698  /*!
1699  * This number includes the inerts.
1700  * -> Don't use this except for scaling purposes
1701  */
1703 
1704  //! Total kmols of species in each phase
1705  /*!
1706  * This contains the total number of moles of species in each phase
1707  *
1708  * Length = number of phases
1709  */
1710  std::vector<double> m_tPhaseMoles_old;
1711 
1712  //! total kmols of species in each phase in the tentative soln vector
1713  /*!
1714  * This contains the total number of moles of species in each phase
1715  * in the tentative solution vector
1716  *
1717  * Length = number of phases
1718  */
1719  std::vector<double> m_tPhaseMoles_new;
1720 
1721  //! Temporary vector of length NPhase
1722  mutable std::vector<double> m_TmpPhase;
1723 
1724  //! Temporary vector of length NPhase
1725  mutable std::vector<double> m_TmpPhase2;
1726 
1727  //! Change in the total moles in each phase
1728  /*!
1729  * Length number of phases.
1730  */
1731  std::vector<double> m_deltaPhaseMoles;
1732 
1733  //! Temperature (Kelvin)
1735 
1736  //! Pressure (units are determined by m_VCS_UnitsFormat
1737  /*!
1738  * | Values | units |
1739  * | ------ | -----
1740  * | -1: | atm |
1741  * | 0: | atm |
1742  * | 1: | atm |
1743  * | 2: | atm |
1744  * | 3: | Pa |
1745  *
1746  * Units being changed to Pa
1747  */
1749 
1750  //! Total kmoles of inert to add to each phase
1751  /*!
1752  * TPhInertMoles[iph] = Total kmoles of inert to add to each phase
1753  * length = number of phases
1754  */
1755  std::vector<double> TPhInertMoles;
1756 
1757  //! Tolerance requirement for major species
1758  double m_tolmaj;
1759 
1760  //! Tolerance requirements for minor species
1761  double m_tolmin;
1762 
1763  //! Below this, major species aren't refined any more
1764  double m_tolmaj2;
1765 
1766  //! Below this, minor species aren't refined any more
1767  double m_tolmin2;
1768 
1769  //! Index vector that keeps track of the species vector rearrangement
1770  /*!
1771  * At the end of each run, the species vector and associated data gets put back
1772  * in the original order.
1773  *
1774  * Example
1775  *
1776  * k = m_speciesMapIndex[kspec]
1777  *
1778  * kspec = current order in the vcs_solve object
1779  * k = original order in the vcs_prob object and in the MultiPhase object
1780  */
1781  std::vector<size_t> m_speciesMapIndex;
1782 
1783  //! Index that keeps track of the index of the species within the local
1784  //! phase
1785  /*!
1786  * This returns the local index of the species within the phase. Its argument
1787  * is the global species index within the VCS problem.
1788  *
1789  * k = m_speciesLocalPhaseIndex[kspec]
1790  *
1791  * k varies between 0 and the nSpecies in the phase
1792  *
1793  * Length = number of species
1794  */
1795  std::vector<size_t> m_speciesLocalPhaseIndex;
1796 
1797  //! Index vector that keeps track of the rearrangement of the elements
1798  /*!
1799  * At the end of each run, the element vector and associated data gets
1800  * put back in the original order.
1801  *
1802  * Example
1803  *
1804  * e = m_elementMapIndex[eNum]
1805  * eNum = current order in the vcs_solve object
1806  * e = original order in the vcs_prob object and in the MultiPhase object
1807  */
1808  std::vector<size_t> m_elementMapIndex;
1809 
1810  //! Mapping between the species index for noncomponent species and the
1811  //! full species index.
1812  /*!
1813  * ir[irxn] = Mapping between the reaction index for noncomponent
1814  * formation reaction of a species and the full species
1815  * index.
1816  *
1817  * Initially set to a value of K = NC + I This vector has length equal to
1818  * number of noncomponent species in the mechanism. It starts with the
1819  * first current noncomponent species in the mechanism. kspec = ir[irxn]
1820  */
1821  std::vector<size_t> m_indexRxnToSpecies;
1822 
1823  //! Major -Minor status vector for the species in the problem
1824  /*!
1825  * The index for this is species. The reaction that this is referring to
1826  * is `kspec = irxn + m_numComponents`. For possible values and their
1827  * meanings, see vcs_evaluate_speciesType().
1828  */
1829  std::vector<int> m_speciesStatus;
1830 
1831  //! Mapping from the species number to the phase number
1832  std::vector<size_t> m_phaseID;
1833 
1834  //! Boolean indicating whether a species belongs to a single-species phase
1835  // vector<bool> can't be used here because it doesn't work with std::swap
1836  std::vector<char> m_SSPhase;
1837 
1838  //! Species string name for the kth species
1839  std::vector<std::string> m_speciesName;
1840 
1841  //! Vector of strings containing the element names
1842  /*!
1843  * ElName[j] = String containing element names
1844  */
1845  std::vector<std::string> m_elementName;
1846 
1847  //! Type of the element constraint
1848  /*!
1849  * * 0 - #VCS_ELEM_TYPE_ABSPOS Normal element that is positive or zero in
1850  * all species.
1851  * * 1 - #VCS_ELEM_TYPE_ELECTRONCHARGE element dof that corresponds to
1852  * the electronic charge DOF.
1853  * * 2 - #VCS_ELEM_TYPE_CHARGENEUTRALITY element dof that corresponds to
1854  * a required charge neutrality constraint on the phase. The element
1855  * abundance is always exactly zero.
1856  * * 3 - #VCS_ELEM_TYPE_OTHERCONSTRAINT Other constraint which may mean
1857  * that a species has neg 0 or pos value of that constraint (other than
1858  * charge)
1859  */
1860  std::vector<int> m_elType;
1861 
1862  //! Specifies whether an element constraint is active
1863  /*!
1864  * The default is true. Length = nelements
1865  */
1866  std::vector<int> m_elementActive;
1867 
1868  //! Array of Phase Structures. Length = number of phases.
1869  std::vector<vcs_VolPhase*> m_VolPhaseList;
1870 
1871  //! String containing the title of the run
1872  std::string m_title;
1873 
1874  //! This specifies the current state of units for the Gibbs free energy
1875  //! properties in the program.
1876  /*!
1877  * The default is to have this unitless
1878  */
1880 
1881  //! Multiplier for the mole numbers within the nondimensionless formulation
1882  /*!
1883  * All numbers within the main routine are on an absolute basis. This
1884  * presents some problems wrt very large and very small mole numbers.
1885  * We get around this by using a multiplier coming into and coming
1886  * out of the equilibrium routines
1887  */
1889 
1890  //! specifies the activity convention of the phase containing the species
1891  /*!
1892  * * 0 = molar based
1893  * * 1 = molality based
1894  *
1895  * length = number of species
1896  */
1897  std::vector<int> m_actConventionSpecies;
1898 
1899  //! specifies the activity convention of the phase.
1900  /*!
1901  * * 0 = molar based
1902  * * 1 = molality based
1903  *
1904  * length = number of phases
1905  */
1906  std::vector<int> m_phaseActConvention;
1907 
1908  //! specifies the ln(Mnaught) used to calculate the chemical potentials
1909  /*!
1910  * For molar based activity conventions this will be equal to 0.0.
1911  * length = number of species.
1912  */
1913  std::vector<double> m_lnMnaughtSpecies;
1914 
1915  //! Molar-based Activity Coefficients for Species.
1916  //! Length = number of species
1917  std::vector<double> m_actCoeffSpecies_new;
1918 
1919  //! Molar-based Activity Coefficients for Species based on old mole numbers
1920  /*!
1921  * These activity coefficients are based on the m_molNumSpecies_old
1922  * values Molar based activity coeffients. Length = number of species
1923  */
1924  std::vector<double> m_actCoeffSpecies_old;
1925 
1926  //! Change in the log of the activity coefficient with respect to the mole number
1927  //! multiplied by the phase mole number
1928  /*!
1929  * size = nspecies x nspecies
1930  *
1931  * This is a temporary array that gets regenerated every time it's
1932  * needed. It is not swapped wrt species.
1933  */
1935 
1936  //! Molecular weight of each species
1937  /*!
1938  * units = kg/kmol. length = number of species.
1939  *
1940  * note: this is a candidate for removal. I don't think we use it.
1941  */
1942  std::vector<double> m_wtSpecies;
1943 
1944  //! Charge of each species. Length = number of species.
1945  std::vector<double> m_chargeSpecies;
1946 
1947  std::vector<std::vector<size_t> > phasePopProblemLists_;
1948 
1949  //! Vector of pointers to thermostructures which identify the model
1950  //! and parameters for evaluating the thermodynamic functions for that
1951  //! particular species.
1952  /*!
1953  * SpeciesThermo[k] pointer to the thermo information for the kth species
1954  */
1955  std::vector<VCS_SPECIES_THERMO*> m_speciesThermoList;
1956 
1957  //! Choice of Hessians
1958  /*!
1959  * If this is true, then we will use a better approximation to the
1960  * Hessian based on Jacobian of the ln(ActCoeff) with respect to mole
1961  * numbers
1962  */
1964 
1965  //! Total volume of all phases. Units are m^3
1966  double m_totalVol;
1967 
1968  //! Partial molar volumes of the species
1969  /*!
1970  * units = mks (m^3/kmol) -determined by m_VCS_UnitsFormat
1971  * Length = number of species
1972  */
1973  std::vector<double> m_PMVolumeSpecies;
1974 
1975  //! dimensionless value of Faraday's constant, F / RT (1/volt)
1977 
1978  //! Timing and iteration counters for the vcs object
1980 
1981  //! Debug printing lvl
1982  /*!
1983  * Levels correspond to the following guidlines
1984  * * 0 No printing at all
1985  * * 1 Serious warnings or fatal errors get one line
1986  * * 2 one line per eacdh successful vcs package call
1987  * * 3 one line per every successful solve_TP calculation
1988  * * 4 one line for every successful operation -> solve_TP gets a summary report
1989  * * 5 each iteration in solve_TP gets a report with one line per species
1990  * * 6 Each decision in solve_TP gets a line per species in addition to 4
1991  * * 10 Additionally Hessian matrix is printed out
1992  *
1993  * Levels of printing above 4 are only accessible when DEBUG_MODE is turned on
1994  */
1996 
1997  //! printing level of timing information
1998  /*!
1999  * * 1 allowing printing of timing
2000  * * 0 do not allow printing of timing -> everything is printed as a NA.
2001  */
2003 
2004  //! Units for the chemical potential data
2005  /*!
2006  * | Value | chemical potential units | pressure units |
2007  * | ----- | ------------------------ | -------------- |
2008  * | -1 | kcal/mol | Pa |
2009  * | 0 | MU/RT | Pa |
2010  * | 1 | kJ/mol | Pa |
2011  * | 2 | Kelvin | Pa |
2012  * | 3 | J / kmol | Pa |
2013  */
2015 
2016  friend class vcs_phaseStabilitySolve;
2017 };
2018 
2019 }
2020 #endif
size_t vcs_RxnStepSizes(int &forceComponentCalc, size_t &kSpecial)
Calculates formation reaction step sizes.
Definition: vcs_rxnadj.cpp:20
int vcs_zero_species(const size_t kspec)
Zero out the concentration of a species.
void vcs_inest(double *const aw, double *const sa, double *const sm, double *const ss, double test)
Estimate equilibrium compositions.
Definition: vcs_inest.cpp:21
int vcs_elem_rearrange(double *const aw, double *const sa, double *const sm, double *const ss)
Rearrange the constraint equations represented by the Formula Matrix so that the operational ones are...
std::vector< size_t > m_phaseID
Mapping from the species number to the phase number.
Definition: vcs_solve.h:1832
std::vector< double > m_TmpPhase2
Temporary vector of length NPhase.
Definition: vcs_solve.h:1725
double vcs_Hessian_actCoeff_diag(size_t irxn)
Calculates the diagonal contribution to the Hessian due to the dependence of the activity coefficient...
Definition: vcs_rxnadj.cpp:594
double vcs_nondimMult_TP(int mu_units, double TKelvin) const
Returns the multiplier for the nondimensionalization of the equations.
Definition: vcs_nondim.cpp:38
double vcs_line_search(const size_t irxn, const double dx_orig, char *const ANOTE=0)
A line search algorithm is carried out on one reaction.
Definition: vcs_rxnadj.cpp:668
double vcs_minor_alt_calc(size_t kspec, size_t irxn, bool *do_delete, char *ANOTE=0) const
Minor species alternative calculation.
std::vector< int > m_elType
Type of the element constraint.
Definition: vcs_solve.h:1860
double vcs_Total_Gibbs(double *w, double *fe, double *tPhMoles)
Calculate the total dimensionless Gibbs free energy.
Definition: vcs_Gibbs.cpp:16
double vcs_nondim_Farad(int mu_units, double TKelvin) const
Returns the multiplier for electric charge terms.
Definition: vcs_nondim.cpp:18
double m_totalMolNum
Total number of kmoles in all phases.
Definition: vcs_solve.h:1702
double m_Faraday_dim
dimensionless value of Faraday's constant, F / RT (1/volt)
Definition: vcs_solve.h:1976
VCS_COUNTERS * m_VCount
Timing and iteration counters for the vcs object.
Definition: vcs_solve.h:1979
void vcs_chemPotPhase(const int stateCalc, const size_t iph, const double *const molNum, double *const ac, double *const mu_i, const bool do_deleted=false)
We calculate the dimensionless chemical potentials of all species in a single phase.
std::vector< VCS_SPECIES_THERMO * > m_speciesThermoList
Vector of pointers to thermostructures which identify the model and parameters for evaluating the the...
Definition: vcs_solve.h:1955
double vcs_Hessian_diag_adj(size_t irxn, double hessianDiag_Ideal)
Calculates the diagonal contribution to the Hessian due to the dependence of the activity coefficient...
Definition: vcs_rxnadj.cpp:576
std::vector< double > m_wtSpecies
Molecular weight of each species.
Definition: vcs_solve.h:1942
Array2D m_stoichCoeffRxnMatrix
Stoichiometric coefficient matrix for the reaction mechanism expressed in Reduced Canonical Form...
Definition: vcs_solve.h:1552
double vcs_birthGuess(const int kspec)
Birth guess returns the number of moles of a species that is coming back to life. ...
int m_useActCoeffJac
Choice of Hessians.
Definition: vcs_solve.h:1963
std::vector< double > m_feSpecies_new
Dimensionless new free energy for all the species in the mechanism at the new tentatite T...
Definition: vcs_solve.h:1593
int vcs_elcorr(double aa[], double x[])
Definition: vcs_elem.cpp:103
int vcs_TP(int ipr, int ip1, int maxit, double T, double pres)
Solve an equilibrium problem at a particular fixed temperature and pressure.
Definition: vcs_TP.cpp:7
size_t vcs_basisOptMax(const double *const molNum, const size_t j, const size_t n)
Choose a species to test for the next component.
void vcs_initSizes(const size_t nspecies0, const size_t nelements, const size_t nphase0)
Initialize the sizes within the VCS_SOLVE object.
Definition: vcs_solve.cpp:59
void vcs_dfe(const int stateCalc, const int ll, const size_t lbot, const size_t ltop)
Calculate the dimensionless chemical potentials of all species or of certain groups of species...
void vcs_elabPhase(size_t iphase, double *const elemAbundPhase)
Definition: vcs_elem.cpp:89
std::vector< vcs_VolPhase * > m_VolPhaseList
Array of Phase Structures. Length = number of phases.
Definition: vcs_solve.h:1869
std::vector< double > m_phasePhi
electric potential of the iph phase
Definition: vcs_solve.h:1639
void vcs_nondim_TP()
Nondimensionalize the problem data.
Definition: vcs_nondim.cpp:60
double vcs_tmoles()
Calculates the total number of moles of species in all phases.
std::vector< double > m_lnMnaughtSpecies
specifies the ln(Mnaught) used to calculate the chemical potentials
Definition: vcs_solve.h:1913
This file contains definitions of terms that are used in internal routines and are unlikely to need m...
std::vector< double > m_spSize
total size of the species
Definition: vcs_solve.h:1569
std::vector< size_t > m_elementMapIndex
Index vector that keeps track of the rearrangement of the elements.
Definition: vcs_solve.h:1808
bool vcs_wellPosed(VCS_PROB *vprob)
In this routine, we check for things that will cause the algorithm to fail.
Definition: vcs_prep.cpp:225
int vcs_prob_update(VCS_PROB *pub)
Transfer the results of the equilibrium calculation back to VCS_PROB.
Definition: vcs_solve.cpp:880
void vcs_deltag_Phase(const size_t iphase, const bool doDeleted, const int stateCalc, const bool alterZeroedPhases=true)
Calculate deltag of formation for all species in a single phase.
double m_tolmin
Tolerance requirements for minor species.
Definition: vcs_solve.h:1761
void vcs_redim_TP()
Redimensionalize the problem data.
Definition: vcs_nondim.cpp:145
void vcs_printChemPotUnits(int unitsFormat) const
Print the string representing the Chemical potential units.
Definition: vcs_nondim.cpp:189
int delta_species(const size_t kspec, double *const delta_ptr)
Change the concentration of a species by delta moles.
int m_VCS_UnitsFormat
Units for the chemical potential data.
Definition: vcs_solve.h:2014
size_t m_numRxnRdc
Current number of non-component species in the problem.
Definition: vcs_solve.h:1519
bool vcs_globStepDamp()
This routine optimizes the minimization of the total Gibbs free energy by making sure the slope of th...
std::vector< double > m_chargeSpecies
Charge of each species. Length = number of species.
Definition: vcs_solve.h:1945
size_t m_numSpeciesRdc
Current number of species in the problems.
Definition: vcs_solve.h:1512
std::vector< int > m_phaseActConvention
specifies the activity convention of the phase.
Definition: vcs_solve.h:1906
int vcs_prep()
Prepare the object for solution.
Definition: vcs_prep.cpp:206
void vcs_deltag(const int l, const bool doDeleted, const int vcsState, const bool alterZeroedPhases=true)
This subroutine calculates reaction free energy changes for all noncomponent formation reactions...
std::vector< double > m_molNumSpecies_old
Total moles of the species.
Definition: vcs_solve.h:1612
A class for 2D arrays stored in column-major (Fortran-compatible) form.
Definition: Array.h:29
Header file for class Cantera::Array2D.
double vcs_GibbsPhase(size_t iphase, const double *const w, const double *const fe)
Calculate the total dimensionless Gibbs free energy of a single phase.
Definition: vcs_Gibbs.cpp:41
std::vector< std::string > m_elementName
Vector of strings containing the element names.
Definition: vcs_solve.h:1845
double vcs_phaseStabilityTest(const size_t iph)
Main program to test whether a deleted phase should be brought back into existence.
std::string m_title
String containing the title of the run.
Definition: vcs_solve.h:1872
int vcs_report(int iconv)
Print out a report on the state of the equilibrium problem to standard output.
Definition: vcs_report.cpp:14
Defines and definitions within the vcs package.
int m_doEstimateEquil
Setting for whether to do an initial estimate.
Definition: vcs_solve.h:1605
int vcs_inest_TP()
Create an initial estimate of the solution to the thermodynamic equilibrium problem.
Definition: vcs_inest.cpp:320
std::vector< size_t > m_speciesMapIndex
Index vector that keeps track of the species vector rearrangement.
Definition: vcs_solve.h:1781
Internal declarations for the VCSnonideal package.
double vcs_VolTotal(const double tkelvin, const double pres, const double w[], double volPM[])
Calculation of the total volume and the partial molar volumes.
Definition: vcs_solve.cpp:983
size_t vcs_add_all_deleted()
Provide an estimate for the deleted species in phases that are not zeroed out.
std::vector< double > m_tPhaseMoles_new
total kmols of species in each phase in the tentative soln vector
Definition: vcs_solve.h:1719
void vcs_printSpeciesChemPot(const int stateCalc) const
Print out a table of chemical potentials.
void vcs_updateVP(const int stateCalc)
This routine uploads the state of the system into all of the vcs_VolumePhase objects in the current p...
std::vector< double > m_elemAbundances
Element abundances vector.
Definition: vcs_solve.h:1686
std::vector< double > m_actCoeffSpecies_new
Molar-based Activity Coefficients for Species.
Definition: vcs_solve.h:1917
size_t m_numElemConstraints
Number of element constraints in the problem.
Definition: vcs_solve.h:1499
void vcs_TCounters_report(int timing_print_lvl=1)
Create a report on the plog file containing timing and its information.
Definition: vcs_report.cpp:373
std::vector< int > m_actConventionSpecies
specifies the activity convention of the phase containing the species
Definition: vcs_solve.h:1897
double m_tolmin2
Below this, minor species aren't refined any more.
Definition: vcs_solve.h:1767
size_t m_numSpeciesTot
Total number of species in the problems.
Definition: vcs_solve.h:1493
std::vector< int > m_speciesUnknownType
Specifies the species unknown type.
Definition: vcs_solve.h:1625
int vcs_recheck_deleted()
Recheck deleted species in multispecies phases.
double l2normdg(double dg[]) const
Calculate the norm of a deltaGibbs free energy vector.
std::vector< double > m_scSize
Absolute size of the stoichiometric coefficients.
Definition: vcs_solve.h:1562
bool recheck_deleted_phase(const int iphase)
Recheck deletion condition for multispecies phases.
int vcs_solve_phaseStability(const int iphase, int ifunc, double &funcval, int print_lvl)
Routine that independently determines whether a phase should be popped under the current conditions...
int vcs(VCS_PROB *vprob, int ifunc, int ipr, int ip1, int maxit)
Solve an equilibrium problem.
Definition: vcs_solve.cpp:244
std::vector< double > m_tPhaseMoles_old
Total kmols of species in each phase.
Definition: vcs_solve.h:1710
std::vector< double > m_deltaPhaseMoles
Change in the total moles in each phase.
Definition: vcs_solve.h:1731
int vcs_species_type(const size_t kspec) const
Evaluate the species category for the indicated species.
int m_debug_print_lvl
Debug printing lvl.
Definition: vcs_solve.h:1995
int m_timing_print_lvl
printing level of timing information
Definition: vcs_solve.h:2002
int vcs_evalSS_TP(int ipr, int ip1, double Temp, double pres)
Definition: vcs_TP.cpp:58
size_t m_numComponents
Number of components calculated for the problem.
Definition: vcs_solve.h:1502
std::vector< double > m_molNumSpecies_new
Tentative value of the mole number vector.
Definition: vcs_solve.h:1644
std::vector< double > m_SSfeSpecies
Standard state chemical potentials for species K at the current temperature and pressure.
Definition: vcs_solve.h:1577
size_t vcs_popPhaseID(std::vector< size_t > &phasePopPhaseIDs)
Decision as to whether a phase pops back into existence.
std::vector< double > m_actCoeffSpecies_old
Molar-based Activity Coefficients for Species based on old mole numbers.
Definition: vcs_solve.h:1924
int vcs_rxn_adj_cg(void)
Calculates reaction adjustments using a full Hessian approximation.
Definition: vcs_rxnadj.cpp:368
size_t NPHASE0
value of the number of phases used to malloc data structures
Definition: vcs_solve.h:1490
int vcs_prep_oneTime(int printLvl)
This routine is mostly concerned with changing the private data to be consistent with what's needed f...
Definition: vcs_prep.cpp:60
size_t m_numRxnMinorZeroed
Number of active species which are currently either treated as minor species.
Definition: vcs_solve.h:1523
int vcs_delete_species(const size_t kspec)
Change a single species from active to inactive status.
double m_totalVol
Total volume of all phases. Units are m^3.
Definition: vcs_solve.h:1966
double m_temperature
Temperature (Kelvin)
Definition: vcs_solve.h:1734
void vcs_reinsert_deleted(size_t kspec)
std::vector< char > m_SSPhase
Boolean indicating whether a species belongs to a single-species phase.
Definition: vcs_solve.h:1836
int vcs_rearrange()
Switch all species data back to the original order.
int vcs_phasePopDeterminePossibleList()
Determine the list of problems that need to be checked to see if there are any phases pops...
std::vector< std::string > m_speciesName
Species string name for the kth species.
Definition: vcs_solve.h:1839
std::vector< double > m_TmpPhase
Temporary vector of length NPhase.
Definition: vcs_solve.h:1722
Array2D m_deltaMolNumPhase
Change in the number of moles of phase, iphase, due to the noncomponent formation reaction...
Definition: vcs_solve.h:1632
Array2D m_formulaMatrix
Formula matrix for the problem.
Definition: vcs_solve.h:1534
std::vector< int > m_elementActive
Specifies whether an element constraint is active.
Definition: vcs_solve.h:1866
double m_pressurePA
Pressure (units are determined by m_VCS_UnitsFormat.
Definition: vcs_solve.h:1748
double m_tolmaj
Tolerance requirement for major species.
Definition: vcs_solve.h:1758
void vcs_elab()
Computes the current elemental abundances vector.
Definition: vcs_elem.cpp:13
int vcs_setMolesLinProg()
Estimate the initial mole numbers by constrained linear programming.
bool vcs_popPhasePossible(const size_t iphasePop) const
Utility function that evaluates whether a phase can be popped into existence.
std::vector< int > m_speciesStatus
Major -Minor status vector for the species in the problem.
Definition: vcs_solve.h:1829
double deltaG_Recalc_Rxn(const int stateCalc, const size_t irxn, const double *const molNum, double *const ac, double *const mu_i)
This function recalculates the deltaG for reaction, irxn.
Definition: vcs_rxnadj.cpp:652
std::vector< double > m_PMVolumeSpecies
Partial molar volumes of the species.
Definition: vcs_solve.h:1973
void vcs_SSPhase(void)
Calculate the status of single species phases.
Definition: vcs_prep.cpp:18
int vcs_basopt(const bool doJustComponents, double aw[], double sa[], double sm[], double ss[], double test, bool *const usedZeroedSpecies)
Choose the optimum species basis for the calculations.
int vcs_prob_specify(const VCS_PROB *pub)
Specify the problem to be solved using VCS_PROB, incrementally.
Definition: vcs_solve.cpp:756
void vcs_switch_pos(const bool ifunc, const size_t k1, const size_t k2)
Swaps the indices for all of the global data for two species, k1 and k2.
void prneav() const
Print out and check the elemental abundance vector.
void vcs_updateMolNumVolPhases(const int stateCalc)
Update all underlying vcs_VolPhase objects.
int vcs_popPhaseRxnStepSizes(const size_t iphasePop)
Calculates the deltas of the reactions due to phases popping into existence.
std::vector< double > m_deltaMolNumSpecies
Reaction Adjustments for each species during the current step.
Definition: vcs_solve.h:1674
std::vector< double > m_deltaGRxn_tmp
Temporary vector of Rxn DeltaG's.
Definition: vcs_solve.h:1667
Class to keep track of time and iterations.
Definition: vcs_internal.h:49
int vcs_rank(const double *awtmp, size_t numSpecies, const double *matrix, size_t numElemConstraints, std::vector< size_t > &compRes, std::vector< size_t > &elemComp, int *const usedZeroedSpecies) const
Calculate the rank of a matrix and return the rows and columns that will generate an independent basi...
Definition: vcs_rank.cpp:38
void vcs_delete_memory()
Delete memory that isn't just resizable STL containers.
Definition: vcs_solve.cpp:221
Array2D m_phaseParticipation
This is 1 if the phase, iphase, participates in the formation reaction irxn, and zero otherwise...
Definition: vcs_solve.h:1636
size_t NSPECIES0
value of the number of species used to malloc data structures
Definition: vcs_solve.h:1487
char m_unitsState
This specifies the current state of units for the Gibbs free energy properties in the program...
Definition: vcs_solve.h:1879
std::vector< size_t > m_indexRxnToSpecies
Mapping between the species index for noncomponent species and the full species index.
Definition: vcs_solve.h:1821
int vcs_solve_TP(int print_lvl, int printDetails, int maxit)
Main routine that solves for equilibrium at constant T and P using a variant of the VCS method...
size_t m_numPhases
Number of Phases in the problem.
Definition: vcs_solve.h:1526
Interface class for the vcs thermo equilibrium solver package, which generally describes the problem ...
Definition: vcs_prob.h:24
double m_totalMoleScale
Multiplier for the mole numbers within the nondimensionless formulation.
Definition: vcs_solve.h:1888
size_t m_numRxnTot
Total number of non-component species in the problem.
Definition: vcs_solve.h:1505
void vcs_fePrep_TP()
Initialize the chemical potential of single species phases.
Definition: vcs_TP.cpp:85
double m_tolmaj2
Below this, major species aren't refined any more.
Definition: vcs_solve.h:1764
bool vcs_evaluate_speciesType()
This routine evaluates the species type for all species.
void vcs_switch_elem_pos(size_t ipos, size_t jpos)
Swaps the indices for all of the global data for two elements, ipos and jpos.
std::vector< double > TPhInertMoles
Total kmoles of inert to add to each phase.
Definition: vcs_solve.h:1755
Array2D m_np_dLnActCoeffdMolNum
Change in the log of the activity coefficient with respect to the mole number multiplied by the phase...
Definition: vcs_solve.h:1934
bool vcs_delete_multiphase(const size_t iph)
This routine handles the bookkeeping involved with the deletion of multiphase phases from the problem...
std::vector< double > m_elemAbundancesGoal
Element abundances vector Goals.
Definition: vcs_solve.h:1695
std::vector< double > m_feSpecies_old
Free energy vector from the start of the current iteration.
Definition: vcs_solve.h:1584
int vcs_prob_specifyFully(const VCS_PROB *pub)
Fully specify the problem to be solved using VCS_PROB.
Definition: vcs_solve.cpp:375
void vcs_counters_init(int ifunc)
Initialize the internal counters.
Definition: vcs_solve.cpp:965
std::vector< double > m_deltaGRxn_Deficient
Last deltag[irxn] from the previous step with additions for possible births of zeroed phases...
Definition: vcs_solve.h:1661
std::vector< double > m_deltaGRxn_new
Delta G(irxn) for the noncomponent species in the mechanism.
Definition: vcs_solve.h:1654
void vcs_CalcLnActCoeffJac(const double *const moleSpeciesVCS)
Recalculate all of the activity coefficients in all of the phases based on input mole numbers...
Definition: vcs_rxnadj.cpp:627
This is the main structure used to hold the internal data used in vcs_solve_TP(), and to solve TP sys...
Definition: vcs_solve.h:49
bool vcs_elabcheck(int ibound)
Definition: vcs_elem.cpp:25
std::vector< size_t > m_speciesLocalPhaseIndex
Index that keeps track of the index of the species within the local phase.
Definition: vcs_solve.h:1795
std::vector< double > m_deltaGRxn_old
Last deltag[irxn] from the previous step.
Definition: vcs_solve.h:1657