Cantera  4.0.0a1
Loading...
Searching...
No Matches
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 equilGroup ).
5 */
6
7// This file is part of Cantera. See License.txt in the top-level directory or
8// at https://cantera.org/license.txt for license and copyright information.
9
10#ifndef _VCS_SOLVE_H
11#define _VCS_SOLVE_H
12
13/*
14* Index of Symbols
15* -------------------
16* irxn -> refers to the species or rxn between the species and
17* the components in the problem
18* k -> refers to the species
19* j -> refers to the element or component
20*
21* ### -> to be eliminated
22*/
23
27#include "cantera/base/Array.h"
28
29namespace Cantera
30{
31
32class vcs_VolPhase;
33class VCS_COUNTERS;
34class MultiPhase;
35
36//! This is the main structure used to hold the internal data
37//! used in vcs_solve_TP(), and to solve TP systems.
38/*!
39 * The indices of information in this structure may change when the species
40 * basis changes or when phases pop in and out of existence. Both of these
41 * operations change the species ordering.
42 */
44{
45public:
46 //! Initialize the sizes within the VCS_SOLVE object
47 /*!
48 * This resizes all of the internal arrays within the object.
49 */
50 VCS_SOLVE(MultiPhase* mphase, int printLvl=0);
51
52 ~VCS_SOLVE();
53
54 //! Main routine that solves for equilibrium at constant T and P using a
55 //! variant of the VCS method
56 /*!
57 * Any number of single-species phases and multi-species phases, including non-ideal
58 * phases, can be handled by the present version.
59 *
60 * @param print_lvl 1 -> Print results to standard output;
61 * 0 -> don't report on anything
62 * @param printDetails 1 -> Print intermediate results.
63 * @param maxit Maximum number of iterations for the algorithm
64 * @return
65 * * 0 = Equilibrium Achieved
66 * * 1 = Range space error encountered. The element abundance criteria
67 * are only partially satisfied. Specifically, the first NC= (number
68 * of components) conditions are satisfied. However, the full NE
69 * (number of elements) conditions are not satisfied. The
70 * equilibrium condition is returned.
71 * * -1 = Maximum number of iterations is exceeded. Convergence was
72 * not found.
73 */
74 int solve_TP(int print_lvl, int printDetails, int maxit);
75
76 /**
77 * We make decisions on the initial mole number, and major-minor status
78 * here. We also fix up the total moles in a phase.
79 *
80 * irxn = id of the noncomponent species formation reaction for the
81 * species to be added in.
82 *
83 * The algorithm proceeds to implement these decisions in the previous
84 * position of the species. Then, vcs_switch_pos is called to move the
85 * species into the last active species slot, incrementing the number of
86 * active species at the same time.
87 *
88 * This routine is responsible for the global data manipulation only.
89 */
90 void vcs_reinsert_deleted(size_t kspec);
91
92 //! Choose the optimum species basis for the calculations
93 /*!
94 * This is done by choosing the species with the largest mole fraction not
95 * currently a linear combination of the previous components. Then,
96 * calculate the stoichiometric coefficient matrix for that basis.
97 *
98 * Rearranges the solution data to put the component data at the
99 * front of the species list.
100 *
101 * Then, calculates `m_stoichCoeffRxnMatrix(jcomp,irxn)` the formation
102 * reactions for all noncomponent species in the mechanism. Also
103 * calculates DNG(I) and DNL(I), the net mole change for each formation
104 * reaction. Also, initializes IR(I) to the default state.
105 *
106 * @param[in] doJustComponents If true, the #m_stoichCoeffRxnMatrix and
107 * #m_deltaMolNumPhase are not calculated.
108 * @param[in] test This is a small negative number dependent upon whether
109 * an estimate is supplied or not.
110 *
111 * ### Internal Variables calculated by this routine:
112 *
113 * - #m_numComponents: Number of component species. This routine
114 * calculates the #m_numComponents species. It switches their positions
115 * in the species vector so that they occupy the first #m_numComponents
116 * spots in the species vector.
117 * - `m_stoichCoeffRxnMatrix(jcomp,irxn)` Stoichiometric coefficient
118 * matrix for the reaction mechanism expressed in Reduced Canonical
119 * Form. jcomp refers to the component number, and irxn refers to the
120 * irxn_th non-component species.
121 * - `m_deltaMolNumPhase(iphase,irxn)`: Change in the number of moles in
122 * phase, iphase, due to the noncomponent formation reaction, irxn.
123 * - `m_phaseParticipation(iphase,irxn)`: This is 1 if the phase, iphase,
124 * participates in the formation reaction, irxn, and zero otherwise.
125 */
126 void vcs_basopt(const bool doJustComponents, double test);
127
128 //! Evaluate the species category for the indicated species
129 /*!
130 * All evaluations are done using the "old" version of the solution.
131 *
132 * @param kspec Species to be evaluated
133 * @returns the calculated species type
134 */
135 int vcs_species_type(const size_t kspec) const;
136
137 //! This routine evaluates the species type for all species
138 /*!
139 * - #VCS_SPECIES_MAJOR: Major species
140 * - #VCS_SPECIES_MINOR: Minor species
141 * - #VCS_SPECIES_ZEROEDMS: The species lies in a multicomponent phase
142 * which currently doesn't exist. Its concentration is currently zero.
143 * - #VCS_SPECIES_ZEROEDSS: Species lies in a single-species phase which
144 * is currently zeroed out.
145 * - #VCS_SPECIES_DELETED: Species has such a small mole fraction it is
146 * deleted even though its phase may possibly exist. The species is
147 * believed to have such a small mole fraction that it best to throw
148 * the calculation of it out. It will be added back in at the end of
149 * the calculation.
150 * - #VCS_SPECIES_INTERFACIALVOLTAGE: Species refers to an electron in
151 * the metal The unknown is equal to the interfacial voltage drop
152 * across the interface on the SHE (standard hydrogen electrode) scale
153 * (volts).
154 * - #VCS_SPECIES_ZEROEDPHASE: Species lies in a multicomponent phase
155 * that is zeroed atm and will stay deleted due to a choice from a
156 * higher level. These species will formally always have zero mole
157 * numbers in the solution vector.
158 * - #VCS_SPECIES_ACTIVEBUTZERO: The species lies in a multicomponent
159 * phase which currently does exist. Its concentration is currently
160 * identically zero, though the phase exists. Note, this is a temporary
161 * condition that exists at the start of an equilibrium problem. The
162 * species is soon "birthed" or "deleted".
163 * - #VCS_SPECIES_STOICHZERO: The species lies in a multicomponent phase
164 * which currently does exist. Its concentration is currently
165 * identically zero, though the phase exists. This is a permanent
166 * condition due to stoich constraints
167 */
169
170 //! Calculate the dimensionless chemical potentials of all species or
171 //! of certain groups of species, at a fixed temperature and pressure.
172 /*!
173 * We calculate the dimensionless chemical potentials of all species or
174 * certain groups of species here, at a fixed temperature and pressure, for
175 * the input mole vector z[] in the parameter list. Nondimensionalization is
176 * achieved by division by RT.
177 *
178 * Note, for multispecies phases which are currently zeroed out, the
179 * chemical potential is filled out with the standard chemical potential.
180 *
181 * For species in multispecies phases whose concentration is zero, we need
182 * to set the mole fraction to a very low value. Its chemical potential is
183 * then calculated using the VCS_DELETE_MINORSPECIES_CUTOFF concentration
184 * to keep numbers positive.
185 *
186 * For formulas, see vcs_chemPotPhase().
187 *
188 * Handling of Small Species:
189 * ------------------------------
190 * As per the discussion above, for small species where the mole fraction
191 *
192 * z(i) < VCS_DELETE_MINORSPECIES_CUTOFF
193 *
194 * The chemical potential is calculated as:
195 *
196 * m_feSpecies(I)(I) = m_SSfeSpecies(I) + ln(ActCoeff[i](VCS_DELETE_MINORSPECIES_CUTOFF))
197 *
198 * Species in the following categories are treated as "small species"
199 * - #VCS_SPECIES_DELETED
200 * - #VCS_SPECIES_ACTIVEBUTZERO
201 *
202 * For species in multispecies phases which are currently not active, the
203 * treatment is different. These species are in the following species
204 * categories:
205 * - #VCS_SPECIES_ZEROEDMS
206 * - #VCS_SPECIES_ZEROEDPHASE
207 *
208 * For these species, the `ln( ActCoeff[I] X[I])` term is dropped
209 * altogether. The following equation is used:
210 *
211 * m_feSpecies(I) = m_SSfeSpecies(I)
212 * + Charge[I] * Faraday_dim * phasePhi[iphase];
213 *
214 * Handling of "Species" Representing Interfacial Voltages
215 * ---------------------------------------------------------
216 *
217 * These species have species types of
218 * #VCS_SPECIES_TYPE_INTERFACIALVOLTAGE The chemical potentials for these
219 * "species" refer to electrons in metal electrodes. They have the
220 * following formula
221 *
222 * m_feSpecies(I) = m_SSfeSpecies(I) - F z[I] / RT
223 *
224 * - `F` is Faraday's constant.
225 * - `R` = gas constant
226 * - `T` = temperature
227 * - `V` = potential of the interface = phi_electrode - phi_solution
228 *
229 * For these species, the solution vector unknown, z[I], is V, the phase
230 * voltage, in volts.
231 *
232 * @param ll Determine which group of species gets updated
233 * - `ll = 0`: Calculate for all species
234 * - `ll < 0`: calculate for components and for major non-components
235 * - `ll = 1`: calculate for components and for minor non-components
236 *
237 * @param lbot Restricts the calculation of the chemical potential
238 * to the species between LBOT <= i < LTOP. Usually
239 * LBOT and LTOP will be equal to 0 and MR, respectively.
240 * @param ltop Top value of the loops
241 *
242 * @param stateCalc Determines whether z is old or new or tentative:
243 * - 1: Use the tentative values for the total number of
244 * moles in the phases, that is, use TG1 instead of TG etc.
245 * - 0: Use the base values of the total number of
246 * moles in each system.
247 *
248 * Also needed:
249 * ff : standard state chemical potentials. These are the
250 * chemical potentials of the standard states at
251 * the same T and P as the solution.
252 * tg : Total Number of moles in the phase.
253 */
254 void vcs_dfe(const int stateCalc, const int ll, const size_t lbot, const size_t ltop);
255
256 //! This routine uploads the state of the system into all of the
257 //! vcs_VolumePhase objects in the current problem.
258 /*!
259 * @param stateCalc Determines where to get the mole numbers from.
260 * - VCS_STATECALC_OLD -> from m_molNumSpecies_old
261 * - VCS_STATECALC_NEW -> from m_molNumSpecies_new
262 */
263 void vcs_updateVP(const int stateCalc);
264
265 //! Utility function that evaluates whether a phase can be popped
266 //! into existence
267 /*!
268 * A phase can be popped iff the stoichiometric coefficients for the
269 * component species, whose concentrations will be lowered during the
270 * process, are positive by at least a small degree.
271 *
272 * If one of the phase species is a zeroed component, then the phase can
273 * be popped if the component increases in mole number as the phase moles
274 * are increased.
275 *
276 * @param iphasePop id of the phase, which is currently zeroed,
277 * @returns true if the phase can come into existence and false otherwise.
278 */
279 bool vcs_popPhasePossible(const size_t iphasePop) const;
280
281 //! Decision as to whether a phase pops back into existence
282 /*!
283 * @returns the phase id of the phase that pops back into existence. Returns
284 * -1 if there are no phases
285 */
286 size_t vcs_popPhaseID();
287
288 //! Calculates the deltas of the reactions due to phases popping
289 //! into existence
290 /*!
291 * Updates #m_deltaMolNumSpecies[irxn] : reaction adjustments, where irxn
292 * refers to the irxn'th species formation reaction. This adjustment is
293 * for species irxn + M, where M is the number of components.
294 *
295 * @param iphasePop Phase id of the phase that will come into existence
296 * @returns an int representing the status of the step
297 * - 0 : normal return
298 * - 1 : A single species phase species has been zeroed out
299 * in this routine. The species is a noncomponent
300 * - 2 : Same as one but, the zeroed species is a component.
301 * - 3 : Nothing was done because the phase couldn't be birthed
302 * because a needed component is zero.
303 */
304 int vcs_popPhaseRxnStepSizes(const size_t iphasePop);
305
306 //! Calculates formation reaction step sizes.
307 /*!
308 * This is equation 6.4-16, p. 143 in Smith and Missen @cite smith1982.
309 *
310 * Output
311 * -------
312 * m_deltaMolNumSpecies(irxn) : reaction adjustments, where irxn refers to
313 * the irxn'th species formation reaction.
314 * This adjustment is for species irxn + M,
315 * where M is the number of components.
316 *
317 * Special branching occurs sometimes. This causes the component basis
318 * to be reevaluated
319 *
320 * @param forceComponentCalc integer flagging whether a component
321 * recalculation needs to be carried out.
322 * @param kSpecial species number of phase being zeroed.
323 * @returns an int representing which phase may need to be zeroed
324 */
325 size_t vcs_RxnStepSizes(int& forceComponentCalc, size_t& kSpecial);
326
327 //! Calculates the total number of moles of species in all phases.
328 /*!
329 * Also updates the variable m_totalMolNum and Reconciles Phase existence
330 * flags with total moles in each phase.
331 */
332 double vcs_tmoles();
333
334 void check_tmoles() const;
335
336 //! This subroutine calculates reaction free energy changes for
337 //! all noncomponent formation reactions.
338 /*!
339 * Formation reactions are
340 * reactions which create each noncomponent species from the component
341 * species. `m_stoichCoeffRxnMatrix(jcomp,irxn)` are the stoichiometric
342 * coefficients for these reactions. A stoichiometric coefficient of
343 * one is assumed for species irxn in this reaction.
344 *
345 * @param L
346 * - `L < 0`: Calculate reactions corresponding to major noncomponent
347 * and zeroed species only
348 * - `L = 0`: Do all noncomponent reactions, i, between
349 * 0 <= i < irxnl
350 * - `L > 0`: Calculate reactions corresponding to minor noncomponent
351 * and zeroed species only
352 *
353 * @param doDeleted Do deleted species
354 * @param vcsState Calculate deltaG corresponding to either old or new
355 * free energies
356 * @param alterZeroedPhases boolean indicating whether we should
357 * add in a special section for zeroed phases.
358 *
359 * Note we special case one important issue. If the component has zero
360 * moles, then we do not allow deltaG < 0.0 for formation reactions which
361 * would lead to the loss of more of that same component. This dG < 0.0
362 * condition feeds back into the algorithm in several places, and leads
363 * to a infinite loop in at least one case.
364 */
365 void vcs_deltag(const int L, const bool doDeleted, const int vcsState,
366 const bool alterZeroedPhases = true);
367
368 //! Swaps the indices for all of the global data for two species, k1
369 //! and k2.
370 /*!
371 * @param ifunc If true, switch the species data and the noncomponent
372 * reaction data. This must be called for a non-component species only.
373 * If false, switch the species data only. Typically, we use this option
374 * when determining the component species and at the end of the
375 * calculation, when we want to return unscrambled results. All rxn data
376 * will be out-of-date.
377 * @param k1 First species index
378 * @param k2 Second species index
379 */
380 void vcs_switch_pos(const bool ifunc, const size_t k1, const size_t k2);
381
382 //! Main program to test whether a deleted phase should be brought
383 //! back into existence
384 /*!
385 * @param iph Phase id of the deleted phase
386 */
387 double vcs_phaseStabilityTest(const size_t iph);
388
389 /**
390 * Evaluate the standard state free energies at the current temperature
391 * and pressure. Ideal gas pressure contribution is added in here.
392 *
393 * @param ipr 1 -> Print results to standard output; 0 -> don't report
394 * on anything
395 * @param ip1 1 -> Print intermediate results; 0 -> don't.
396 * @param Temp Temperature (Kelvin)
397 * @param pres Pressure (Pascal)
398 */
399 int vcs_evalSS_TP(int ipr, int ip1, double Temp, double pres);
400
401 //! This routine is mostly concerned with changing the private data to be
402 //! consistent with what's needed for solution. It is called one time for
403 //! each new problem structure definition.
404 /*!
405 * The problem structure refers to:
406 *
407 * - the number and identity of the species.
408 * - the formula matrix and thus the number of components.
409 * - the number and identity of the phases.
410 * - the equation of state
411 * - the method and parameters for determining the standard state
412 * - The method and parameters for determining the activity coefficients.
413 *
414 * Tasks:
415 * 1. Fill in the SSPhase[] array.
416 * 2. Check to see if any multispecies phases actually have only one
417 * species in that phase. If true, reassign that phase and species
418 * to be a single-species phase.
419 * 3. Determine the number of components in the problem if not already
420 * done so. During this process the order of the species is changed
421 * in the private data structure. All references to the species
422 * properties must employ the ind[] index vector.
423 * 4. Initialization of arrays to zero.
424 * 5. Check to see if the problem is well posed (If all the element
425 * abundances are zero, the algorithm will fail)
426 *
427 * @param printLvl Print level of the routine
428 */
429 void vcs_prep(int printLvl);
430
431 //! Rearrange the constraint equations represented by the Formula
432 //! Matrix so that the operational ones are in the front
433 /*!
434 * This subroutine handles the rearrangement of the constraint equations
435 * represented by the Formula Matrix. Rearrangement is only necessary when
436 * the number of components is less than the number of elements. For this
437 * case, some constraints can never be satisfied exactly, because the
438 * range space represented by the Formula Matrix of the components can't
439 * span the extra space. These constraints, which are out of the range
440 * space of the component Formula matrix entries, are migrated to the back
441 * of the Formula matrix.
442 *
443 * A prototypical example is an extra element column in FormulaMatrix[],
444 * which is identically zero. For example, let's say that argon is has an
445 * element column in FormulaMatrix[], but no species in the mechanism
446 * actually contains argon. Then, nc < ne. Also, without perturbation of
447 * FormulaMatrix[] vcs_basopt[] would produce a zero pivot because the
448 * matrix would be singular (unless the argon element column was already
449 * the last column of FormulaMatrix[].
450 *
451 * This routine borrows heavily from vcs_basopt's algorithm. It finds nc
452 * constraints which span the range space of the Component Formula matrix,
453 * and assigns them as the first nc components in the formula matrix. This
454 * guarantees that vcs_basopt[] has a nonsingular matrix to invert.
455 */
456 void vcs_elem_rearrange();
457
458 //! Swaps the indices for all of the global data for two elements, ipos
459 //! and jpos.
460 /*!
461 * This function knows all of the element information with VCS_SOLVE, and
462 * can therefore switch element positions
463 *
464 * @param ipos first global element index
465 * @param jpos second global element index
466 */
467 void vcs_switch_elem_pos(size_t ipos, size_t jpos);
468
469 //! Calculates the diagonal contribution to the Hessian due to
470 //! the dependence of the activity coefficients on the mole numbers.
471 /*!
472 * (See framemaker notes, Eqn. 20 - VCS Equations document)
473 *
474 * We allow the diagonal to be increased positively to any degree.
475 * We allow the diagonal to be decreased to 1/3 of the ideal solution
476 * value, but no more -> it must remain positive.
477 */
478 double vcs_Hessian_diag_adj(size_t irxn, double hessianDiag_Ideal);
479
480 //! Calculates the diagonal contribution to the Hessian due to
481 //! the dependence of the activity coefficients on the mole numbers.
482 /*!
483 * (See framemaker notes, Eqn. 20 - VCS Equations document)
484 */
485 double vcs_Hessian_actCoeff_diag(size_t irxn);
486
487 //! Recalculate all of the activity coefficients in all of the phases
488 //! based on input mole numbers
489 /**
490 * @param moleSpeciesVCS kmol of species to be used in the update.
491 *
492 * NOTE: This routine needs to be regulated.
493 */
494 void vcs_CalcLnActCoeffJac(span<const double> moleSpeciesVCS);
495
496 //! Print out a report on the state of the equilibrium problem to
497 //! standard output.
498 /*!
499 * @param iconv Indicator of convergence, to be printed out in the report:
500 * - 0 converged
501 * - 1 range space error
502 * - -1 not converged
503 */
504 void vcs_report(int iconv);
505
506 //! Computes the current elemental abundances vector
507 /*!
508 * Computes the elemental abundances vector, m_elemAbundances[], and stores
509 * it back into the global structure
510 */
511 void vcs_elab();
512
513 /**
514 * Checks to see if the element abundances are in compliance. If they are,
515 * then TRUE is returned. If not, FALSE is returned. Note the number of
516 * constraints checked is usually equal to the number of components in the
517 * problem. This routine can check satisfaction of all of the constraints
518 * in the problem, which is equal to ne. However, the solver can't fix
519 * breakage of constraints above nc, because that nc is the range space by
520 * definition. Satisfaction of extra constraints would have had to occur
521 * in the problem specification.
522 *
523 * The constraints should be broken up into 2 sections. If a constraint
524 * involves a formula matrix with positive and negative signs, and eaSet =
525 * 0.0, then you can't expect that the sum will be zero. There may be
526 * roundoff that inhibits this. However, if the formula matrix is all of
527 * one sign, then this requires that all species with nonzero entries in
528 * the formula matrix be identically zero. We put this into the logic
529 * below.
530 *
531 * @param ibound 1: Checks constraints up to the number of elements;
532 * 0: Checks constraints up to the number of components.
533 */
534 bool vcs_elabcheck(int ibound);
535
536
537 /**
538 * This subroutine corrects for element abundances. At the end of the
539 * subroutine, the total moles in all phases are recalculated again,
540 * because we have changed the number of moles in this routine.
541 *
542 * @param aa temporary work vector, length ne*ne
543 * @param x temporary work vector, length ne
544 *
545 * @returns
546 * - 0 = Nothing of significance happened, Element abundances were and
547 * still are good.
548 * - 1 = The solution changed significantly; The element abundances are
549 * now good.
550 * - 2 = The solution changed significantly, The element abundances are
551 * still bad.
552 * - 3 = The solution changed significantly, The element abundances are
553 * still bad and a component species got zeroed out.
554 *
555 * Internal data to be worked on::
556 * - ga Current element abundances
557 * - m_elemAbundancesGoal Required elemental abundances
558 * - m_molNumSpecies_old Current mole number of species.
559 * - m_formulaMatrix[][] Formula matrix of the species
560 * - ne Number of elements
561 * - nc Number of components.
562 *
563 * NOTES:
564 * This routine is turning out to be very problematic. There are
565 * lots of special cases and problems with zeroing out species.
566 *
567 * Still need to check out when we do loops over nc vs. ne.
568 */
569 int vcs_elcorr(span<double> aa, span<double> x);
570
571 //! Create an initial estimate of the solution to the equilibrium problem.
572 void vcs_inest();
573
574 //! Estimate the initial mole numbers by constrained linear programming
575 /*!
576 * This is done by running each reaction as far forward or backward as
577 * possible, subject to the constraint that all mole numbers remain non-
578 * negative. Reactions for which @f$ \Delta \mu^0 @f$ are positive are run
579 * in reverse, and ones for which it is negative are run in the forward
580 * direction. The end result is equivalent to solving the linear
581 * programming problem of minimizing the linear Gibbs function subject to
582 * the element and non-negativity constraints.
583 */
585
586 //! Calculate the total dimensionless Gibbs free energy
587 /*!
588 * Note, for this algorithm this function should be monotonically decreasing.
589 */
590 double vcs_Total_Gibbs(span<const double> w, span<const double> fe,
591 span<const double> tPhMoles);
592
593 //! Fully specify the problem to be solved
595
596private:
597 //! Zero out the concentration of a species.
598 /*!
599 * Make sure to conserve elements and keep track of the total moles in all
600 * phases.
601 * - w[]
602 * - m_tPhaseMoles_old[]
603 *
604 * @param kspec Species index
605 * @return:
606 * 1: succeeded
607 * 0: failed.
608 */
609 int vcs_zero_species(const size_t kspec);
610
611 //! Change a single species from active to inactive status
612 /*!
613 * Rearrange data when species is added or removed. The Lth species is
614 * moved to the back of the species vector. The back of the species
615 * vector is indicated by the value of MR, the current number of
616 * active species in the mechanism.
617 *
618 * @param kspec Species Index
619 * @return
620 * Returns 0 unless.
621 * The return is 1 when the current number of
622 * noncomponent species is equal to zero. A recheck of deleted species
623 * is carried out in the main code.
624 */
625 int vcs_delete_species(const size_t kspec);
626
627 //! This routine handles the bookkeeping involved with the deletion of
628 //! multiphase phases from the problem.
629 /*!
630 * When they are deleted, all of their species become active species, even
631 * though their mole numbers are set to zero. The routine does not make the
632 * decision to eliminate multiphases.
633 *
634 * Note, species in phases with zero mole numbers are still considered
635 * active. Whether the phase pops back into existence or not is checked as
636 * part of the main iteration loop.
637 *
638 * @param iph Phase to be deleted
639 * @returns whether the operation was successful or not
640 */
641 bool vcs_delete_multiphase(const size_t iph);
642
643 //! Change the concentration of a species by delta moles.
644 /*!
645 * Make sure to conserve elements and keep track of the total kmoles in
646 * all phases.
647 *
648 * @param kspec The species index
649 * @param delta The delta for the species. This may change during the calculation.
650 * @return
651 * 1: succeeded without change of dx
652 * 0: Had to adjust dx, perhaps to zero, in order to do the delta.
653 */
654 int delta_species(const size_t kspec, double& delta);
655
656 //! Provide an estimate for the deleted species in phases that are not
657 //! zeroed out
658 /*!
659 * Try to add back in all deleted species. An estimate of the kmol numbers
660 * are obtained and the species is added back into the equation system, into
661 * the old state vector.
662 *
663 * This routine is called at the end of the calculation, just before
664 * returning to the user.
665 */
666 size_t vcs_add_all_deleted();
667
668 //! Recheck deleted species in multispecies phases.
669 /*!
670 * We are checking the equation:
671 *
672 * sum_u = sum_j_comp [ sigma_i_j * u_j ]
673 * = u_i_O + \ln((AC_i * W_i)/m_tPhaseMoles_old)
674 *
675 * by first evaluating:
676 *
677 * DG_i_O = u_i_O - sum_u.
678 *
679 * Then, if TL is zero, the phase pops into existence if DG_i_O < 0. Also,
680 * if the phase exists, then we check to see if the species can have a mole
681 * number larger than VCS_DELETE_SPECIES_CUTOFF (default value = 1.0E-32).
682 */
684
685 //! Minor species alternative calculation
686 /*!
687 * This is based upon the following approximation: The mole fraction
688 * changes due to these reactions don't affect the mole numbers of the
689 * component species. Therefore the following approximation is valid for
690 * a small component of an ideal phase:
691 *
692 * 0 = m_deltaGRxn_old(I) + log(molNum_new(I)/molNum_old(I))
693 *
694 * `m_deltaGRxn_old` contains the contribution from
695 *
696 * m_feSpecies_old(I) =
697 * m_SSfeSpecies(I) +
698 * log(ActCoeff[i] * molNum_old(I) / m_tPhaseMoles_old(iph))
699 * Thus,
700 *
701 * molNum_new(I)= molNum_old(I) * EXP(-m_deltaGRxn_old(I))
702 *
703 * Most of this section is mainly restricting the update to reasonable
704 * values. We restrict the update a factor of 1.0E10 up and 1.0E-10 down
705 * because we run into trouble with the addition operator due to roundoff if
706 * we go larger than ~1.0E15. Roundoff will then sometimes produce zero mole
707 * fractions.
708 *
709 * Note: This routine was generalized to incorporate nonideal phases and
710 * phases on the molality basis
711 *
712 * @param[in] kspec The current species and corresponding formation
713 * reaction number.
714 * @param[in] irxn The current species and corresponding formation
715 * reaction number.
716 * @param[out] do_delete: BOOLEAN which if true on return, then we
717 * branch to the section that deletes a species
718 * from the current set of active species.
719 */
720 double vcs_minor_alt_calc(size_t kspec, size_t irxn, bool& do_delete) const;
721
722 //! This routine optimizes the minimization of the total Gibbs free energy
723 //! by making sure the slope of the following functional stays negative:
724 /*!
725 * The slope of the following functional is equivalent to the slope of the
726 * total Gibbs free energy of the system:
727 *
728 * d_Gibbs/ds = sum_k( m_deltaGRxn * m_deltaMolNumSpecies[k] )
729 *
730 * along the current direction m_deltaMolNumSpecies[], by choosing a value, al: (0<al<1)
731 * such that the a parabola approximation to Gibbs(al) fit to the
732 * end points al = 0 and al = 1 is minimized.
733 * - s1 = slope of Gibbs function at al = 0, which is the previous
734 * solution = d(Gibbs)/d(al).
735 * - s2 = slope of Gibbs function at al = 1, which is the current
736 * solution = d(Gibbs)/d(al).
737 *
738 * Only if there has been an inflection point (that is, s1 < 0 and s2 > 0),
739 * does this code section kick in. It finds the point on the parabola
740 * where the slope is equal to zero.
741 */
742 bool vcs_globStepDamp();
743
744 //! Calculate the norm of a deltaGibbs free energy vector
745 /*!
746 * Positive DG for species which don't exist are ignored.
747 *
748 * @param dg Vector of local delta G's.
749 */
750 double l2normdg(span<const double> dg) const;
751
752 void checkDelta1(span<const double> ds, span<const double> delTPhMoles,
753 size_t kspec);
754
755 //! Calculate the status of single species phases.
756 void vcs_SSPhase();
757
758 //! Initialize the internal counters
759 /*!
760 * Initialize the internal counters tracking subroutine calls.
761 *
762 * ifunc = 0 Initialize only those counters appropriate for the top of
763 * vcs_solve_TP().
764 * = 1 Initialize all counters.
765 */
766 void vcs_counters_init(int ifunc);
767
768 //! Create a report on the plog file containing iteration counters
770
771 void vcs_setFlagsVolPhases(const bool upToDate, const int stateCalc);
772
773 //! Update all underlying vcs_VolPhase objects
774 /*!
775 * Update the mole numbers and the phase voltages of all phases in the vcs
776 * problem
777 *
778 * @param stateCalc Location of the update (either VCS_STATECALC_NEW or
779 * VCS_STATECALC_OLD).
780 */
781 void vcs_updateMolNumVolPhases(const int stateCalc);
782
783 // Helper functions used internally by vcs_solve_TP
784 void solve_tp_component_calc(bool& allMinorZeroedSpecies);
785 void solve_tp_inner(size_t& iti, size_t& it1, bool& uptodate_minors,
786 bool& allMinorZeroedSpecies, int& forceComponentCalc,
787 int& stage, bool printDetails);
788 void solve_tp_equilib_check(bool& allMinorZeroedSpecies, bool& uptodate_minors,
789 bool& giveUpOnElemAbund, int& solveFail,
790 size_t& iti, size_t& it1, int maxit,
791 int& stage, bool& lec);
792 void solve_tp_elem_abund_check(size_t& iti, int& stage, bool& lec,
793 bool& giveUpOnElemAbund,
794 int& finalElemAbundAttempts,
795 int& rangeErrorFound);
796
797 vector<double> m_sm; //!< QR matrix work space used in vcs_basopt()
798 vector<double> m_ss; //!< Gram-Schmidt work space used in vcs_basopt()
799 vector<double> m_sa; //!< Gram-Schmidt work space used in vcs_basopt()
800 vector<double> m_aw; //!< work array of mole fractions used in vcs_basopt()
801 vector<double> m_wx;
802
803public:
804 //! Print level for print routines
806
807 MultiPhase* m_mix;
808
809 //! Print out the problem specification in all generality as it currently
810 //! exists in the VCS_SOLVE object
811 /*!
812 * @param print_lvl Parameter lvl for printing
813 * * 0 - no printing
814 * * 1 - all printing
815 */
816 void prob_report(int print_lvl);
817
818 //! Add elements to the local element list
819 //! This routines adds entries for the formula matrix for one species
820 /*!
821 * This routines adds entries for the formula matrix for this object for one
822 * species
823 *
824 * This object also fills in the index filed, IndSpecies, within the
825 * volPhase object.
826 *
827 * @param volPhase object containing the species
828 * @param k Species number within the volPhase k
829 * @param kT global Species number within this object
830 *
831 */
832 size_t addOnePhaseSpecies(vcs_VolPhase* volPhase, size_t k, size_t kT);
833
834 //! This routine resizes the number of elements in the VCS_SOLVE object by
835 //! adding a new element to the end of the element list
836 /*!
837 * The element name is added. Formula vector entries ang element abundances
838 * for the new element are set to zero.
839 *
840 * @param elNameNew New name of the element
841 * @param elType Type of the element
842 * @param elactive boolean indicating whether the element is active
843 * @returns the index number of the new element
844 */
845 size_t addElement(const char* elNameNew, int elType, int elactive);
846
847 //! Return the index of an element by name, or npos if not found
848 size_t elementIndex(const string& elementName) const;
849
850 //! Total number of species in the problems
851 size_t m_nsp;
852
853 //! Number of element constraints in the problem
854 size_t m_nelem = 0;
855
856 //! Number of components calculated for the problem
857 size_t m_numComponents = 0;
858
859 //! Total number of non-component species in the problem
860 size_t m_numRxnTot = 0;
861
862 //! Current number of species in the problems. Species can be deleted if
863 //! they aren't stable under the current conditions
865
866 //! Current number of non-component species in the problem. Species can be
867 //! deleted if they aren't stable under the current conditions
868 size_t m_numRxnRdc = 0;
869
870 //! Number of active species which are currently either treated as
871 //! minor species
873
874 //! Number of Phases in the problem
876
877 //! Formula matrix for the problem
878 /*!
879 * FormulaMatrix(kspec,j) = Number of elements, j, in the kspec species
880 *
881 * Both element and species indices are swapped.
882 */
884
885 //! Stoichiometric coefficient matrix for the reaction mechanism expressed
886 //! in Reduced Canonical Form.
887 /*!
888 * This is the stoichiometric coefficient matrix for the reaction which
889 * forms species kspec from the component species. A stoichiometric
890 * coefficient of one is assumed for the species kspec in this mechanism.
891 *
892 * NOTE: kspec = irxn + m_numComponents
893 *
894 * `m_stoichCoeffRxnMatrix(j,irxn)` : j refers to the component number, and
895 * irxn refers to the irxn_th non-component species. The stoichiometric
896 * coefficients multiplied by the Formula coefficients of the component
897 * species add up to the negative value of the number of elements in the
898 * species kspec.
899 *
900 * size = nelements0 x nspecies0
901 */
903
904 //! Absolute size of the stoichiometric coefficients
905 /*!
906 * scSize[irxn] = abs(Size) of the stoichiometric coefficients. These are
907 * used to determine whether a given species should be
908 * handled by the alt_min treatment or should be handled as a
909 * major species.
910 */
911 vector<double> m_scSize;
912
913 //! total size of the species
914 /*!
915 * This is used as a multiplier to the mole number in figuring out which
916 * species should be components.
917 */
918 vector<double> m_spSize;
919
920 //! Standard state chemical potentials for species K at the current
921 //! temperature and pressure.
922 /*!
923 * The first NC entries are for components. The following NR entries are
924 * for the current non-component species in the mechanism.
925 */
926 vector<double> m_SSfeSpecies;
927
928 //! Free energy vector from the start of the current iteration
929 /*!
930 * The free energies are saved at the start of the current iteration.
931 * Length = number of species
932 */
933 vector<double> m_feSpecies_old;
934
935 //! Dimensionless new free energy for all the species in the mechanism
936 //! at the new tentative T, P, and mole numbers.
937 /*!
938 * The first NC entries are for components. The following
939 * NR entries are for the current non-component species in the mechanism.
940 * Length = number of species
941 */
942 vector<double> m_feSpecies_new;
943
944 //! Setting for whether to do an initial estimate
945 /*!
946 * Initial estimate:
947 * 0 Do not estimate the solution at all. Use the
948 * supplied mole numbers as is.
949 * 1 Only do an estimate if the element abundances aren't satisfied.
950 * -1 Force an estimate of the soln. Throw out the input mole numbers.
951 */
953
954 //! Total moles of the species
955 /*!
956 * Total number of moles of the kth species.
957 * Length = Total number of species = m
958 */
959 vector<double> m_molNumSpecies_old;
960
961 //! Specifies the species unknown type
962 /*!
963 * There are two types. One is the straightforward species, with the mole
964 * number w[k], as the unknown. The second is the an interfacial voltage
965 * where w[k] refers to the interfacial voltage in volts.
966 *
967 * These species types correspond to metallic electrons corresponding to
968 * electrodes. The voltage and other interfacial conditions sets up an
969 * interfacial current, which is set to zero in this initial treatment.
970 * Later we may have non-zero interfacial currents.
971 */
973
974 //! Change in the number of moles of phase, iphase, due to the
975 //! noncomponent formation reaction, irxn, for species, k:
976 /*!
977 * m_deltaMolNumPhase(iphase,irxn) = k = nc + irxn
978 */
980
981 //! This is 1 if the phase, iphase, participates in the formation reaction
982 //! irxn, and zero otherwise. PhaseParticipation(iphase,irxn)
984
985 //! electric potential of the iph phase
986 vector<double> m_phasePhi;
987
988 //! Tentative value of the mole number vector. It's also used to store the
989 //! mole fraction vector.
990 vector<double> m_molNumSpecies_new;
991
992 //! Delta G(irxn) for the noncomponent species in the mechanism.
993 /*!
994 * Computed by the subroutine deltaG. m_deltaGRxn is the free energy change
995 * for the reaction which forms species K from the component species. This
996 * vector has length equal to the number of noncomponent species in the
997 * mechanism. It starts with the first current noncomponent species in the
998 * mechanism.
999 */
1000 vector<double> m_deltaGRxn_new;
1001
1002 //! Last deltag[irxn] from the previous step
1003 vector<double> m_deltaGRxn_old;
1004
1005 //! Last deltag[irxn] from the previous step with additions for
1006 //! possible births of zeroed phases.
1008
1009 //! Temporary vector of Rxn DeltaG's
1010 /*!
1011 * This is used from time to time, for printing purposes
1012 */
1013 vector<double> m_deltaGRxn_tmp;
1014
1015 //! Reaction Adjustments for each species during the current step
1016 /*!
1017 * delta Moles for each species during the current step.
1018 * Length = number of species
1019 */
1020 vector<double> m_deltaMolNumSpecies;
1021
1022 //! Element abundances vector
1023 /*!
1024 * Vector of moles of each element actually in the solution vector. Except
1025 * for certain parts of the algorithm, this is a constant. Note other
1026 * constraint conditions are added to this vector. This is input from the
1027 * input file and is considered a constant from thereon. units = kmoles
1028 */
1029 vector<double> m_elemAbundances;
1030
1031 //! Element abundances vector Goals
1032 /*!
1033 * Vector of moles of each element that are the goals of the simulation.
1034 * This is a constant in the problem. Note other constraint conditions are
1035 * added to this vector. This is input from the input file and is considered
1036 * a constant from thereon. units = kmoles
1037 */
1038 vector<double> m_elemAbundancesGoal;
1039
1040 //! Total number of kmoles in all phases
1041 /*!
1042 * Don't use this except for scaling purposes
1043 */
1044 double m_totalMolNum = 0.0;
1045
1046 //! Total kmols of species in each phase
1047 /*!
1048 * This contains the total number of moles of species in each phase
1049 *
1050 * Length = number of phases
1051 */
1052 vector<double> m_tPhaseMoles_old;
1053
1054 //! total kmols of species in each phase in the tentative soln vector
1055 /*!
1056 * This contains the total number of moles of species in each phase
1057 * in the tentative solution vector
1058 *
1059 * Length = number of phases
1060 */
1061 vector<double> m_tPhaseMoles_new;
1062
1063 //! Temporary vector of length NPhase
1064 mutable vector<double> m_TmpPhase;
1065
1066 //! Temporary vector of length NPhase
1067 mutable vector<double> m_TmpPhase2;
1068
1069 //! Change in the total moles in each phase
1070 /*!
1071 * Length number of phases.
1072 */
1073 vector<double> m_deltaPhaseMoles;
1074
1075 //! Temperature (Kelvin)
1077
1078 //! Pressure
1080
1081 //! Tolerance requirement for major species
1082 double m_tolmaj= 1e-8;
1083
1084 //! Tolerance requirements for minor species
1085 double m_tolmin = 1e-6;
1086
1087 //! Below this, major species aren't refined any more
1088 double m_tolmaj2 = 1e-10;
1089
1090 //! Below this, minor species aren't refined any more
1091 double m_tolmin2 = 1e-8;
1092
1093 //! Index vector that keeps track of the species vector rearrangement
1094 /*!
1095 * At the end of each run, the species vector and associated data gets put
1096 * back in the original order.
1097 *
1098 * Example
1099 *
1100 * k = m_speciesMapIndex[kspec]
1101 *
1102 * kspec = current order in the vcs_solve object
1103 * k = original order in the MultiPhase object
1104 */
1105 vector<size_t> m_speciesMapIndex;
1106
1107 //! Index that keeps track of the index of the species within the local
1108 //! phase
1109 /*!
1110 * This returns the local index of the species within the phase. Its
1111 * argument is the global species index within the VCS problem.
1112 *
1113 * k = m_speciesLocalPhaseIndex[kspec]
1114 *
1115 * k varies between 0 and the nSpecies in the phase
1116 *
1117 * Length = number of species
1118 */
1120
1121 //! Index vector that keeps track of the rearrangement of the elements
1122 /*!
1123 * At the end of each run, the element vector and associated data gets
1124 * put back in the original order.
1125 *
1126 * Example
1127 *
1128 * e = m_elementMapIndex[eNum]
1129 * eNum = current order in the vcs_solve object
1130 * e = original order in the MultiPhase object
1131 */
1132 vector<size_t> m_elementMapIndex;
1133
1134 //! Mapping between the species index for noncomponent species and the
1135 //! full species index.
1136 /*!
1137 * ir[irxn] = Mapping between the reaction index for noncomponent
1138 * formation reaction of a species and the full species
1139 * index.
1140 *
1141 * Initially set to a value of K = NC + I This vector has length equal to
1142 * number of noncomponent species in the mechanism. It starts with the
1143 * first current noncomponent species in the mechanism. kspec = ir[irxn]
1144 */
1145 vector<size_t> m_indexRxnToSpecies;
1146
1147 //! Major -Minor status vector for the species in the problem
1148 /*!
1149 * The index for this is species. The reaction that this is referring to
1150 * is `kspec = irxn + m_numComponents`. For possible values and their
1151 * meanings, see vcs_evaluate_speciesType().
1152 */
1153 vector<int> m_speciesStatus;
1154
1155 //! Mapping from the species number to the phase number
1156 vector<size_t> m_phaseID;
1157
1158 //! Boolean indicating whether a species belongs to a single-species phase
1159 // vector<bool> can't be used here because it doesn't work with std::swap
1160 vector<char> m_SSPhase;
1161
1162 //! Species string name for the kth species
1163 vector<string> m_speciesName;
1164
1165 //! Vector of strings containing the element names
1166 /*!
1167 * ElName[j] = String containing element names
1168 */
1169 vector<string> m_elementName;
1170
1171 //! Type of the element constraint
1172 /*!
1173 * * 0 - #VCS_ELEM_TYPE_ABSPOS Normal element that is positive or zero in
1174 * all species.
1175 * * 1 - #VCS_ELEM_TYPE_ELECTRONCHARGE element dof that corresponds to
1176 * the electronic charge DOF.
1177 * * 2 - #VCS_ELEM_TYPE_CHARGENEUTRALITY element dof that corresponds to
1178 * a required charge neutrality constraint on the phase. The element
1179 * abundance is always exactly zero.
1180 */
1181 vector<int> m_elType;
1182
1183 //! Specifies whether an element constraint is active
1184 /*!
1185 * The default is true. Length = nelements
1186 */
1187 vector<int> m_elementActive;
1188
1189 //! Array of Phase Structures. Length = number of phases.
1190 vector<unique_ptr<vcs_VolPhase>> m_VolPhaseList;
1191
1192 //! specifies the activity convention of the phase containing the species
1193 /*!
1194 * * 0 = molar based
1195 * * 1 = molality based
1196 *
1197 * length = number of species
1198 */
1200
1201 //! specifies the activity convention of the phase.
1202 /*!
1203 * * 0 = molar based
1204 * * 1 = molality based
1205 *
1206 * length = number of phases
1207 */
1209
1210 //! specifies the ln(Mnaught) used to calculate the chemical potentials
1211 /*!
1212 * For molar based activity conventions this will be equal to 0.0.
1213 * length = number of species.
1214 */
1215 vector<double> m_lnMnaughtSpecies;
1216
1217 //! Molar-based Activity Coefficients for Species.
1218 //! Length = number of species
1220
1221 //! Molar-based Activity Coefficients for Species based on old mole numbers
1222 /*!
1223 * These activity coefficients are based on the m_molNumSpecies_old
1224 * values Molar based activity coefficients. Length = number of species
1225 */
1227
1228 //! Change in the log of the activity coefficient with respect to the mole
1229 //! number multiplied by the phase mole number
1230 /*!
1231 * size = nspecies x nspecies
1232 *
1233 * This is a temporary array that gets regenerated every time it's needed.
1234 * It is not swapped wrt species.
1235 */
1237
1238 //! Molecular weight of each species
1239 /*!
1240 * units = kg/kmol. length = number of species.
1241 *
1242 * note: this is a candidate for removal. I don't think we use it.
1243 */
1244 vector<double> m_wtSpecies;
1245
1246 //! Charge of each species. Length = number of species.
1247 vector<double> m_chargeSpecies;
1248
1249 vector<vector<size_t>> phasePopProblemLists_;
1250
1251 //! Choice of Hessians
1252 /*!
1253 * If this is true, then we will use a better approximation to the Hessian
1254 * based on Jacobian of the ln(ActCoeff) with respect to mole numbers
1255 */
1257
1258 //! Partial molar volumes of the species
1259 /*!
1260 * units = mks (m^3/kmol)
1261 * Length = number of species
1262 */
1263 vector<double> m_PMVolumeSpecies;
1264
1265 //! dimensionless value of Faraday's constant, F / RT (1/volt)
1267
1268 //! Timing and iteration counters for the vcs object
1270
1271 //! Debug printing lvl
1272 /*!
1273 * Levels correspond to the following guidelines
1274 * * 0 No printing at all
1275 * * 1 Serious warnings or fatal errors get one line
1276 * * 2 one line per each successful vcs package call
1277 * * 3 one line per every successful solve_TP calculation
1278 * * 4 one line for every successful operation -> solve_TP gets a summary report
1279 * * 5 each iteration in solve_TP gets a report with one line per species
1280 * * 6 Each decision in solve_TP gets a line per species in addition to 4
1281 * * 10 Additionally Hessian matrix is printed out
1282 */
1284
1285 friend class vcs_phaseStabilitySolve;
1286};
1287
1288}
1289#endif
Header file for class Cantera::Array2D.
A class for 2D arrays stored in column-major (Fortran-compatible) form.
Definition Array.h:32
A class for multiphase mixtures.
Definition MultiPhase.h:62
Class to keep track of time and iterations.
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:44
vector< double > m_deltaPhaseMoles
Change in the total moles in each phase.
Definition vcs_solve.h:1073
size_t m_numPhases
Number of Phases in the problem.
Definition vcs_solve.h:875
int m_doEstimateEquil
Setting for whether to do an initial estimate.
Definition vcs_solve.h:952
vector< double > m_deltaMolNumSpecies
Reaction Adjustments for each species during the current step.
Definition vcs_solve.h:1020
vector< double > m_spSize
total size of the species
Definition vcs_solve.h:918
vector< double > m_scSize
Absolute size of the stoichiometric coefficients.
Definition vcs_solve.h:911
int vcs_evalSS_TP(int ipr, int ip1, double Temp, double pres)
Evaluate the standard state free energies at the current temperature and pressure.
vector< double > m_sm
QR matrix work space used in vcs_basopt()
Definition vcs_solve.h:797
vector< double > m_wtSpecies
Molecular weight of each species.
Definition vcs_solve.h:1244
vector< size_t > m_indexRxnToSpecies
Mapping between the species index for noncomponent species and the full species index.
Definition vcs_solve.h:1145
void vcs_prep(int printLvl)
This routine is mostly concerned with changing the private data to be consistent with what's needed f...
vector< double > m_sa
Gram-Schmidt work space used in vcs_basopt()
Definition vcs_solve.h:799
double l2normdg(span< const double > dg) const
Calculate the norm of a deltaGibbs free energy vector.
vector< size_t > m_phaseID
Mapping from the species number to the phase number.
Definition vcs_solve.h:1156
vector< string > m_speciesName
Species string name for the kth species.
Definition vcs_solve.h:1163
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.
double m_tolmaj2
Below this, major species aren't refined any more.
Definition vcs_solve.h:1088
vector< double > m_deltaGRxn_new
Delta G(irxn) for the noncomponent species in the mechanism.
Definition vcs_solve.h:1000
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:1236
void vcs_updateMolNumVolPhases(const int stateCalc)
Update all underlying vcs_VolPhase objects.
void vcs_report(int iconv)
Print out a report on the state of the equilibrium problem to standard output.
size_t m_nelem
Number of element constraints in the problem.
Definition vcs_solve.h:854
vector< double > m_PMVolumeSpecies
Partial molar volumes of the species.
Definition vcs_solve.h:1263
vector< int > m_elementActive
Specifies whether an element constraint is active.
Definition vcs_solve.h:1187
int m_useActCoeffJac
Choice of Hessians.
Definition vcs_solve.h:1256
void vcs_CalcLnActCoeffJac(span< const double > moleSpeciesVCS)
Recalculate all of the activity coefficients in all of the phases based on input mole numbers.
vector< size_t > m_speciesLocalPhaseIndex
Index that keeps track of the index of the species within the local phase.
Definition vcs_solve.h:1119
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...
size_t vcs_add_all_deleted()
Provide an estimate for the deleted species in phases that are not zeroed out.
bool vcs_globStepDamp()
This routine optimizes the minimization of the total Gibbs free energy by making sure the slope of th...
void vcs_basopt(const bool doJustComponents, double test)
Choose the optimum species basis for the calculations.
vector< string > m_elementName
Vector of strings containing the element names.
Definition vcs_solve.h:1169
vector< int > m_speciesUnknownType
Specifies the species unknown type.
Definition vcs_solve.h:972
vector< double > m_TmpPhase
Temporary vector of length NPhase.
Definition vcs_solve.h:1064
vector< size_t > m_speciesMapIndex
Index vector that keeps track of the species vector rearrangement.
Definition vcs_solve.h:1105
Array2D m_phaseParticipation
This is 1 if the phase, iphase, participates in the formation reaction irxn, and zero otherwise.
Definition vcs_solve.h:983
vector< int > m_actConventionSpecies
specifies the activity convention of the phase containing the species
Definition vcs_solve.h:1199
double vcs_phaseStabilityTest(const size_t iph)
Main program to test whether a deleted phase should be brought back into existence.
void vcs_inest()
Create an initial estimate of the solution to the equilibrium problem.
vector< double > m_elemAbundances
Element abundances vector.
Definition vcs_solve.h:1029
vector< double > m_tPhaseMoles_new
total kmols of species in each phase in the tentative soln vector
Definition vcs_solve.h:1061
vector< double > m_actCoeffSpecies_old
Molar-based Activity Coefficients for Species based on old mole numbers.
Definition vcs_solve.h:1226
vector< double > m_SSfeSpecies
Standard state chemical potentials for species K at the current temperature and pressure.
Definition vcs_solve.h:926
bool vcs_popPhasePossible(const size_t iphasePop) const
Utility function that evaluates whether a phase can be popped into existence.
int 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.
double m_tolmin2
Below this, minor species aren't refined any more.
Definition vcs_solve.h:1091
vector< double > m_lnMnaughtSpecies
specifies the ln(Mnaught) used to calculate the chemical potentials
Definition vcs_solve.h:1215
vector< int > m_elType
Type of the element constraint.
Definition vcs_solve.h:1181
void vcs_elab()
Computes the current elemental abundances vector.
void prob_report(int print_lvl)
Print out the problem specification in all generality as it currently exists in the VCS_SOLVE object.
double vcs_minor_alt_calc(size_t kspec, size_t irxn, bool &do_delete) const
Minor species alternative calculation.
int vcs_setMolesLinProg()
Estimate the initial mole numbers by constrained linear programming.
vector< double > m_molNumSpecies_old
Total moles of the species.
Definition vcs_solve.h:959
size_t m_numRxnMinorZeroed
Number of active species which are currently either treated as minor species.
Definition vcs_solve.h:872
vector< size_t > m_elementMapIndex
Index vector that keeps track of the rearrangement of the elements.
Definition vcs_solve.h:1132
vector< double > m_aw
work array of mole fractions used in vcs_basopt()
Definition vcs_solve.h:800
size_t m_numSpeciesRdc
Current number of species in the problems.
Definition vcs_solve.h:864
int vcs_delete_species(const size_t kspec)
Change a single species from active to inactive status.
vector< double > m_feSpecies_old
Free energy vector from the start of the current iteration.
Definition vcs_solve.h:933
int m_debug_print_lvl
Debug printing lvl.
Definition vcs_solve.h:1283
void vcs_counters_init(int ifunc)
Initialize the internal counters.
void vcs_reinsert_deleted(size_t kspec)
We make decisions on the initial mole number, and major-minor status here.
vector< double > m_ss
Gram-Schmidt work space used in vcs_basopt()
Definition vcs_solve.h:798
vector< double > m_chargeSpecies
Charge of each species. Length = number of species.
Definition vcs_solve.h:1247
Array2D m_formulaMatrix
Formula matrix for the problem.
Definition vcs_solve.h:883
size_t addOnePhaseSpecies(vcs_VolPhase *volPhase, size_t k, size_t kT)
Add elements to the local element list This routines adds entries for the formula matrix for one spec...
Array2D m_stoichCoeffRxnMatrix
Stoichiometric coefficient matrix for the reaction mechanism expressed in Reduced Canonical Form.
Definition vcs_solve.h:902
size_t m_numRxnRdc
Current number of non-component species in the problem.
Definition vcs_solve.h:868
VCS_COUNTERS * m_VCount
Timing and iteration counters for the vcs object.
Definition vcs_solve.h:1269
vector< double > m_deltaGRxn_Deficient
Last deltag[irxn] from the previous step with additions for possible births of zeroed phases.
Definition vcs_solve.h:1007
vector< unique_ptr< vcs_VolPhase > > m_VolPhaseList
Array of Phase Structures. Length = number of phases.
Definition vcs_solve.h:1190
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...
size_t m_numComponents
Number of components calculated for the problem.
Definition vcs_solve.h:857
void vcs_prob_specifyFully()
Fully specify the problem to be solved.
int vcs_species_type(const size_t kspec) const
Evaluate the species category for the indicated species.
double m_pressurePA
Pressure.
Definition vcs_solve.h:1079
void vcs_SSPhase()
Calculate the status of single species phases.
vector< double > m_feSpecies_new
Dimensionless new free energy for all the species in the mechanism at the new tentative T,...
Definition vcs_solve.h:942
double m_Faraday_dim
dimensionless value of Faraday's constant, F / RT (1/volt)
Definition vcs_solve.h:1266
int vcs_recheck_deleted()
Recheck deleted species in multispecies phases.
void vcs_elem_rearrange()
Rearrange the constraint equations represented by the Formula Matrix so that the operational ones are...
double m_tolmin
Tolerance requirements for minor species.
Definition vcs_solve.h:1085
size_t m_nsp
Total number of species in the problems.
Definition vcs_solve.h:851
vector< double > m_actCoeffSpecies_new
Molar-based Activity Coefficients for Species.
Definition vcs_solve.h:1219
bool vcs_elabcheck(int ibound)
Checks to see if the element abundances are in compliance.
vector< char > m_SSPhase
Boolean indicating whether a species belongs to a single-species phase.
Definition vcs_solve.h:1160
bool vcs_delete_multiphase(const size_t iph)
This routine handles the bookkeeping involved with the deletion of multiphase phases from the problem...
vector< double > m_phasePhi
electric potential of the iph phase
Definition vcs_solve.h:986
double vcs_Hessian_actCoeff_diag(size_t irxn)
Calculates the diagonal contribution to the Hessian due to the dependence of the activity coefficient...
int delta_species(const size_t kspec, double &delta)
Change the concentration of a species by delta moles.
int vcs_zero_species(const size_t kspec)
Zero out the concentration of a species.
vector< double > m_elemAbundancesGoal
Element abundances vector Goals.
Definition vcs_solve.h:1038
size_t vcs_RxnStepSizes(int &forceComponentCalc, size_t &kSpecial)
Calculates formation reaction step sizes.
size_t vcs_popPhaseID()
Decision as to whether a phase pops back into existence.
Array2D m_deltaMolNumPhase
Change in the number of moles of phase, iphase, due to the noncomponent formation reaction,...
Definition vcs_solve.h:979
vector< int > m_phaseActConvention
specifies the activity convention of the phase.
Definition vcs_solve.h:1208
vector< int > m_speciesStatus
Major -Minor status vector for the species in the problem.
Definition vcs_solve.h:1153
vector< double > m_deltaGRxn_tmp
Temporary vector of Rxn DeltaG's.
Definition vcs_solve.h:1013
vector< double > m_TmpPhase2
Temporary vector of length NPhase.
Definition vcs_solve.h:1067
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.
double vcs_tmoles()
Calculates the total number of moles of species in all phases.
size_t elementIndex(const string &elementName) const
Return the index of an element by name, or npos if not found.
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.
int m_printLvl
Print level for print routines.
Definition vcs_solve.h:805
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,...
vector< double > m_tPhaseMoles_old
Total kmols of species in each phase.
Definition vcs_solve.h:1052
double vcs_Total_Gibbs(span< const double > w, span< const double > fe, span< const double > tPhMoles)
Calculate the total dimensionless Gibbs free energy.
int vcs_popPhaseRxnStepSizes(const size_t iphasePop)
Calculates the deltas of the reactions due to phases popping into existence.
size_t m_numRxnTot
Total number of non-component species in the problem.
Definition vcs_solve.h:860
vector< double > m_molNumSpecies_new
Tentative value of the mole number vector.
Definition vcs_solve.h:990
double m_temperature
Temperature (Kelvin)
Definition vcs_solve.h:1076
double m_tolmaj
Tolerance requirement for major species.
Definition vcs_solve.h:1082
vector< double > m_deltaGRxn_old
Last deltag[irxn] from the previous step.
Definition vcs_solve.h:1003
bool vcs_evaluate_speciesType()
This routine evaluates the species type for all species.
void vcs_TCounters_report()
Create a report on the plog file containing iteration counters.
int vcs_elcorr(span< double > aa, span< double > x)
This subroutine corrects for element abundances.
size_t addElement(const char *elNameNew, int elType, int elactive)
This routine resizes the number of elements in the VCS_SOLVE object by adding a new element to the en...
double m_totalMolNum
Total number of kmoles in all phases.
Definition vcs_solve.h:1044
Phase information and Phase calculations for vcs.
This file contains definitions of constants, types and terms that are used in internal routines and a...
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595
Defines and definitions within the vcs package.
Internal declarations for the VCSnonideal package.