Cantera  4.0.0a1
Loading...
Searching...
No Matches
vcs_solve.cpp
Go to the documentation of this file.
1/*!
2 * @file vcs_solve.cpp Implementation file for the internal class that holds
3 * the problem.
4 */
5
6// This file is part of Cantera. See License.txt in the top-level directory or
7// at https://cantera.org/license.txt for license and copyright information.
8
16
17using namespace std;
18
19namespace Cantera
20{
21
22namespace {
23
24void printProgress(const vector<string>& spName, const vector<double>& soln,
25 const vector<double>& ff)
26{
27 double sum = 0.0;
28 plogf(" --- Summary of current progress:\n");
29 plogf(" --- Name Moles - SSGibbs \n");
30 plogf(" -------------------------------------------------------------------------------------\n");
31 for (size_t k = 0; k < soln.size(); k++) {
32 plogf(" --- %20s %12.4g - %12.4g\n", spName[k], soln[k], ff[k]);
33 sum += soln[k] * ff[k];
34 }
35 plogf(" --- Total sum to be minimized = %g\n", sum);
36}
37
38} // anonymous namespace
39
40VCS_SOLVE::VCS_SOLVE(MultiPhase* mphase, int printLvl) :
41 m_printLvl(printLvl),
42 m_mix(mphase),
43 m_nsp(mphase->nSpecies()),
44 m_numSpeciesRdc(mphase->nSpecies()),
45 m_numPhases(mphase->nPhases()),
46 m_temperature(mphase->temperature()),
47 m_pressurePA(mphase->pressure()),
48 m_Faraday_dim(Faraday / (m_temperature * GasConstant))
49{
50 string ser = "VCS_SOLVE: ERROR:\n\t";
51 if (m_nsp <= 0) {
52 plogf("%s Number of species is nonpositive\n", ser);
53 throw CanteraError("VCS_SOLVE::VCS_SOLVE", ser +
54 " Number of species is nonpositive\n");
55 }
56 if (m_numPhases <= 0) {
57 plogf("%s Number of phases is nonpositive\n", ser);
58 throw CanteraError("VCS_SOLVE::VCS_SOLVE", ser +
59 " Number of phases is nonpositive\n");
60 }
61
62 /*
63 * We will initialize sc[] to note the fact that it needs to be
64 * filled with meaningful information.
65 */
66 m_scSize.resize(m_nsp, 0.0);
67 m_spSize.resize(m_nsp, 1.0);
68 m_SSfeSpecies.resize(m_nsp, 0.0);
69 m_feSpecies_new.resize(m_nsp, 0.0);
70 m_molNumSpecies_old.resize(m_nsp, 0.0);
74 m_phasePhi.resize(m_numPhases, 0.0);
75 m_molNumSpecies_new.resize(m_nsp, 0.0);
76 m_deltaGRxn_new.resize(m_nsp, 0.0);
77 m_deltaGRxn_old.resize(m_nsp, 0.0);
78 m_deltaGRxn_Deficient.resize(m_nsp, 0.0);
79 m_deltaGRxn_tmp.resize(m_nsp, 0.0);
80 m_deltaMolNumSpecies.resize(m_nsp, 0.0);
81 m_feSpecies_old.resize(m_nsp, 0.0);
82 m_tPhaseMoles_old.resize(m_numPhases, 0.0);
83 m_tPhaseMoles_new.resize(m_numPhases, 0.0);
84 m_deltaPhaseMoles.resize(m_numPhases, 0.0);
85 m_TmpPhase.resize(m_numPhases, 0.0);
86 m_TmpPhase2.resize(m_numPhases, 0.0);
87
88 // ind[] is an index variable that keep track of solution vector rotations.
89 m_speciesMapIndex.resize(m_nsp, 0);
91
92 // ir[] is an index vector that keeps track of the irxn to species mapping.
93 // We can't fill it in until we know the number of c components in the
94 // problem
95 m_indexRxnToSpecies.resize(m_nsp, 0);
96
97 // Initialize all species to be major species
99
100 m_SSPhase.resize(2*m_nsp, 0);
101 m_phaseID.resize(m_nsp, 0);
102 m_speciesName.resize(m_nsp);
103
104 // space for activity coefficients for all species. Set it equal to one.
105 m_actConventionSpecies.resize(m_nsp, 0);
107 m_lnMnaughtSpecies.resize(m_nsp, 0.0);
108 m_actCoeffSpecies_new.resize(m_nsp, 1.0);
109 m_actCoeffSpecies_old.resize(m_nsp, 1.0);
110 m_wtSpecies.resize(m_nsp, 0.0);
111 m_chargeSpecies.resize(m_nsp, 0.0);
112
113 // Phase Info
115
116 // For Future expansion
117 m_useActCoeffJac = true;
118 if (m_useActCoeffJac) {
120 }
121
122 m_PMVolumeSpecies.resize(m_nsp, 0.0);
123
124 // counters kept within vcs
125 m_VCount = new VCS_COUNTERS();
127
128 // Loop over the phases, transferring pertinent information
129 int kT = 0;
130 for (size_t iphase = 0; iphase < m_numPhases; iphase++) {
131 // Get the ThermoPhase object - assume volume phase
132 ThermoPhase* tPhase = &mphase->phase(iphase);
133
134 // Assign the value of eqn of state. Handle conflicts here.
135 if (tPhase->nDim() != 3) {
136 throw CanteraError("VCS_SOLVE::VCS_SOLVE",
137 "Surface/edge phase not handled yet.");
138 }
139
140 m_VolPhaseList.emplace_back(
141 make_unique<vcs_VolPhase>(this, tPhase, iphase));
142 vcs_VolPhase* VolPhase = m_VolPhaseList[iphase].get();
144
145 // Loop through each species in the current phase
146 size_t nSpPhase = tPhase->nSpecies();
147 for (size_t k = 0; k < nSpPhase; k++) {
148 // Obtain the molecular weight of the species from the
149 // ThermoPhase object
150 m_wtSpecies[kT] = tPhase->molecularWeight(k);
151
152 // Obtain the charges of the species from the ThermoPhase object
153 m_chargeSpecies[kT] = tPhase->charge(k);
154
155 // Set the phaseid of the species
156 m_phaseID[kT] = iphase;
157
158 // Transfer the type of unknown
159 m_speciesUnknownType[kT] = VolPhase->speciesUnknownType(k);
161 // Set the initial number of kmoles of the species
162 m_molNumSpecies_old[kT] = mphase->speciesMoles(kT);
165 } else {
166 throw CanteraError("VCS_SOLVE::VCS_SOLVE",
167 "Unknown species type: {}", m_speciesUnknownType[kT]);
168 }
169
170 // Transfer the species information from the
171 // volPhase structure to the VPROB structure
172 // This includes:
173 // FormulaMatrix[][]
174 // VolPhase->IndSpecies[]
175 addOnePhaseSpecies(VolPhase, k, kT);
176
177 kT++;
178 }
179
181 }
182
183 // Work arrays used by vcs_basopt
184 m_sm.assign(m_nelem * m_nelem, 0.0);
185 m_ss.assign(m_nelem, 0.0);
186 m_sa.assign(m_nelem, 0.0);
187 m_aw.assign(m_nsp, 0.0);
188 m_wx.assign(m_nelem, 0.0);
189
190 // Transfer initial element abundances based on the species mole numbers
191 for (size_t j = 0; j < m_nelem; j++) {
192 for (size_t kspec = 0; kspec < m_nsp; kspec++) {
195 }
196 }
197 }
198
199 // Printout the species information: PhaseID's and mole nums
200 if (m_printLvl > 1) {
201 writeline('=', 80, true, true);
202 writeline('=', 16, false);
203 plogf(" Cantera_to_vprob: START OF PROBLEM STATEMENT ");
204 writeline('=', 20);
205 writeline('=', 80);
206 plogf(" Phase IDs of species\n");
207 plogf(" species phaseID phaseName ");
208 plogf(" Initial_Estimated_kMols\n");
209 for (size_t i = 0; i < m_nsp; i++) {
210 size_t iphase = m_phaseID[i];
211
212 vcs_VolPhase* VolPhase = m_VolPhaseList[iphase].get();
213 plogf("%16s %5d %16s", mphase->speciesName(i).c_str(), iphase,
214 VolPhase->PhaseName.c_str());
216 plogf(" Volts = %-10.5g\n", m_molNumSpecies_old[i]);
217 } else {
218 plogf(" %-10.5g\n", m_molNumSpecies_old[i]);
219 }
220 }
221
222 // Printout of the Phase structure information
223 writeline('-', 80, true, true);
224 plogf(" Information about phases\n");
225 plogf(" PhaseName PhaseNum SingSpec EqnState NumSpec "
226 "Tmoles(kmol)\n");
227
228 for (size_t iphase = 0; iphase < m_numPhases; iphase++) {
229 vcs_VolPhase* VolPhase = m_VolPhaseList[iphase].get();
230 plogf("%16s %8d %8d %16s %8d ", VolPhase->PhaseName.c_str(),
231 VolPhase->VP_ID_, VolPhase->m_singleSpecies,
232 VolPhase->eos_name(), VolPhase->nSpecies());
233 plogf("%16e\n", VolPhase->totalMoles());
234 }
235
236 writeline('=', 80, true, true);
237 writeline('=', 16, false);
238 plogf(" Cantera_to_vprob: END OF PROBLEM STATEMENT ");
239 writeline('=', 20);
240 writeline('=', 80);
241 plogf("\n");
242 }
243
244 // m_speciesIndexVector[] is an index variable that keep track of solution
245 // vector rotations.
246 for (size_t i = 0; i < m_nsp; i++) {
247 m_speciesMapIndex[i] = i;
248 }
249
250 // IndEl[] is an index variable that keep track of element vector rotations.
251 for (size_t i = 0; i < m_nelem; i++) {
252 m_elementMapIndex[i] = i;
253 }
254
255 // Fill in the species to phase mapping. Check for bad values at the same
256 // time.
257 vector<size_t> numPhSp(m_numPhases, 0);
258 for (size_t kspec = 0; kspec < m_nsp; kspec++) {
259 size_t iph = m_phaseID[kspec];
260 if (iph >= m_numPhases) {
261 throw CanteraError("VCS_SOLVE::VCS_SOLVE",
262 "Species to Phase Mapping, PhaseID, has a bad value\n"
263 "\tm_phaseID[{}] = {}\n"
264 "Allowed values: 0 to {}", kspec, iph, m_numPhases - 1);
265 }
266 m_phaseID[kspec] = m_phaseID[kspec];
267 m_speciesLocalPhaseIndex[kspec] = numPhSp[iph];
268 numPhSp[iph]++;
269 }
270 for (size_t iph = 0; iph < m_numPhases; iph++) {
271 vcs_VolPhase* Vphase = m_VolPhaseList[iph].get();
272 if (numPhSp[iph] != Vphase->nSpecies()) {
273 throw CanteraError("VCS_SOLVE::VCS_SOLVE",
274 "Number of species in phase {}, {}, doesn't match ({} != {}) [vphase = {}]",
275 ser, iph, Vphase->PhaseName, numPhSp[iph], Vphase->nSpecies(), (size_t) Vphase);
276 }
277 }
278
279 for (size_t i = 0; i < m_nelem; i++) {
281 if (m_elemAbundancesGoal[i] != 0.0) {
282 if (fabs(m_elemAbundancesGoal[i]) > 1.0E-9) {
283 throw CanteraError("VCS_SOLVE::VCS_SOLVE",
284 "Charge neutrality condition {} is signicantly "
285 "nonzero, {}. Giving up",
287 } else {
288 if (m_debug_print_lvl >= 2) {
289 plogf("Charge neutrality condition %s not zero, %g. Setting it zero\n",
291 }
292 m_elemAbundancesGoal[i] = 0.0;
293 }
294 }
295 }
296 }
297
298 // Copy over the species names
299 for (size_t i = 0; i < m_nsp; i++) {
300 m_speciesName[i] = m_mix->speciesName(i);
301 }
302
303 // Specify the Activity Convention information
304 for (size_t iph = 0; iph < m_numPhases; iph++) {
305 vcs_VolPhase* Vphase = m_VolPhaseList[iph].get();
307 if (Vphase->p_activityConvention != 0) {
308 // We assume here that species 0 is the solvent. The solvent isn't
309 // on a unity activity basis The activity for the solvent assumes
310 // that the it goes to one as the species mole fraction goes to one;
311 // that is, it's really on a molarity framework. So
312 // SpecLnMnaught[iSolvent] = 0.0, and the loop below starts at 1,
313 // not 0.
314 size_t iSolvent = Vphase->spGlobalIndexVCS(0);
315 double mnaught = m_wtSpecies[iSolvent] / 1000.;
316 for (size_t k = 1; k < Vphase->nSpecies(); k++) {
317 size_t kspec = Vphase->spGlobalIndexVCS(k);
319 m_lnMnaughtSpecies[kspec] = log(mnaught);
320 }
321 }
322 }
323
324}
325
326VCS_SOLVE::~VCS_SOLVE()
327{
328 delete m_VCount;
329}
330
331bool VCS_SOLVE::vcs_popPhasePossible(const size_t iphasePop) const
332{
333 vcs_VolPhase* Vphase = m_VolPhaseList[iphasePop].get();
334 AssertThrowMsg(!Vphase->exists(), "VCS_SOLVE::vcs_popPhasePossible",
335 "called for a phase that exists!");
336
337 // Loop through all of the species in the phase. We say the phase can be
338 // popped, if there is one species in the phase that can be popped. This
339 // does not mean that the phase will be popped or that it leads to a lower
340 // Gibbs free energy.
341 for (size_t k = 0; k < Vphase->nSpecies(); k++) {
342 size_t kspec = Vphase->spGlobalIndexVCS(k);
344 "VCS_SOLVE::vcs_popPhasePossible",
345 "we shouldn't be here {}: {} > 0.0", kspec,
346 m_molNumSpecies_old[kspec]);
347 size_t irxn = kspec - m_numComponents;
348 if (kspec >= m_numComponents) {
349 bool iPopPossible = true;
350
351 // Note one case is if the component is a member of the popping
352 // phase. This component will be zeroed and the logic here will
353 // negate the current species from causing a positive if this
354 // component is consumed.
355 for (size_t j = 0; j < m_numComponents; ++j) {
356 if (m_elType[j] == VCS_ELEM_TYPE_ABSPOS) {
357 double stoicC = m_stoichCoeffRxnMatrix(j,irxn);
358 if (stoicC != 0.0) {
359 double negChangeComp = - stoicC;
360 if (negChangeComp > 0.0) {
361 // If there is no component to give, then the
362 // species can't be created
364 iPopPossible = false;
365 }
366 }
367 }
368 }
369 }
370 // We are here when the species can be popped because all its needed
371 // components have positive mole numbers
372 if (iPopPossible) {
373 return true;
374 }
375 } else {
376 // We are here when the species, k, in the phase is a component. Its
377 // mole number is zero. We loop through the regular reaction looking
378 // for a reaction that can pop the component.
379 for (size_t jrxn = 0; jrxn < m_numRxnRdc; jrxn++) {
380 bool foundJrxn = false;
381 // First, if the component is a product of the reaction
382 if (m_stoichCoeffRxnMatrix(kspec,jrxn) > 0.0) {
383 foundJrxn = true;
384 // We can do the reaction if all other reactant components
385 // have positive mole fractions
386 for (size_t kcomp = 0; kcomp < m_numComponents; kcomp++) {
387 if (m_stoichCoeffRxnMatrix(kcomp,jrxn) < 0.0 && m_molNumSpecies_old[kcomp] <= VCS_DELETE_ELEMENTABS_CUTOFF*0.5) {
388 foundJrxn = false;
389 }
390 }
391 if (foundJrxn) {
392 return true;
393 }
394 } else if (m_stoichCoeffRxnMatrix(kspec,jrxn) < 0.0) {
395 // Second we are here if the component is a reactant in the
396 // reaction, and the reaction goes backwards.
397 foundJrxn = true;
398 size_t jspec = jrxn + m_numComponents;
400 foundJrxn = false;
401 continue;
402 }
403 // We can do the backwards reaction if all of the product
404 // components species are positive
405 for (size_t kcomp = 0; kcomp < m_numComponents; kcomp++) {
406 if (m_stoichCoeffRxnMatrix(kcomp,jrxn) > 0.0 && m_molNumSpecies_old[kcomp] <= VCS_DELETE_ELEMENTABS_CUTOFF*0.5) {
407 foundJrxn = false;
408 }
409 }
410 if (foundJrxn) {
411 return true;
412 }
413 }
414 }
415 }
416 }
417 return false;
418}
419
421{
422 size_t iphasePop = npos;
423 double FephaseMax = -1.0E30;
424 double Fephase = -1.0E30;
425
426 string note;
427 if (m_debug_print_lvl >= 2) {
428 plogf(" --- vcs_popPhaseID() called\n");
429 plogf(" --- Phase Status F_e MoleNum\n");
430 plogf(" --------------------------------------------------------------------------\n");
431 }
432 for (size_t iph = 0; iph < m_numPhases; iph++) {
433 vcs_VolPhase* Vphase = m_VolPhaseList[iph].get();
434 int existence = Vphase->exists();
435 note = "";
436 if (existence > 0) {
437 if (m_debug_print_lvl >= 2) {
438 plogf(" --- %18s %5d NA %11.3e\n",
439 Vphase->PhaseName, existence, m_tPhaseMoles_old[iph]);
440 }
441 } else {
442 if (Vphase->m_singleSpecies) {
443 // Single Phase Stability Resolution
444 size_t kspec = Vphase->spGlobalIndexVCS(0);
445 // Skip if the species is a component: a zero-abundance element may
446 // cause its only containing species to be selected as a component with
447 // zero moles, making its phase dead while still in the component basis.
448 // Such phases can never be popped.
449 if (kspec < m_numComponents) {
450 continue;
451 }
452 size_t irxn = kspec - m_numComponents;
453 if (irxn > m_deltaGRxn_old.size()) {
454 throw CanteraError("VCS_SOLVE::vcs_popPhaseID",
455 "Index out of bounds due to logic error.");
456 }
457 double deltaGRxn = m_deltaGRxn_old[irxn];
458 Fephase = exp(-deltaGRxn) - 1.0;
459 if (Fephase > 0.0) {
460 note = "(ready to be birthed)";
461 if (Fephase > FephaseMax) {
462 iphasePop = iph;
463 FephaseMax = Fephase;
464 note = "(chosen to be birthed)";
465 }
466 }
467 if (Fephase < 0.0) {
468 note = "(not stable)";
470 "VCS_SOLVE::vcs_popPhaseID", "shouldn't be here");
471 }
472
473 if (m_debug_print_lvl >= 2) {
474 plogf(" --- %18s %5d %10.3g %10.3g %s\n",
475 Vphase->PhaseName, existence, Fephase,
476 m_tPhaseMoles_old[iph], note.c_str());
477 }
478 } else if (vcs_popPhasePossible(iph)) {
479 // MultiSpecies Phase Stability Resolution
480 Fephase = vcs_phaseStabilityTest(iph);
481 if (Fephase > 0.0) {
482 if (Fephase > FephaseMax) {
483 iphasePop = iph;
484 FephaseMax = Fephase;
485 }
486 } else {
487 FephaseMax = std::max(FephaseMax, Fephase);
488 }
489 if (m_debug_print_lvl >= 2) {
490 plogf(" --- %18s %5d %11.3g %11.3g\n",
491 Vphase->PhaseName, existence, Fephase,
492 m_tPhaseMoles_old[iph]);
493 }
494 } else if (m_debug_print_lvl >= 2) {
495 plogf(" --- %18s %5d blocked %11.3g\n",
496 Vphase->PhaseName,
497 existence, m_tPhaseMoles_old[iph]);
498 }
499 }
500 }
501 // Insert logic here to figure out if phase pops are linked together. Only
502 // do one linked pop at a time.
503 if (m_debug_print_lvl >= 2) {
504 plogf(" ---------------------------------------------------------------------\n");
505 }
506 return iphasePop;
507}
508
509int VCS_SOLVE::vcs_popPhaseRxnStepSizes(const size_t iphasePop)
510{
511 vcs_VolPhase* Vphase = m_VolPhaseList[iphasePop].get();
512 // Identify the first species in the phase
513 size_t kspec = Vphase->spGlobalIndexVCS(0);
514 AssertThrowMsg(kspec >= m_numComponents, "VCS_SOLVE::vcs_popPhaseRxnStepSizes",
515 "called for a phase whose species is a component");
516 // Identify the formation reaction for that species
517 size_t irxn = kspec - m_numComponents;
518 vector<size_t> creationGlobalRxnNumbers(Vphase->nSpecies(), npos);
519
520 // Calculate the initial moles of the phase being born.
521 // Here we set it to 10x of the value which would cause the phase to be
522 // zeroed out within the algorithm. We may later adjust the value.
523 double tPhaseMoles = 10. * m_totalMolNum * VCS_DELETE_PHASE_CUTOFF;
524
525 AssertThrowMsg(!Vphase->exists(), "VCS_SOLVE::vcs_popPhaseRxnStepSizes",
526 "called for a phase that exists!");
527 if (m_debug_print_lvl >= 2) {
528 plogf(" --- vcs_popPhaseRxnStepSizes() called to pop phase %s %d into existence\n",
529 Vphase->PhaseName, iphasePop);
530 }
531 // Section for a single-species phase
532 if (Vphase->m_singleSpecies) {
533 double s = 0.0;
534 for (size_t j = 0; j < m_numComponents; ++j) {
535 if (!m_SSPhase[j] && m_molNumSpecies_old[j] > 0.0) {
536 s += pow(m_stoichCoeffRxnMatrix(j,irxn), 2) / m_molNumSpecies_old[j];
537 }
538 }
539 for (size_t j = 0; j < m_numPhases; j++) {
540 Vphase = m_VolPhaseList[j].get();
541 if (! Vphase->m_singleSpecies && m_tPhaseMoles_old[j] > 0.0) {
542 s -= pow(m_deltaMolNumPhase(j,irxn), 2) / m_tPhaseMoles_old[j];
543 }
544 }
545 if (s != 0.0) {
546 s = vcs_Hessian_diag_adj(irxn, s);
547 m_deltaMolNumSpecies[kspec] = -m_deltaGRxn_new[irxn] / s;
548 } else {
549 // Ok, s is equal to zero. We can not apply a sophisticated theory
550 // to birth the phase. Just pick a small delta and go with it.
551 m_deltaMolNumSpecies[kspec] = tPhaseMoles;
552 }
553
554 // section to do damping of the m_deltaMolNumSpecies[]
555 for (size_t j = 0; j < m_numComponents; ++j) {
556 double stoicC = m_stoichCoeffRxnMatrix(j,irxn);
557 if (stoicC != 0.0 && m_elType[j] == VCS_ELEM_TYPE_ABSPOS) {
558 double negChangeComp = - stoicC * m_deltaMolNumSpecies[kspec];
559 if (negChangeComp > m_molNumSpecies_old[j]) {
560 if (m_molNumSpecies_old[j] > 0.0) {
561 m_deltaMolNumSpecies[kspec] = - 0.5 * m_molNumSpecies_old[j] / stoicC;
562 } else {
563 m_deltaMolNumSpecies[kspec] = 0.0;
564 }
565 }
566 }
567 }
568 // Implement a damping term that limits m_deltaMolNumSpecies to the size
569 // of the mole number
570 if (-m_deltaMolNumSpecies[kspec] > m_molNumSpecies_old[kspec]) {
572 }
573 } else {
574 vector<double> fracDelta(Vphase->nSpecies());
575 vector<double> X_est(Vphase->nSpecies());
576 auto storedDelta = Vphase->creationMoleNumbers(creationGlobalRxnNumbers);
577 copy(storedDelta.begin(), storedDelta.end(), fracDelta.begin());
578
579 double sumFrac = 0.0;
580 for (size_t k = 0; k < Vphase->nSpecies(); k++) {
581 sumFrac += fracDelta[k];
582 }
583 for (size_t k = 0; k < Vphase->nSpecies(); k++) {
584 X_est[k] = fracDelta[k] / sumFrac;
585 }
586
587 double deltaMolNumPhase = tPhaseMoles;
588 double damp = 1.0;
590 double* molNumSpecies_tmp = m_deltaGRxn_tmp.data();
591
592 for (size_t k = 0; k < Vphase->nSpecies(); k++) {
593 kspec = Vphase->spGlobalIndexVCS(k);
594 double delmol = deltaMolNumPhase * X_est[k];
595 if (kspec >= m_numComponents) {
596 irxn = kspec - m_numComponents;
597 if (irxn > m_stoichCoeffRxnMatrix.nColumns()) {
598 throw CanteraError("VCS_SOLVE::vcs_popPhaseRxnStepSizes",
599 "Index out of bounds due to logic error.");
600 }
601 for (size_t j = 0; j < m_numComponents; ++j) {
602 double stoicC = m_stoichCoeffRxnMatrix(j,irxn);
603 if (stoicC != 0.0 && m_elType[j] == VCS_ELEM_TYPE_ABSPOS) {
604 molNumSpecies_tmp[j] += stoicC * delmol;
605 }
606 }
607 }
608 }
609
610 double ratioComp = 0.0;
611 for (size_t j = 0; j < m_numComponents; ++j) {
612 double deltaJ = m_molNumSpecies_old[j] - molNumSpecies_tmp[j];
613 if (molNumSpecies_tmp[j] < 0.0) {
614 ratioComp = 1.0;
615 if (deltaJ > 0.0) {
616 double delta0 = m_molNumSpecies_old[j];
617 damp = std::min(damp, delta0 / deltaJ * 0.9);
618 }
619 } else {
620 if (m_elType[j] == VCS_ELEM_TYPE_ABSPOS) {
621 size_t jph = m_phaseID[j];
622 if ((jph != iphasePop) && (!m_SSPhase[j])) {
623 double fdeltaJ = fabs(deltaJ);
624 if (m_molNumSpecies_old[j] > 0.0) {
625 ratioComp = std::max(ratioComp, fdeltaJ/ m_molNumSpecies_old[j]);
626 }
627 }
628 }
629 }
630 }
631
632 // We may have greatly underestimated the deltaMoles for the phase pop
633 // Here we create a damp > 1 to account for this possibility. We adjust
634 // upwards to make sure that a component in an existing multispecies
635 // phase is modified by a factor of 1/1000.
636 if (ratioComp > 1.0E-30 && ratioComp < 0.001) {
637 damp = 0.001 / ratioComp;
638 }
639 if (damp <= 1.0E-6) {
640 return 3;
641 }
642
643 for (size_t k = 0; k < Vphase->nSpecies(); k++) {
644 kspec = Vphase->spGlobalIndexVCS(k);
645 if (kspec < m_numComponents) {
647 } else {
648 m_deltaMolNumSpecies[kspec] = deltaMolNumPhase * X_est[k] * damp;
649 m_speciesStatus[kspec] = (X_est[k] > 1.0E-3) ? VCS_SPECIES_MAJOR
651 }
652 }
653 }
654 return 0;
655}
656
657size_t VCS_SOLVE::vcs_RxnStepSizes(int& forceComponentCalc, size_t& kSpecial)
658{
659 size_t iphDel = npos;
660 size_t k = 0;
661 if (m_debug_print_lvl >= 2) {
662 plogf(" ");
663 for (int j = 0; j < 82; j++) {
664 plogf("-");
665 }
666 plogf("\n");
667 plogf(" --- Subroutine vcs_RxnStepSizes called - Details:\n");
668 plogf(" ");
669 for (int j = 0; j < 82; j++) {
670 plogf("-");
671 }
672 plogf("\n");
673 plogf(" --- Species KMoles Rxn_Adjustment DeltaG\n");
674 }
675
676 // We update the matrix dlnActCoeffdmolNumber[][] at the top of the loop,
677 // when necessary
678 if (m_useActCoeffJac) {
680 }
681
682 // LOOP OVER THE FORMATION REACTIONS
683 for (size_t irxn = 0; irxn < m_numRxnRdc; ++irxn) {
684
685 size_t kspec = m_indexRxnToSpecies[irxn];
687 m_deltaMolNumSpecies[kspec] = 0.0;
689 if (m_molNumSpecies_old[kspec] == 0.0 && (!m_SSPhase[kspec])) {
690 // MULTISPECIES PHASE WITH total moles equal to zero
691 //
692 // If dg[irxn] is negative, then the multispecies phase should
693 // come alive again. Add a small positive step size to make it
694 // come alive.
695 if (m_deltaGRxn_new[irxn] < -1.0e-4) {
696 // First decide if this species is part of a multiphase that
697 // is nontrivial in size.
698 size_t iph = m_phaseID[kspec];
699 double tphmoles = m_tPhaseMoles_old[iph];
700 double trphmoles = tphmoles / m_totalMolNum;
701 vcs_VolPhase* Vphase = m_VolPhaseList[iph].get();
702 if (Vphase->exists() && (trphmoles > VCS_DELETE_PHASE_CUTOFF)) {
705 m_deltaMolNumSpecies[kspec] = 0.0;
706 } else {
708 }
709 } else {
710 m_deltaMolNumSpecies[kspec] = 0.0;
711 if (Vphase->exists() > 0 && trphmoles > 0.0) {
713 }
714 }
715 } else {
716 m_deltaMolNumSpecies[kspec] = 0.0;
717 }
718 } else {
719 // REGULAR PROCESSING
720 //
721 // First take care of cases where we want to bail out. Don't
722 // bother if superconvergence has already been achieved in this
723 // mode.
724 if (fabs(m_deltaGRxn_new[irxn]) <= m_tolmaj2) {
725 if (m_debug_print_lvl >= 2) {
726 plogf(" --- %-12.12s", m_speciesName[kspec]);
727 plogf(" %12.4E %12.4E %12.4E\n",
729 m_deltaGRxn_new[irxn]);
730 }
731 continue;
732 }
733
734 // Don't calculate for minor or nonexistent species if their
735 // values are to be decreasing anyway.
736 if ((m_speciesStatus[kspec] != VCS_SPECIES_MAJOR) && (m_deltaGRxn_new[irxn] >= 0.0)) {
737 if (m_debug_print_lvl >= 2) {
738 plogf(" --- %-12.12s", m_speciesName[kspec]);
739 plogf(" %12.4E %12.4E %12.4E\n",
741 m_deltaGRxn_new[irxn]);
742 }
743 continue;
744 }
745
746 // Start of the regular processing
747 double s = (m_SSPhase[kspec]) ? 0.0 : 1.0 / m_molNumSpecies_old[kspec];
748 for (size_t j = 0; j < m_numComponents; ++j) {
749 if (!m_SSPhase[j] && m_molNumSpecies_old[j] > 0.0) {
750 s += pow(m_stoichCoeffRxnMatrix(j,irxn), 2) / m_molNumSpecies_old[j];
751 }
752 }
753 for (size_t j = 0; j < m_numPhases; j++) {
754 vcs_VolPhase* Vphase = m_VolPhaseList[j].get();
755 if (!Vphase->m_singleSpecies && m_tPhaseMoles_old[j] > 0.0) {
756 s -= pow(m_deltaMolNumPhase(j,irxn), 2) / m_tPhaseMoles_old[j];
757 }
758 }
759 if (s != 0.0) {
760 // Take into account of the derivatives of the activity
761 // coefficients with respect to the mole numbers, even in
762 // our diagonal approximation.
763 if (m_useActCoeffJac) {
764 double s_old = s;
765 s = vcs_Hessian_diag_adj(irxn, s_old);
766 }
767
768 m_deltaMolNumSpecies[kspec] = -m_deltaGRxn_new[irxn] / s;
769 // New section to do damping of the m_deltaMolNumSpecies[]
770 for (size_t j = 0; j < m_numComponents; ++j) {
771 double stoicC = m_stoichCoeffRxnMatrix(j,irxn);
772 if (stoicC != 0.0) {
773 double negChangeComp = -stoicC * m_deltaMolNumSpecies[kspec];
774 if (negChangeComp > m_molNumSpecies_old[j]) {
775 if (m_molNumSpecies_old[j] > 0.0) {
776 m_deltaMolNumSpecies[kspec] = -m_molNumSpecies_old[j] / stoicC;
777 } else {
778 m_deltaMolNumSpecies[kspec] = 0.0;
779 }
780 }
781 }
782 }
783 // Implement a damping term that limits m_deltaMolNumSpecies
784 // to the size of the mole number
785 if (-m_deltaMolNumSpecies[kspec] > m_molNumSpecies_old[kspec]) {
787 }
788 } else {
789 // REACTION IS ENTIRELY AMONGST SINGLE SPECIES PHASES.
790 // DELETE ONE OF THE PHASES AND RECOMPUTE BASIS.
791 //
792 // Either the species L will disappear or one of the
793 // component single species phases will disappear. The sign
794 // of DG(I) will indicate which way the reaction will go.
795 // Then, we need to follow the reaction to see which species
796 // will zero out first. The species to be zeroed out will be
797 // "k".
798 double dss;
799 if (m_deltaGRxn_new[irxn] > 0.0) {
800 dss = m_molNumSpecies_old[kspec];
801 k = kspec;
802 for (size_t j = 0; j < m_numComponents; ++j) {
803 if (m_stoichCoeffRxnMatrix(j,irxn) > 0.0) {
804 double xx = m_molNumSpecies_old[j] / m_stoichCoeffRxnMatrix(j,irxn);
805 if (xx < dss) {
806 dss = xx;
807 k = j;
808 }
809 }
810 }
811 dss = -dss;
812 } else {
813 dss = 1.0e10;
814 for (size_t j = 0; j < m_numComponents; ++j) {
815 if (m_stoichCoeffRxnMatrix(j,irxn) < 0.0) {
816 double xx = -m_molNumSpecies_old[j] / m_stoichCoeffRxnMatrix(j,irxn);
817 if (xx < dss) {
818 dss = xx;
819 k = j;
820 }
821 }
822 }
823 }
824
825 // Here we adjust the mole fractions according to DSS and
826 // the stoichiometric array to take into account that we are
827 // eliminating the kth species. DSS contains the amount of
828 // moles of the kth species that needs to be added back into
829 // the component species.
830 if (dss != 0.0) {
831 if ((k == kspec) && (m_SSPhase[kspec] != 1)) {
832 // Found out that we can be in this spot, when
833 // components of multispecies phases are zeroed,
834 // leaving noncomponent species of the same phase
835 // having all of the mole numbers of that phases. it
836 // seems that we can suggest a zero of the species
837 // and the code will recover.
839 if (m_debug_print_lvl >= 2) {
840 plogf(" --- %-12.12s", m_speciesName[kspec]);
841 plogf(" %12.4E %12.4E %12.4E\n",
843 m_deltaGRxn_new[irxn]);
844 }
845 continue;
846 }
847
848 // Delete the single species phase
849 for (size_t j = 0; j < m_nsp; j++) {
850 m_deltaMolNumSpecies[j] = 0.0;
851 }
852 m_deltaMolNumSpecies[kspec] = dss;
853 for (size_t j = 0; j < m_numComponents; ++j) {
855 }
856
857 iphDel = m_phaseID[k];
858 kSpecial = k;
859
860 if (k != kspec) {
861 } else {
862 }
863 if (m_debug_print_lvl >= 2) {
864 plogf(" --- %-12.12s", m_speciesName[kspec]);
865 plogf(" %12.4E %12.4E %12.4E\n",
867 m_deltaGRxn_new[irxn]);
868 plogf(" --- vcs_RxnStepSizes Special section to set up to delete %s\n",
869 m_speciesName[k]);
870 }
871 if (k != kspec) {
872 forceComponentCalc = 1;
873 debuglog(" --- Force a component recalculation\n\n", m_debug_print_lvl >= 2);
874 }
875 if (m_debug_print_lvl >= 2) {
876 plogf(" ");
877 writeline('-', 82);
878 }
879 return iphDel;
880 }
881 }
882 } // End of regular processing
883 if (m_debug_print_lvl >= 2) {
884 plogf(" --- %-12.12s", m_speciesName[kspec]);
885 plogf(" %12.4E %12.4E %12.4E\n",
887 m_deltaGRxn_new[irxn]);
888 }
889 } // End of loop over m_speciesUnknownType
890 } // End of loop over non-component stoichiometric formation reactions
891 if (m_debug_print_lvl >= 2) {
892 plogf(" ");
893 writeline('-', 82);
894 }
895 return iphDel;
896}
897
898double VCS_SOLVE::vcs_phaseStabilityTest(const size_t iph)
899{
900 // We will use the _new state calc here
901 vcs_VolPhase* Vphase = m_VolPhaseList[iph].get();
902 const size_t nsp = Vphase->nSpecies();
903 int minNumberIterations = 3;
904 if (nsp <= 1) {
905 minNumberIterations = 1;
906 }
907
908 // We will do a full Newton calculation later, but for now, ...
909 bool doSuccessiveSubstitution = true;
910 double funcPhaseStability;
911 vector<double> X_est(nsp, 0.0);
912 vector<double> delFrac(nsp, 0.0);
913 vector<double> E_phi(nsp, 0.0);
914 vector<double> fracDelta_old(nsp, 0.0);
915 vector<double> fracDelta_raw(nsp, 0.0);
916 vector<size_t> creationGlobalRxnNumbers(nsp, npos);
918 vector<double> feSpecies_Deficient = m_feSpecies_old;
919
920 // get the activity coefficients
922
923 // Get the stored estimate for the composition of the phase if
924 // it gets created
925 vector<double> fracDelta_new(nsp, 0.0);
926 const auto& storedDelta = Vphase->creationMoleNumbers(creationGlobalRxnNumbers);
927 copy(storedDelta.begin(), storedDelta.end(), fracDelta_new.begin());
928
929 vector<size_t> componentList;
930 for (size_t k = 0; k < nsp; k++) {
931 size_t kspec = Vphase->spGlobalIndexVCS(k);
932 if (kspec < m_numComponents) {
933 componentList.push_back(k);
934 }
935 }
936
937 double normUpdate = 0.1 * vcs_l2norm(fracDelta_new);
938 double damp = 1.0E-2;
939
940 if (doSuccessiveSubstitution) {
941 int KP = 0;
942 if (m_debug_print_lvl >= 2) {
943 plogf(" --- vcs_phaseStabilityTest() called\n");
944 plogf(" --- Its X_old[%2d] FracDel_old[%2d] deltaF[%2d] FracDel_new[%2d]"
945 " normUpdate damp FuncPhaseStability\n", KP, KP, KP, KP);
946 plogf(" --------------------------------------------------------------"
947 "--------------------------------------------------------\n");
948 } else if (m_debug_print_lvl == 1) {
949 plogf(" --- vcs_phaseStabilityTest() called for phase %d\n", iph);
950 }
951
952 for (size_t k = 0; k < nsp; k++) {
953 if (fracDelta_new[k] < 1.0E-13) {
954 fracDelta_new[k] =1.0E-13;
955 }
956 }
957 bool converged = false;
958 double dirProd = 0.0;
959 for (int its = 0; its < 200 && (!converged); its++) {
960 double dampOld = damp;
961 double normUpdateOld = normUpdate;
962 fracDelta_old = fracDelta_new;
963 double dirProdOld = dirProd;
964
965 // Given a set of fracDelta's, we calculate the fracDelta's
966 // for the component species, if any
967 for (size_t i = 0; i < componentList.size(); i++) {
968 size_t kc = componentList[i];
969 size_t kc_spec = Vphase->spGlobalIndexVCS(kc);
970 fracDelta_old[kc] = 0.0;
971 for (size_t k = 0; k < nsp; k++) {
972 size_t kspec = Vphase->spGlobalIndexVCS(k);
973 if (kspec >= m_numComponents) {
974 size_t irxn = kspec - m_numComponents;
975 fracDelta_old[kc] += m_stoichCoeffRxnMatrix(kc_spec,irxn) * fracDelta_old[k];
976 }
977 }
978 }
979
980 // Now, calculate the predicted mole fractions, X_est[k]
981 double sumFrac = 0.0;
982 for (size_t k = 0; k < nsp; k++) {
983 sumFrac += fracDelta_old[k];
984 }
985 // Necessary because this can be identically zero. -> we need to fix
986 // this algorithm!
987 if (sumFrac <= 0.0) {
988 sumFrac = 1.0;
989 }
990 double sum_Xcomp = 0.0;
991 for (size_t k = 0; k < nsp; k++) {
992 X_est[k] = fracDelta_old[k] / sumFrac;
993 if (Vphase->spGlobalIndexVCS(k) < m_numComponents) {
994 sum_Xcomp += X_est[k];
995 }
996 }
997
998 // Feed the newly formed estimate of the mole fractions back into the
999 // ThermoPhase object
1001
1002 // get the activity coefficients
1004
1005 // First calculate altered chemical potentials for component species
1006 // belonging to this phase.
1007 for (size_t i = 0; i < componentList.size(); i++) {
1008 size_t kc = componentList[i];
1009 size_t kc_spec = Vphase->spGlobalIndexVCS(kc);
1010 double moleFrac = std::max(X_est[kc], VCS_DELETE_MINORSPECIES_CUTOFF);
1011 feSpecies_Deficient[kc_spec] = m_feSpecies_old[kc_spec]
1012 + log(m_actCoeffSpecies_new[kc_spec] * moleFrac);
1013 }
1014
1015 for (size_t i = 0; i < componentList.size(); i++) {
1016 size_t kc_spec = Vphase->spGlobalIndexVCS(componentList[i]);
1017 for (size_t k = 0; k < Vphase->nSpecies(); k++) {
1018 size_t kspec = Vphase->spGlobalIndexVCS(k);
1019 if (kspec >= m_numComponents) {
1020 size_t irxn = kspec - m_numComponents;
1021 if (i == 0) {
1023 }
1024 if (m_stoichCoeffRxnMatrix(kc_spec,irxn) != 0.0) {
1025 m_deltaGRxn_Deficient[irxn] +=
1026 m_stoichCoeffRxnMatrix(kc_spec,irxn) * (feSpecies_Deficient[kc_spec]- m_feSpecies_old[kc_spec]);
1027 }
1028 }
1029 }
1030 }
1031
1032 // Calculate the E_phi's
1033 double sum = 0.0;
1034 funcPhaseStability = sum_Xcomp - 1.0;
1035 for (size_t k = 0; k < nsp; k++) {
1036 size_t kspec = Vphase->spGlobalIndexVCS(k);
1037 if (kspec >= m_numComponents) {
1038 size_t irxn = kspec - m_numComponents;
1039 double deltaGRxn = clip(m_deltaGRxn_Deficient[irxn], -50.0, 50.0);
1040 E_phi[k] = std::exp(-deltaGRxn) / m_actCoeffSpecies_new[kspec];
1041 sum += E_phi[k];
1042 funcPhaseStability += E_phi[k];
1043 } else {
1044 E_phi[k] = 0.0;
1045 }
1046 }
1047
1048 // Calculate the raw estimate of the new fracs
1049 for (size_t k = 0; k < nsp; k++) {
1050 size_t kspec = Vphase->spGlobalIndexVCS(k);
1051 double b = E_phi[k] / sum * (1.0 - sum_Xcomp);
1052 if (kspec >= m_numComponents) {
1053 fracDelta_raw[k] = b;
1054 }
1055 }
1056
1057 // Given a set of fracDelta's, we calculate the fracDelta's
1058 // for the component species, if any
1059 for (size_t i = 0; i < componentList.size(); i++) {
1060 size_t kc = componentList[i];
1061 size_t kc_spec = Vphase->spGlobalIndexVCS(kc);
1062 fracDelta_raw[kc] = 0.0;
1063 for (size_t k = 0; k < nsp; k++) {
1064 size_t kspec = Vphase->spGlobalIndexVCS(k);
1065 if (kspec >= m_numComponents) {
1066 size_t irxn = kspec - m_numComponents;
1067 fracDelta_raw[kc] += m_stoichCoeffRxnMatrix(kc_spec,irxn) * fracDelta_raw[k];
1068 }
1069 }
1070 }
1071
1072 // Now possibly dampen the estimate.
1073 for (size_t k = 0; k < nsp; k++) {
1074 delFrac[k] = fracDelta_raw[k] - fracDelta_old[k];
1075 }
1076 normUpdate = vcs_l2norm(delFrac);
1077
1078 dirProd = 0.0;
1079 for (size_t k = 0; k < nsp; k++) {
1080 dirProd += fracDelta_old[k] * delFrac[k];
1081 }
1082 bool crossedSign = false;
1083 if (dirProd * dirProdOld < 0.0) {
1084 crossedSign = true;
1085 }
1086
1087 damp = 0.5;
1088 if (dampOld < 0.25) {
1089 damp = 2.0 * dampOld;
1090 }
1091 if (crossedSign) {
1092 if (normUpdate *1.5 > normUpdateOld) {
1093 damp = 0.5 * dampOld;
1094 } else if (normUpdate *2.0 > normUpdateOld) {
1095 damp = 0.8 * dampOld;
1096 }
1097 } else {
1098 if (normUpdate > normUpdateOld * 2.0) {
1099 damp = 0.6 * dampOld;
1100 } else if (normUpdate > normUpdateOld * 1.2) {
1101 damp = 0.9 * dampOld;
1102 }
1103 }
1104
1105 for (size_t k = 0; k < nsp; k++) {
1106 if (fabs(damp * delFrac[k]) > 0.3*fabs(fracDelta_old[k])) {
1107 damp = std::max(0.3*fabs(fracDelta_old[k]) / fabs(delFrac[k]),
1108 1.0E-8/fabs(delFrac[k]));
1109 }
1110 if (delFrac[k] < 0.0 && 2.0 * damp * (-delFrac[k]) > fracDelta_old[k]) {
1111 damp = fracDelta_old[k] / (2.0 * -delFrac[k]);
1112 }
1113 if (delFrac[k] > 0.0 && 2.0 * damp * delFrac[k] > fracDelta_old[k]) {
1114 damp = fracDelta_old[k] / (2.0 * delFrac[k]);
1115 }
1116 }
1117 damp = std::max(damp, 0.000001);
1118 for (size_t k = 0; k < nsp; k++) {
1119 fracDelta_new[k] = fracDelta_old[k] + damp * delFrac[k];
1120 }
1121
1122 if (m_debug_print_lvl >= 2) {
1123 plogf(" --- %3d %12g %12g %12g %12g %12g %12g %12g\n", its, X_est[KP], fracDelta_old[KP],
1124 delFrac[KP], fracDelta_new[KP], normUpdate, damp, funcPhaseStability);
1125 }
1126
1127 if (normUpdate < 1.0E-5 * damp) {
1128 converged = true;
1129 if (its < minNumberIterations) {
1130 converged = false;
1131 }
1132 }
1133 }
1134
1135 if (converged) {
1136 // Save the final optimized stated back into the VolPhase object for later use
1138
1139 // Save fracDelta for later use to initialize the problem better
1140 // @todo creationGlobalRxnNumbers needs to be calculated here and stored.
1141 Vphase->setCreationMoleNumbers(fracDelta_new, creationGlobalRxnNumbers);
1142 }
1143 } else {
1144 throw CanteraError("VCS_SOLVE::vcs_phaseStabilityTest", "not done yet");
1145 }
1146 if (m_debug_print_lvl >= 2) {
1147 plogf(" ------------------------------------------------------------"
1148 "-------------------------------------------------------------\n");
1149 } else if (m_debug_print_lvl == 1) {
1150 if (funcPhaseStability > 0.0) {
1151 plogf(" --- phase %d with func = %g is to be born\n", iph, funcPhaseStability);
1152 } else {
1153 plogf(" --- phase %d with func = %g stays dead\n", iph, funcPhaseStability);
1154 }
1155 }
1156 return funcPhaseStability;
1157}
1158
1159int VCS_SOLVE::vcs_evalSS_TP(int ipr, int ip1, double Temp, double pres)
1160{
1161 for (size_t iph = 0; iph < m_numPhases; iph++) {
1162 vcs_VolPhase* vph = m_VolPhaseList[iph].get();
1165 }
1166 for (size_t k = 0; k < m_nsp; k++) {
1168 }
1169
1170 return VCS_SUCCESS;
1171}
1172
1173void VCS_SOLVE::vcs_prep(int printLvl)
1174{
1175 m_debug_print_lvl = printLvl;
1176
1177 // Calculate the Single Species status of phases
1178 // Also calculate the number of species per phase
1179 vcs_SSPhase();
1180
1181 // Set an initial estimate for the number of noncomponent species equal to
1182 // nspecies - nelements. This may be changed below
1183 m_numRxnTot = (m_nelem > m_nsp) ? 0 : m_nsp - m_nelem;
1186 for (size_t i = 0; i < m_numRxnRdc; ++i) {
1188 }
1189
1190 for (size_t kspec = 0; kspec < m_nsp; ++kspec) {
1191 double sz = 0.0;
1192 for (size_t e = 0; e < m_nelem; e++) {
1193 sz += fabs(m_formulaMatrix(kspec, e));
1194 }
1195 m_spSize[kspec] = (sz > 0.0) ? sz : 1.0;
1196 }
1197
1198 // DETERMINE THE NUMBER OF COMPONENTS
1199 //
1200 // Obtain a valid estimate of the mole fraction. This will be used as an
1201 // initial ordering vector for prioritizing which species are defined as
1202 // components.
1203 //
1204 // If a mole number estimate was supplied from the input file, use that mole
1205 // number estimate.
1206 //
1207 // If a solution estimate wasn't supplied from the input file, supply an
1208 // initial estimate for the mole fractions based on the relative reverse
1209 // ordering of the chemical potentials.
1210 //
1211 // For voltage unknowns, set these to zero for the moment.
1212 double test = -1.0e-10;
1213 bool modifiedSoln = false;
1214 if (m_doEstimateEquil < 0) {
1215 double sum = 0.0;
1216 for (size_t kspec = 0; kspec < m_nsp; ++kspec) {
1218 sum += fabs(m_molNumSpecies_old[kspec]);
1219 }
1220 }
1221 if (fabs(sum) < 1.0E-6) {
1222 modifiedSoln = true;
1223 double pres = (m_pressurePA <= 0.0) ? 1.01325E5 : m_pressurePA;
1224 vcs_evalSS_TP(0, 0, m_temperature, pres);
1225 for (size_t kspec = 0; kspec < m_nsp; ++kspec) {
1227 m_molNumSpecies_old[kspec] = - m_SSfeSpecies[kspec];
1228 } else {
1229 m_molNumSpecies_old[kspec] = 0.0;
1230 }
1231 }
1232 }
1233 test = -1.0e20;
1234 }
1235
1236 // NC = number of components is in the vcs.h common block. This call to
1237 // BASOPT doesn't calculate the stoichiometric reaction matrix.
1238 vcs_basopt(true, test);
1239
1240 if (m_nsp >= m_numComponents) {
1242 for (size_t i = 0; i < m_numRxnRdc; ++i) {
1244 }
1245 } else {
1247 }
1248
1249 // The elements might need to be rearranged.
1251
1252 // If we mucked up the solution unknowns because they were all
1253 // zero to start with, set them back to zero here
1254 if (modifiedSoln) {
1255 for (size_t kspec = 0; kspec < m_nsp; ++kspec) {
1256 m_molNumSpecies_old[kspec] = 0.0;
1257 }
1258 }
1259
1260 // Initialize various arrays in the data to zero
1261 m_feSpecies_old.assign(m_feSpecies_old.size(), 0.0);
1262 m_feSpecies_new.assign(m_feSpecies_new.size(), 0.0);
1263 m_molNumSpecies_new.assign(m_molNumSpecies_new.size(), 0.0);
1266 m_deltaPhaseMoles.assign(m_deltaPhaseMoles.size(), 0.0);
1267 m_tPhaseMoles_new.assign(m_tPhaseMoles_new.size(), 0.0);
1268
1269 // Calculate the total number of moles in all phases.
1270 vcs_tmoles();
1271
1272 // Check to see if the current problem is well posed.
1273 double sum = 0.0;
1274 for (size_t e = 0; e < m_mix->nElements(); e++) {
1275 sum += m_mix->elementMoles(e);
1276 }
1277 if (sum < 1.0E-20) {
1278 throw CanteraError("VCS_SOLVE::vcs_prep", "The problem is not well posed"
1279 " because the total number of element moles is zero.");
1280 }
1281}
1282
1284{
1285 vector<double> awSpace(m_nsp + (m_nelem + 2)*(m_nelem), 0.0);
1286 vector<double> aw(m_nelem), sa(m_nelem), ss(m_nelem);
1287 vector<double> sm(m_nelem*m_nelem);
1288
1289 size_t ncomponents = m_numComponents;
1290 if (m_debug_print_lvl >= 2) {
1291 plogf(" ");
1292 writeline('-', 77);
1293 plogf(" --- Subroutine elem_rearrange() called to ");
1294 plogf("check stoich. coefficient matrix\n");
1295 plogf(" --- and to rearrange the element ordering once\n");
1296 }
1297
1298 // Use a temporary work array for the element numbers
1299 // Also make sure the value of test is unique.
1300 bool lindep = true;
1301 double test = -1.0E10;
1302 while (lindep) {
1303 lindep = false;
1304 for (size_t i = 0; i < m_nelem; ++i) {
1305 test -= 1.0;
1306 aw[i] = m_elemAbundancesGoal[i];
1307 if (test == aw[i]) {
1308 lindep = true;
1309 }
1310 }
1311 }
1312
1313 // Top of a loop of some sort based on the index JR. JR is the current
1314 // number independent elements found.
1315 size_t jr = 0;
1316 while (jr < ncomponents) {
1317 size_t k;
1318
1319 // Top of another loop point based on finding a linearly independent
1320 // species
1321 while (true) {
1322 // Search the remaining part of the mole fraction vector, AW, for
1323 // the largest remaining species. Return its identity in K.
1324 k = m_nelem;
1325 for (size_t ielem = jr; ielem < m_nelem; ielem++) {
1326 if (m_elementActive[ielem] && aw[ielem] != test) {
1327 k = ielem;
1328 break;
1329 }
1330 }
1331 if (k == m_nelem) {
1332 throw CanteraError("VCS_SOLVE::vcs_elem_rearrange",
1333 "Shouldn't be here. Algorithm misfired.");
1334 }
1335
1336 // Assign a large negative number to the element that we have just
1337 // found, in order to take it out of further consideration.
1338 aw[k] = test;
1339
1340 // CHECK LINEAR INDEPENDENCE OF CURRENT FORMULA MATRIX LINE WITH
1341 // PREVIOUS LINES OF THE FORMULA MATRIX
1342 //
1343 // Modified Gram-Schmidt Method, p. 202 Dalquist QR factorization of
1344 // a matrix without row pivoting.
1345 size_t jl = jr;
1346
1347 // Fill in the row for the current element, k, under consideration
1348 // The row will contain the Formula matrix value for that element
1349 // from the current component.
1350 for (size_t j = 0; j < ncomponents; ++j) {
1351 sm[j + jr*ncomponents] = m_formulaMatrix(j,k);
1352 }
1353 if (jl > 0) {
1354 // Compute the coefficients of JA column of the the upper
1355 // triangular R matrix, SS(J) = R_J_JR (this is slightly
1356 // different than Dalquist) R_JA_JA = 1
1357 for (size_t j = 0; j < jl; ++j) {
1358 ss[j] = 0.0;
1359 for (size_t i = 0; i < ncomponents; ++i) {
1360 ss[j] += sm[i + jr*ncomponents] * sm[i + j*ncomponents];
1361 }
1362 ss[j] /= sa[j];
1363 }
1364
1365 // Now make the new column, (*,JR), orthogonal to the previous
1366 // columns
1367 for (size_t j = 0; j < jl; ++j) {
1368 for (size_t i = 0; i < ncomponents; ++i) {
1369 sm[i + jr*ncomponents] -= ss[j] * sm[i + j*ncomponents];
1370 }
1371 }
1372 }
1373
1374 // Find the new length of the new column in Q. It will be used in
1375 // the denominator in future row calcs.
1376 sa[jr] = 0.0;
1377 for (size_t ml = 0; ml < ncomponents; ++ml) {
1378 sa[jr] += pow(sm[ml + jr*ncomponents], 2);
1379 }
1380 // IF NORM OF NEW ROW .LT. 1E-6 REJECT
1381 if (sa[jr] > 1.0e-6) {
1382 break;
1383 }
1384 }
1385 // REARRANGE THE DATA
1386 if (jr != k) {
1387 if (m_debug_print_lvl >= 2) {
1388 plogf(" --- %-2.2s(%9.2g) replaces %-2.2s(%9.2g) as element %3d\n",
1390 m_elementName[jr], m_elemAbundancesGoal[jr], jr);
1391 }
1392 vcs_switch_elem_pos(jr, k);
1393 std::swap(aw[jr], aw[k]);
1394 }
1395
1396 // If we haven't found enough components, go back and find some more.
1397 jr++;
1398 }
1399}
1400
1401void VCS_SOLVE::vcs_switch_elem_pos(size_t ipos, size_t jpos)
1402{
1403 if (ipos == jpos) {
1404 return;
1405 }
1406 AssertThrowMsg(ipos < m_nelem && jpos < m_nelem,
1407 "vcs_switch_elem_pos",
1408 "inappropriate args: {} {}", ipos, jpos);
1409
1410 // Change the element Global Index list in each vcs_VolPhase object
1411 // to reflect the switch in the element positions.
1412 for (size_t iph = 0; iph < m_numPhases; iph++) {
1413 vcs_VolPhase* volPhase = m_VolPhaseList[iph].get();
1414 for (size_t e = 0; e < volPhase->nElemConstraints(); e++) {
1415 if (volPhase->elemGlobalIndex(e) == ipos) {
1416 volPhase->setElemGlobalIndex(e, jpos);
1417 }
1418 if (volPhase->elemGlobalIndex(e) == jpos) {
1419 volPhase->setElemGlobalIndex(e, ipos);
1420 }
1421 }
1422 }
1423 std::swap(m_elemAbundancesGoal[ipos], m_elemAbundancesGoal[jpos]);
1424 std::swap(m_elemAbundances[ipos], m_elemAbundances[jpos]);
1425 std::swap(m_elementMapIndex[ipos], m_elementMapIndex[jpos]);
1426 std::swap(m_elType[ipos], m_elType[jpos]);
1427 std::swap(m_elementActive[ipos], m_elementActive[jpos]);
1428 for (size_t j = 0; j < m_nsp; ++j) {
1429 std::swap(m_formulaMatrix(j,ipos), m_formulaMatrix(j,jpos));
1430 }
1431 std::swap(m_elementName[ipos], m_elementName[jpos]);
1432}
1433
1434double VCS_SOLVE::vcs_Hessian_diag_adj(size_t irxn, double hessianDiag_Ideal)
1435{
1436 double diag = hessianDiag_Ideal;
1437 double hessActCoef = vcs_Hessian_actCoeff_diag(irxn);
1438 if (hessianDiag_Ideal <= 0.0) {
1439 throw CanteraError("VCS_SOLVE::vcs_Hessian_diag_adj",
1440 "We shouldn't be here");
1441 }
1442 if (hessActCoef >= 0.0) {
1443 diag += hessActCoef;
1444 } else if (fabs(hessActCoef) < 0.6666 * hessianDiag_Ideal) {
1445 diag += hessActCoef;
1446 } else {
1447 diag -= 0.6666 * hessianDiag_Ideal;
1448 }
1449 return diag;
1450}
1451
1453{
1454 size_t kspec = m_indexRxnToSpecies[irxn];
1455 size_t kph = m_phaseID[kspec];
1456 double np_kspec = std::max(m_tPhaseMoles_old[kph], 1e-13);
1457 auto sc_irxn = m_stoichCoeffRxnMatrix.col(irxn);
1458
1459 // First the diagonal term of the Jacobian
1460 double s = m_np_dLnActCoeffdMolNum(kspec,kspec) / np_kspec;
1461
1462 // Next, the other terms. Note this only a loop over the components So, it's
1463 // not too expensive to calculate.
1464 for (size_t j = 0; j < m_numComponents; j++) {
1465 if (!m_SSPhase[j]) {
1466 for (size_t k = 0; k < m_numComponents; ++k) {
1467 if (m_phaseID[k] == m_phaseID[j]) {
1468 double np = m_tPhaseMoles_old[m_phaseID[k]];
1469 if (np > 0.0) {
1470 s += sc_irxn[k] * sc_irxn[j] * m_np_dLnActCoeffdMolNum(j,k) / np;
1471 }
1472 }
1473 }
1474 if (kph == m_phaseID[j]) {
1475 s += sc_irxn[j] * (m_np_dLnActCoeffdMolNum(j,kspec) + m_np_dLnActCoeffdMolNum(kspec,j)) / np_kspec;
1476 }
1477 }
1478 }
1479 return s;
1480}
1481
1482void VCS_SOLVE::vcs_CalcLnActCoeffJac(span<const double> moleSpeciesVCS)
1483{
1484 for (auto& Vphase : m_VolPhaseList) {
1485 // We don't need to call single species phases;
1486 if (!Vphase->m_singleSpecies && !Vphase->isIdealSoln()) {
1487 // update the mole numbers
1488 Vphase->setMolesFromVCS(VCS_STATECALC_OLD, moleSpeciesVCS);
1489
1490 // Download the resulting calculation into the full vector. This
1491 // scatter calculation is carried out in the vcs_VolPhase object.
1492 Vphase->sendToVCS_LnActCoeffJac(m_np_dLnActCoeffdMolNum);
1493 }
1494 }
1495}
1496
1498{
1499 // SORT DEPENDENT SPECIES IN DECREASING ORDER OF MOLES
1500 vector<pair<double, size_t>> x_order;
1501 for (size_t i = 0; i < m_nsp; i++) {
1502 x_order.push_back({-m_molNumSpecies_old[i], i});
1503 }
1504 std::sort(x_order.begin() + m_numComponents,
1505 x_order.begin() + m_numSpeciesRdc);
1506
1507 vcs_setFlagsVolPhases(false, VCS_STATECALC_OLD);
1509
1510 // PRINT OUT RESULTS
1511 plogf("\n\n\n\n");
1512 writeline('-', 80);
1513 writeline('-', 80);
1514 plogf("\t\t VCS_TP REPORT\n");
1515 writeline('-', 80);
1516 writeline('-', 80);
1517 if (iconv < 0) {
1518 plogf(" ERROR: CONVERGENCE CRITERION NOT SATISFIED.\n");
1519 } else if (iconv == 1) {
1520 plogf(" RANGE SPACE ERROR: Equilibrium Found but not all Element Abundances are Satisfied\n");
1521 }
1522
1523 // Calculate some quantities that may need updating
1524 vcs_tmoles();
1525 double totalVolume = 0.0;
1526 for (size_t iphase = 0; iphase < m_numPhases; iphase++) {
1527 vcs_VolPhase* Vphase = m_VolPhaseList[iphase].get();
1530 double Volp = Vphase->sendToVCS_VolPM(m_PMVolumeSpecies);
1531 totalVolume += Volp;
1532 }
1533
1534 plogf("\t\tTemperature = %15.2g Kelvin\n", m_temperature);
1535 plogf("\t\tPressure = %15.5g Pa \n", m_pressurePA);
1536 plogf("\t\ttotal Volume = %15.5g m**3\n", totalVolume);
1537
1538 // TABLE OF SPECIES IN DECREASING MOLE NUMBERS
1539 plogf("\n\n");
1540 writeline('-', 80);
1541 plogf(" Species Equilibrium kmoles ");
1542 plogf("Mole Fraction ChemPot/RT SpecUnkType\n");
1543 writeline('-', 80);
1544 for (size_t i = 0; i < m_numComponents; ++i) {
1545 plogf(" %-12.12s", m_speciesName[i]);
1546 writeline(' ', 13, false);
1547 plogf("%14.7E %14.7E %12.4E", m_molNumSpecies_old[i],
1549 plogf(" %3d", m_speciesUnknownType[i]);
1550 plogf("\n");
1551 }
1552 for (size_t i = m_numComponents; i < m_numSpeciesRdc; ++i) {
1553 size_t j = x_order[i].second;
1554 plogf(" %-12.12s", m_speciesName[j]);
1555 writeline(' ', 13, false);
1556
1558 plogf("%14.7E %14.7E %12.4E", m_molNumSpecies_old[j],
1560 plogf(" KMolNum ");
1562 plogf(" NA %14.7E %12.4E", 1.0, m_feSpecies_old[j]);
1563 plogf(" Voltage = %14.7E", m_molNumSpecies_old[j]);
1564 } else {
1565 throw CanteraError("VCS_SOLVE::vcs_report", "we have a problem");
1566 }
1567 plogf("\n");
1568 }
1569 if (m_numSpeciesRdc != m_nsp) {
1570 plogf("\n SPECIES WITH LESS THAN 1.0E-32 KMOLES:\n\n");
1571 for (size_t kspec = m_numSpeciesRdc; kspec < m_nsp; ++kspec) {
1572 plogf(" %-12.12s", m_speciesName[kspec]);
1573 // Note m_deltaGRxn_new[] stores in kspec slot not irxn slot, after solve
1574 plogf(" %14.7E %14.7E %12.4E",
1575 m_molNumSpecies_old[kspec],
1576 m_molNumSpecies_new[kspec], m_deltaGRxn_new[kspec]);
1578 plogf(" KMol_Num");
1580 plogf(" Voltage");
1581 } else {
1582 plogf(" Unknown");
1583 }
1584 plogf("\n");
1585 }
1586 }
1587 writeline('-', 80);
1588 plogf("\n");
1589
1590 // TABLE OF SPECIES FORMATION REACTIONS
1591 writeline('-', m_numComponents*10 + 45, true, true);
1592 plogf(" |ComponentID|");
1593 for (size_t j = 0; j < m_numComponents; j++) {
1594 plogf(" %3d", j);
1595 }
1596 plogf(" | |\n");
1597 plogf(" | Components|");
1598 for (size_t j = 0; j < m_numComponents; j++) {
1599 plogf(" %10.10s", m_speciesName[j]);
1600 }
1601 plogf(" | |\n");
1602 plogf(" NonComponent | Moles |");
1603 for (size_t j = 0; j < m_numComponents; j++) {
1604 plogf(" %10.3g", m_molNumSpecies_old[j]);
1605 }
1606 plogf(" | DG/RT Rxn |\n");
1607 writeline('-', m_numComponents*10 + 45);
1608 for (size_t irxn = 0; irxn < m_numRxnTot; irxn++) {
1609 size_t kspec = m_indexRxnToSpecies[irxn];
1610 plogf(" %3d ", kspec);
1611 plogf("%-10.10s", m_speciesName[kspec]);
1612 plogf("|%10.3g |", m_molNumSpecies_old[kspec]);
1613 for (size_t j = 0; j < m_numComponents; j++) {
1614 plogf(" %6.2f", m_stoichCoeffRxnMatrix(j,irxn));
1615 }
1616 plogf(" |%10.3g |", m_deltaGRxn_new[irxn]);
1617 plogf("\n");
1618 }
1619 writeline('-', m_numComponents*10 + 45);
1620 plogf("\n");
1621
1622 // TABLE OF PHASE INFORMATION
1623 vector<double> gaTPhase(m_nelem, 0.0);
1624 double totalMoles = 0.0;
1625 double gibbsPhase = 0.0;
1626 double gibbsTotal = 0.0;
1627 plogf("\n\n");
1628 plogf("\n");
1629 writeline('-', m_nelem*10 + 58);
1630 plogf(" | ElementID |");
1631 for (size_t j = 0; j < m_nelem; j++) {
1632 plogf(" %3d", j);
1633 }
1634 plogf(" | |\n");
1635 plogf(" | Element |");
1636 for (size_t j = 0; j < m_nelem; j++) {
1637 plogf(" %10.10s", m_elementName[j]);
1638 }
1639 plogf(" | |\n");
1640 plogf(" PhaseName |KMolTarget |");
1641 for (size_t j = 0; j < m_nelem; j++) {
1642 plogf(" %10.3g", m_elemAbundancesGoal[j]);
1643 }
1644 plogf(" | Gibbs Total |\n");
1645 writeline('-', m_nelem*10 + 58);
1646 for (size_t iphase = 0; iphase < m_numPhases; iphase++) {
1647 plogf(" %3d ", iphase);
1648 vcs_VolPhase* VPhase = m_VolPhaseList[iphase].get();
1649 plogf("%-12.12s |",VPhase->PhaseName);
1650 plogf("%10.3e |", m_tPhaseMoles_old[iphase]);
1651 totalMoles += m_tPhaseMoles_old[iphase];
1652 if (m_tPhaseMoles_old[iphase] != VPhase->totalMoles() &&
1653 !vcs_doubleEqual(m_tPhaseMoles_old[iphase], VPhase->totalMoles())) {
1654 throw CanteraError("VCS_SOLVE::vcs_report", "we have a problem");
1655 }
1656 // Compute the elemental abundances for each element in this phase
1657 for (size_t j = 0; j < m_nelem; ++j) {
1658 double abundance_j = 0.0;
1659 for (size_t i = 0; i < m_nsp; ++i) {
1661 && m_phaseID[i] == iphase)
1662 {
1663 abundance_j += m_formulaMatrix(i,j) * m_molNumSpecies_old[i];
1664 }
1665 }
1666 plogf(" %10.3g", abundance_j);
1667 gaTPhase[j] += abundance_j;
1668 }
1669 gibbsPhase = 0.0;
1670 for (size_t kspec = 0; kspec < m_numSpeciesRdc; ++kspec) {
1671 if (m_phaseID[kspec] == iphase && m_speciesUnknownType[kspec] != VCS_SPECIES_TYPE_INTERFACIALVOLTAGE) {
1672 gibbsPhase += m_molNumSpecies_old[kspec] * m_feSpecies_old[kspec];
1673 }
1674 }
1675 gibbsTotal += gibbsPhase;
1676 plogf(" | %18.11E |\n", gibbsPhase);
1677 }
1678 writeline('-', m_nelem*10 + 58);
1679 plogf(" TOTAL |%10.3e |", totalMoles);
1680 for (size_t j = 0; j < m_nelem; j++) {
1681 plogf(" %10.3g", gaTPhase[j]);
1682 }
1683 plogf(" | %18.11E |\n", gibbsTotal);
1684
1685 writeline('-', m_nelem*10 + 58);
1686 plogf("\n");
1687
1688 // GLOBAL SATISFACTION INFORMATION
1689
1690 // Calculate the total dimensionless Gibbs Free Energy.
1692 plogf("\n\tTotal Dimensionless Gibbs Free Energy = G/RT = %15.7E\n", g);
1693 plogf("\nElemental Abundances (kmol): ");
1694 plogf(" Actual Target Type ElActive\n");
1695 for (size_t i = 0; i < m_nelem; ++i) {
1696 writeline(' ', 26, false);
1697 plogf("%-2.2s", m_elementName[i]);
1698 plogf("%20.12E %20.12E", m_elemAbundances[i], m_elemAbundancesGoal[i]);
1699 plogf(" %3d %3d\n", m_elType[i], m_elementActive[i]);
1700 }
1701 plogf("\n");
1702
1703 // TABLE OF SPECIES CHEM POTS
1704 writeline('-', 93, true, true);
1705 plogf("Chemical Potentials of the Species: (dimensionless)\n");
1706
1707 plogf("\t\t(RT = %g J/kmol)\n", GasConstant * m_temperature);
1708 plogf(" Name TKMoles StandStateChemPot "
1709 " ln(AC) ln(X_i) | F z_i phi | ChemPot | (-lnMnaught)");
1710 plogf("| (MolNum ChemPot)|");
1711 writeline('-', 147, true, true);
1712 for (size_t i = 0; i < m_nsp; ++i) {
1713 size_t j = x_order[i].second;
1714 size_t pid = m_phaseID[j];
1715 plogf(" %-12.12s", m_speciesName[j]);
1716 plogf(" %14.7E ", m_molNumSpecies_old[j]);
1717 plogf("%14.7E ", m_SSfeSpecies[j]);
1718 plogf("%14.7E ", log(m_actCoeffSpecies_old[j]));
1719 double tpmoles = m_tPhaseMoles_old[pid];
1720 double phi = m_phasePhi[pid];
1721 double eContrib = phi * m_chargeSpecies[j] * m_Faraday_dim;
1722 double lx = 0.0;
1724 lx = 0.0;
1725 } else {
1726 if (tpmoles > 0.0 && m_molNumSpecies_old[j] > 0.0) {
1727 double tmp = std::max(VCS_DELETE_MINORSPECIES_CUTOFF, m_molNumSpecies_old[j]);
1728 lx = log(tmp) - log(tpmoles);
1729 } else {
1730 lx = m_feSpecies_old[j] - m_SSfeSpecies[j]
1732 }
1733 }
1734 plogf("%14.7E |", lx);
1735 plogf("%14.7E | ", eContrib);
1736 double tmp = m_SSfeSpecies[j] + log(m_actCoeffSpecies_old[j])
1737 + lx - m_lnMnaughtSpecies[j] + eContrib;
1738 if (fabs(m_feSpecies_old[j] - tmp) > 1.0E-7) {
1739 throw CanteraError("VCS_SOLVE::vcs_report",
1740 "we have a problem - doesn't add up");
1741 }
1742 plogf(" %12.4E |", m_feSpecies_old[j]);
1743 if (m_lnMnaughtSpecies[j] != 0.0) {
1744 plogf("(%11.5E)", - m_lnMnaughtSpecies[j]);
1745 } else {
1746 plogf(" ");
1747 }
1748
1749 plogf("| %20.9E |", m_feSpecies_old[j] * m_molNumSpecies_old[j]);
1750 plogf("\n");
1751 }
1752 for (size_t i = 0; i < 125; i++) {
1753 plogf(" ");
1754 }
1755 plogf(" %20.9E\n", g);
1756 writeline('-', 147);
1757
1758 // TABLE OF SOLUTION COUNTERS
1759 plogf("\n");
1760 plogf("\nCounters: Iterations\n");
1761 plogf(" vcs_basopt: %5d\n", m_VCount->Basis_Opts);
1762 plogf(" vcs_TP: %5d\n", m_VCount->Its);
1763 writeline('-', 80);
1764 writeline('-', 80);
1765}
1766
1768{
1769 for (size_t j = 0; j < m_nelem; ++j) {
1770 m_elemAbundances[j] = 0.0;
1771 for (size_t i = 0; i < m_nsp; ++i) {
1774 }
1775 }
1776 }
1777}
1778
1780{
1781 size_t top = (ibound) ? m_nelem : m_numComponents;
1782 for (size_t i = 0; i < top; ++i) {
1783 // Require 12 digits of accuracy on non-zero constraints.
1784 if (m_elementActive[i] && fabs(m_elemAbundances[i] - m_elemAbundancesGoal[i]) > fabs(m_elemAbundancesGoal[i]) * 1.0e-12) {
1785 // This logic is for charge neutrality condition
1787 m_elemAbundancesGoal[i] != 0.0) {
1788 throw CanteraError("VCS_SOLVE::vcs_elabcheck",
1789 "Problem with charge neutrality condition");
1790 }
1793
1794 // Find out if the constraint is a multisign constraint. If it
1795 // is, then we have to worry about roundoff error in the
1796 // addition of terms. We are limited to 13 digits of finite
1797 // arithmetic accuracy.
1798 bool multisign = false;
1799 for (size_t kspec = 0; kspec < m_nsp; kspec++) {
1800 double eval = m_formulaMatrix(kspec,i);
1801 if (eval < 0.0) {
1802 multisign = true;
1803 }
1804 if (eval != 0.0) {
1805 scale = std::max(scale, fabs(eval * m_molNumSpecies_old[kspec]));
1806 }
1807 }
1808 if (multisign) {
1809 if (fabs(m_elemAbundances[i] - m_elemAbundancesGoal[i]) > 1e-11 * scale) {
1810 return false;
1811 }
1812 } else {
1814 return false;
1815 }
1816 }
1817 } else {
1818 return false;
1819 }
1820 }
1821 }
1822 return true;
1823}
1824
1825int VCS_SOLVE::vcs_elcorr(span<double> aa, span<double> x)
1826{
1827 checkArraySize("VCS_SOLVE::vcs_elcorr[x]", x.size(), m_nelem);
1828 checkArraySize("VCS_SOLVE::vcs_elcorr[aa]", aa.size(), m_nelem * m_nelem);
1829 int retn = 0;
1830
1831 vector<double> ga_save(m_elemAbundances);
1832 if (m_debug_print_lvl >= 2) {
1833 plogf(" --- vcsc_elcorr: Element abundances correction routine");
1834 if (m_nelem != m_numComponents) {
1835 plogf(" (m_numComponents != m_nelem)");
1836 }
1837 plogf("\n");
1838 }
1839
1840 for (size_t i = 0; i < m_nelem; ++i) {
1842 }
1843 double l2before = 0.0;
1844 for (size_t i = 0; i < m_nelem; ++i) {
1845 l2before += x[i] * x[i];
1846 }
1847 l2before = sqrt(l2before/m_nelem);
1848
1849 // Special section to take out single species, single component,
1850 // moles. These are species which have non-zero entries in the
1851 // formula matrix, and no other species have zero values either.
1852 bool changed = false;
1853 for (size_t i = 0; i < m_nelem; ++i) {
1854 int numNonZero = 0;
1855 bool multisign = false;
1856 for (size_t kspec = 0; kspec < m_nsp; kspec++) {
1858 double eval = m_formulaMatrix(kspec,i);
1859 if (eval < 0.0) {
1860 multisign = true;
1861 }
1862 if (eval != 0.0) {
1863 numNonZero++;
1864 }
1865 }
1866 }
1867 if (!multisign) {
1868 if (numNonZero < 2) {
1869 for (size_t kspec = 0; kspec < m_nsp; kspec++) {
1871 double eval = m_formulaMatrix(kspec,i);
1872 if (eval > 0.0) {
1873 m_molNumSpecies_old[kspec] = m_elemAbundancesGoal[i] / eval;
1874 changed = true;
1875 }
1876 }
1877 }
1878 } else {
1879 int numCompNonZero = 0;
1880 size_t compID = npos;
1881 for (size_t kspec = 0; kspec < m_numComponents; kspec++) {
1883 double eval = m_formulaMatrix(kspec,i);
1884 if (eval > 0.0) {
1885 compID = kspec;
1886 numCompNonZero++;
1887 }
1888 }
1889 }
1890 if (numCompNonZero == 1) {
1891 double diff = m_elemAbundancesGoal[i];
1892 for (size_t kspec = m_numComponents; kspec < m_nsp; kspec++) {
1894 double eval = m_formulaMatrix(kspec,i);
1895 diff -= eval * m_molNumSpecies_old[kspec];
1896 }
1897 m_molNumSpecies_old[compID] = std::max(0.0,diff/m_formulaMatrix(compID,i));
1898 changed = true;
1899 }
1900 }
1901 }
1902 }
1903 }
1904 if (changed) {
1905 vcs_elab();
1906 }
1907
1908 // Section to check for maximum bounds errors on all species due to
1909 // elements. This may only be tried on element types which are
1910 // VCS_ELEM_TYPE_ABSPOS. This is because no other species may have a
1911 // negative number of these.
1912 //
1913 // Note, also we can do this over ne, the number of elements, not just the
1914 // number of components.
1915 changed = false;
1916 for (size_t i = 0; i < m_nelem; ++i) {
1917 int elType = m_elType[i];
1918 if (elType == VCS_ELEM_TYPE_ABSPOS) {
1919 for (size_t kspec = 0; kspec < m_nsp; kspec++) {
1921 double atomComp = m_formulaMatrix(kspec,i);
1922 if (atomComp > 0.0) {
1923 double maxPermissible = m_elemAbundancesGoal[i] / atomComp;
1924 if (m_molNumSpecies_old[kspec] > maxPermissible) {
1925 if (m_debug_print_lvl >= 3) {
1926 plogf(" --- vcs_elcorr: Reduced species %s from %g to %g "
1927 "due to %s max bounds constraint\n",
1928 m_speciesName[kspec], m_molNumSpecies_old[kspec],
1929 maxPermissible, m_elementName[i]);
1930 }
1931 m_molNumSpecies_old[kspec] = maxPermissible;
1932 changed = true;
1934 m_molNumSpecies_old[kspec] = 0.0;
1935 if (m_SSPhase[kspec]) {
1937 } else {
1939 }
1940 if (m_debug_print_lvl >= 2) {
1941 plogf(" --- vcs_elcorr: Zeroed species %s and changed "
1942 "status to %d due to max bounds constraint\n",
1943 m_speciesName[kspec], m_speciesStatus[kspec]);
1944 }
1945 }
1946 }
1947 }
1948 }
1949 }
1950 }
1951 }
1952
1953 // Recalculate the element abundances if something has changed.
1954 if (changed) {
1955 vcs_elab();
1956 }
1957
1958 // Ok, do the general case. Linear algebra problem is of length nc, not ne,
1959 // as there may be degenerate rows when nc .ne. ne.
1961 for (size_t i = 0; i < m_numComponents; ++i) {
1963 if (fabs(x[i]) > 1.0E-13) {
1964 retn = 1;
1965 }
1966 for (size_t j = 0; j < m_numComponents; ++j) {
1967 A(j, i) = - m_formulaMatrix(i,j);
1968 }
1969 }
1970
1971 solve(A, x.first(m_nelem), 1, m_nelem);
1972
1973 // Now apply the new direction without creating negative species.
1974 double par = 0.5;
1975 for (size_t i = 0; i < m_numComponents; ++i) {
1976 if (m_molNumSpecies_old[i] > 0.0) {
1977 par = std::max(par, -x[i] / m_molNumSpecies_old[i]);
1978 }
1979 }
1980 par = std::min(par, 100.0);
1981 par = 1.0 / par;
1982 if (par < 1.0 && par > 0.0) {
1983 retn = 2;
1984 par *= 0.9999;
1985 } else {
1986 par = 1.0;
1987 }
1988 for (size_t i = 0; i < m_numComponents; ++i) {
1989 double tmp = m_molNumSpecies_old[i] + par * x[i];
1990 if (tmp > 0.0) {
1991 m_molNumSpecies_old[i] = tmp;
1992 } else if (m_SSPhase[i]) {
1993 m_molNumSpecies_old[i] = 0.0;
1994 } else {
1995 m_molNumSpecies_old[i] *= 0.0001;
1996 }
1997 }
1998
1999 // We have changed the element abundances. Calculate them again
2000 vcs_elab();
2001
2002 // We have changed the total moles in each phase. Calculate them again
2003 vcs_tmoles();
2004
2005 // cleanup/logging function to be called before returning
2006 auto finalize = [&]() {
2007 vcs_tmoles();
2008 double l2after = 0.0;
2009 for (size_t i = 0; i < m_nelem; ++i) {
2010 l2after += pow(m_elemAbundances[i] - m_elemAbundancesGoal[i], 2);
2011 }
2012 l2after = sqrt(l2after/m_nelem);
2013 if (m_debug_print_lvl >= 2) {
2014 plogf(" --- Elem_Abund: Correct Initial "
2015 " Final\n");
2016 for (size_t i = 0; i < m_nelem; ++i) {
2017 plogf(" --- ");
2018 plogf("%-2.2s", m_elementName[i]);
2019 plogf(" %20.12E %20.12E %20.12E\n", m_elemAbundancesGoal[i], ga_save[i], m_elemAbundances[i]);
2020 }
2021 plogf(" --- Diff_Norm: %20.12E %20.12E\n",
2022 l2before, l2after);
2023 }
2024 };
2025
2026 // Try some ad hoc procedures for fixing the problem
2027 if (retn >= 2) {
2028 // First find a species whose adjustment is a win-win situation.
2029 for (size_t kspec = 0; kspec < m_nsp; kspec++) {
2031 continue;
2032 }
2033 double saveDir = 0.0;
2034 bool goodSpec = true;
2035 for (size_t i = 0; i < m_numComponents; ++i) {
2036 double dir = m_formulaMatrix(kspec,i) * (m_elemAbundancesGoal[i] - m_elemAbundances[i]);
2037 if (fabs(dir) > 1.0E-10) {
2038 if (dir > 0.0) {
2039 if (saveDir < 0.0) {
2040 goodSpec = false;
2041 break;
2042 }
2043 } else if (saveDir > 0.0) {
2044 goodSpec = false;
2045 break;
2046 }
2047 saveDir = dir;
2048 } else if (m_formulaMatrix(kspec,i) != 0.0) {
2049 goodSpec = false;
2050 break;
2051 }
2052 }
2053 if (goodSpec) {
2054 int its = 0;
2055 double xx = 0.0;
2056 for (size_t i = 0; i < m_numComponents; ++i) {
2057 if (m_formulaMatrix(kspec,i) != 0.0) {
2058 xx += (m_elemAbundancesGoal[i] - m_elemAbundances[i]) / m_formulaMatrix(kspec,i);
2059 its++;
2060 }
2061 }
2062 if (its > 0) {
2063 xx /= its;
2064 }
2065 m_molNumSpecies_old[kspec] += xx;
2066 m_molNumSpecies_old[kspec] = std::max(m_molNumSpecies_old[kspec], 1.0E-10);
2067
2068 // If we are dealing with a deleted species, then we need to
2069 // reinsert it into the active list.
2070 if (kspec >= m_numSpeciesRdc) {
2071 vcs_reinsert_deleted(kspec);
2073 vcs_elab();
2074 finalize();
2075 return retn;
2076 }
2077 vcs_elab();
2078 }
2079 }
2080 }
2081 if (vcs_elabcheck(0)) {
2082 finalize();
2083 return 1;
2084 }
2085
2086 for (size_t i = 0; i < m_nelem; ++i) {
2088 (m_elType[i] == VCS_ELEM_TYPE_ABSPOS && m_elemAbundancesGoal[i] == 0.0)) {
2089 for (size_t kspec = 0; kspec < m_numSpeciesRdc; kspec++) {
2090 if (m_elemAbundances[i] > 0.0 && m_formulaMatrix(kspec,i) < 0.0) {
2092 m_molNumSpecies_old[kspec] = std::max(m_molNumSpecies_old[kspec], 0.0);
2093 vcs_elab();
2094 break;
2095 }
2096 if (m_elemAbundances[i] < 0.0 && m_formulaMatrix(kspec,i) > 0.0) {
2098 m_molNumSpecies_old[kspec] = std::max(m_molNumSpecies_old[kspec], 0.0);
2099 vcs_elab();
2100 break;
2101 }
2102 }
2103 }
2104 }
2105 if (vcs_elabcheck(1)) {
2106 finalize();
2107 return 1;
2108 }
2109
2110 // For electron charges element types, we try positive deltas in the species
2111 // concentrations to match the desired electron charge exactly.
2112 for (size_t i = 0; i < m_nelem; ++i) {
2113 double dev = m_elemAbundancesGoal[i] - m_elemAbundances[i];
2114 if (m_elType[i] == VCS_ELEM_TYPE_ELECTRONCHARGE && (fabs(dev) > 1.0E-300)) {
2115 bool useZeroed = true;
2116 for (size_t kspec = 0; kspec < m_numSpeciesRdc; kspec++) {
2117 if (dev < 0.0) {
2118 if (m_formulaMatrix(kspec,i) < 0.0 && m_molNumSpecies_old[kspec] > 0.0) {
2119 useZeroed = false;
2120 }
2121 } else {
2122 if (m_formulaMatrix(kspec,i) > 0.0 && m_molNumSpecies_old[kspec] > 0.0) {
2123 useZeroed = false;
2124 }
2125 }
2126 }
2127 for (size_t kspec = 0; kspec < m_numSpeciesRdc; kspec++) {
2128 if (m_molNumSpecies_old[kspec] > 0.0 || useZeroed) {
2129 if (dev < 0.0 && m_formulaMatrix(kspec,i) < 0.0) {
2130 double delta = dev / m_formulaMatrix(kspec,i);
2131 m_molNumSpecies_old[kspec] += delta;
2132 m_molNumSpecies_old[kspec] = std::max(m_molNumSpecies_old[kspec], 0.0);
2133 vcs_elab();
2134 break;
2135 }
2136 if (dev > 0.0 && m_formulaMatrix(kspec,i) > 0.0) {
2137 double delta = dev / m_formulaMatrix(kspec,i);
2138 m_molNumSpecies_old[kspec] += delta;
2139 m_molNumSpecies_old[kspec] = std::max(m_molNumSpecies_old[kspec], 0.0);
2140 vcs_elab();
2141 break;
2142 }
2143 }
2144 }
2145 }
2146 }
2147 if (vcs_elabcheck(1)) {
2148 finalize();
2149 return 1;
2150 }
2151 return retn;
2152}
2153
2155{
2156 double test = -1.0E-10;
2157
2158 if (m_debug_print_lvl >= 2) {
2159 plogf(" --- call setInitialMoles\n");
2160 }
2161
2162 int iter = 0;
2163 bool abundancesOK = true;
2164 vector<double> sm(m_nelem * m_nelem, 0.0);
2165 vector<double> ss(m_nelem, 0.0);
2166 vector<double> sa(m_nelem, 0.0);
2167 vector<double> wx(m_nelem, 0.0);
2168 vector<double> aw(m_nsp, 0.0);
2169
2170 for (size_t ik = 0; ik < m_nsp; ik++) {
2172 m_molNumSpecies_old[ik] = max(0.0, m_molNumSpecies_old[ik]);
2173 }
2174 }
2175
2176 if (m_debug_print_lvl >= 2) {
2178 }
2179
2180 bool redo = true;
2181 while (redo) {
2182 if (!vcs_elabcheck(0)) {
2183 if (m_debug_print_lvl >= 2) {
2184 plogf(" --- seMolesLinProg Mole numbers failing element abundances\n");
2185 plogf(" --- seMolesLinProg Call vcs_elcorr to attempt fix\n");
2186 }
2187 int retn = vcs_elcorr(sm, wx);
2188 abundancesOK = (retn < 2);
2189 } else {
2190 abundancesOK = true;
2191 }
2192
2193 // Now find the optimized basis that spans the stoichiometric
2194 // coefficient matrix, based on the current composition,
2195 // m_molNumSpecies_old[] We also calculate sc[][], the reaction matrix.
2196 vcs_basopt(false, test);
2197
2198 if (m_debug_print_lvl >= 2) {
2199 plogf("iteration %d\n", iter);
2200 }
2201 redo = false;
2202 iter++;
2203 if (iter > 15) {
2204 break;
2205 }
2206
2207 // loop over all reactions
2208 for (size_t irxn = 0; irxn < m_numRxnTot; irxn++) {
2209 // dg_rt is the Delta_G / RT value for the reaction
2210 size_t ik = m_numComponents + irxn;
2211 double dg_rt = m_SSfeSpecies[ik];
2212 double dxi_min = 1.0e10;
2213 auto sc_irxn = m_stoichCoeffRxnMatrix.col(irxn);
2214 for (size_t jcomp = 0; jcomp < m_nelem; jcomp++) {
2215 dg_rt += m_SSfeSpecies[jcomp] * sc_irxn[jcomp];
2216 }
2217 // fwd or rev direction.
2218 // idir > 0 implies increasing the current species
2219 // idir < 0 implies decreasing the current species
2220 int idir = (dg_rt < 0.0 ? 1 : -1);
2221 if (idir < 0) {
2222 dxi_min = m_molNumSpecies_old[ik];
2223 }
2224
2225 for (size_t jcomp = 0; jcomp < m_numComponents; jcomp++) {
2226 double nu = sc_irxn[jcomp];
2227 // set max change in progress variable by
2228 // non-negativity requirement
2229 if (nu*idir < 0) {
2230 double delta_xi = fabs(m_molNumSpecies_old[jcomp]/nu);
2231 // if a component has nearly zero moles, redo
2232 // with a new set of components
2233 if (!redo && delta_xi < 1.0e-10 && (m_molNumSpecies_old[ik] >= 1.0E-10)) {
2234 if (m_debug_print_lvl >= 2) {
2235 plogf(" --- Component too small: %s\n", m_speciesName[jcomp]);
2236 }
2237 redo = true;
2238 }
2239 dxi_min = std::min(dxi_min, delta_xi);
2240 }
2241 }
2242
2243 // step the composition by dxi_min, check against zero, since
2244 // we are zeroing components and species on every step.
2245 // Redo the iteration, if a component went from positive to zero on this step.
2246 double dsLocal = idir*dxi_min;
2247 m_molNumSpecies_old[ik] += dsLocal;
2248 m_molNumSpecies_old[ik] = max(0.0, m_molNumSpecies_old[ik]);
2249 for (size_t jcomp = 0; jcomp < m_numComponents; jcomp++) {
2250 bool full = false;
2251 if (m_molNumSpecies_old[jcomp] > 1.0E-15) {
2252 full = true;
2253 }
2254 m_molNumSpecies_old[jcomp] += sc_irxn[jcomp] * dsLocal;
2255 m_molNumSpecies_old[jcomp] = max(0.0, m_molNumSpecies_old[jcomp]);
2256 if (full && m_molNumSpecies_old[jcomp] < 1.0E-60) {
2257 redo = true;
2258 }
2259 }
2260 }
2261
2262 if (m_debug_print_lvl >= 2) {
2264 }
2265 }
2266
2267 if (m_debug_print_lvl == 1) {
2269 plogf(" --- setInitialMoles end\n");
2270 }
2271 if (!abundancesOK) {
2272 return -1;
2273 } else if (iter > 15) {
2274 return 1;
2275 }
2276 return 0;
2277}
2278
2279double VCS_SOLVE::vcs_Total_Gibbs(span<const double> molesSp,
2280 span<const double> chemPot,
2281 span<const double> tPhMoles)
2282{
2283 double g = 0.0;
2284
2285 for (size_t kspec = 0; kspec < m_numSpeciesRdc; ++kspec) {
2287 g += molesSp[kspec] * chemPot[kspec];
2288 }
2289 }
2290
2291 return g;
2292}
2293
2295{
2296 // Whether we have an estimate or not gets overwritten on
2297 // the call to the equilibrium solver.
2298 m_temperature = m_mix->temperature();
2299 m_pressurePA = m_mix->pressure();
2301
2302 for (size_t iphase = 0; iphase < m_numPhases; iphase++) {
2303 vcs_VolPhase* volPhase = m_VolPhaseList[iphase].get();
2304
2306
2307 // Loop through each species in the current phase
2308 size_t nSpPhase = m_mix->phase(iphase).nSpecies();
2309 if ((nSpPhase == 1) && (volPhase->phiVarIndex() == 0)) {
2311 } else if (volPhase->totalMoles() > 0.0) {
2313 } else {
2315 }
2316 }
2317
2318 // Printout the species information: PhaseID's and mole nums
2319 if (m_printLvl > 1) {
2320 writeline('=', 80, true, true);
2321 writeline('=', 20, false);
2322 plogf(" Cantera_to_vprob: START OF PROBLEM STATEMENT ");
2323 writeline('=', 20);
2324 writeline('=', 80);
2325 plogf("\n");
2326 plogf(" Phase IDs of species\n");
2327 plogf(" species phaseID phaseName ");
2328 plogf(" Initial_Estimated_kMols\n");
2329 for (size_t i = 0; i < m_nsp; i++) {
2330 size_t iphase = m_phaseID[i];
2331 vcs_VolPhase* VolPhase = m_VolPhaseList[iphase].get();
2332 plogf("%16s %5d %16s", m_speciesName[i].c_str(), iphase,
2333 VolPhase->PhaseName.c_str());
2335 plogf(" Volts = %-10.5g\n", m_molNumSpecies_old[i]);
2336 } else {
2337 plogf(" %-10.5g\n", m_molNumSpecies_old[i]);
2338 }
2339 }
2340
2341 // Printout of the Phase structure information
2342 writeline('-', 80, true, true);
2343 plogf(" Information about phases\n");
2344 plogf(" PhaseName PhaseNum SingSpec EqnState NumSpec "
2345 "Tmoles(kmol)\n");
2346
2347 for (size_t iphase = 0; iphase < m_numPhases; iphase++) {
2348 vcs_VolPhase* VolPhase = m_VolPhaseList[iphase].get();
2349 plogf("%16s %8d %8d %16s %8d ", VolPhase->PhaseName.c_str(),
2350 VolPhase->VP_ID_, VolPhase->m_singleSpecies,
2351 VolPhase->eos_name(), VolPhase->nSpecies());
2352 plogf("%16e\n", VolPhase->totalMoles());
2353 }
2354
2355 writeline('=', 80, true, true);
2356 writeline('=', 20, false);
2357 plogf(" Cantera_to_vprob: END OF PROBLEM STATEMENT ");
2358 writeline('=', 20);
2359 writeline('=', 80);
2360 plogf("\n");
2361 }
2362
2363 // m_numRxnTot = number of noncomponents, also equal to the number of
2364 // reactions. Note, it's possible that the number of elements is greater
2365 // than the number of species. In that case set the number of reactions to
2366 // zero.
2367 m_numRxnTot = (m_nelem > m_nsp) ? 0 : m_nsp - m_nelem;
2369}
2370
2372{
2373 const char* pprefix = " --- vcs_inest: ";
2374 if (m_doEstimateEquil > 0) {
2375 // Calculate the elemental abundances
2376 vcs_elab();
2377 if (vcs_elabcheck(0)) {
2378 if (m_debug_print_lvl >= 2) {
2379 plogf("%s Initial guess passed element abundances on input\n", pprefix);
2380 plogf("%s m_doEstimateEquil = 1 so will use the input mole "
2381 "numbers as estimates\n", pprefix);
2382 }
2383 return;
2384 } else if (m_debug_print_lvl >= 2) {
2385 plogf("%s Initial guess failed element abundances on input\n", pprefix);
2386 plogf("%s m_doEstimateEquil = 1 so will discard input "
2387 "mole numbers and find our own estimate\n", pprefix);
2388 }
2389 }
2390
2391 // temporary space for usage in this routine and in subroutines
2392 vector<double> sm(m_nelem*m_nelem, 0.0);
2393 vector<double> ss(m_nelem, 0.0);
2394 vector<double> sa(m_nelem, 0.0);
2395 vector<double> aw(m_nsp + m_nelem, 0.0);
2396
2397 // Go get the estimate of the solution
2398 if (m_debug_print_lvl >= 2) {
2399 plogf("%sGo find an initial estimate for the equilibrium problem\n",
2400 pprefix);
2401 }
2402 double test = -1.0E20;
2403 size_t nrxn = m_numRxnTot;
2404
2405 // CALL ROUTINE TO SOLVE MAX(CC*molNum) SUCH THAT AX*molNum = BB AND
2406 // molNum(I) .GE. 0.0. Note, both of these programs do this.
2408
2409 if (m_debug_print_lvl >= 2) {
2410 plogf("%s Mole Numbers returned from linear programming (vcs_inest initial guess):\n",
2411 pprefix);
2412 plogf("%s SPECIES MOLE_NUMBER -SS_ChemPotential\n", pprefix);
2413 for (size_t kspec = 0; kspec < m_nsp; ++kspec) {
2414 plogf("%s ", pprefix);
2415 plogf("%-12.12s", m_speciesName[kspec]);
2416 plogf(" %15.5g %12.3g\n", m_molNumSpecies_old[kspec], -m_SSfeSpecies[kspec]);
2417 }
2418 plogf("%s Element Abundance Agreement returned from linear "
2419 "programming (vcs_inest initial guess):\n", pprefix);
2420 plogf("%s Element Goal Actual\n", pprefix);
2421 for (size_t j = 0; j < m_nelem; j++) {
2422 if (m_elementActive[j]) {
2423 double tmp = 0.0;
2424 for (size_t kspec = 0; kspec < m_nsp; ++kspec) {
2425 tmp += m_formulaMatrix(kspec,j) * m_molNumSpecies_old[kspec];
2426 }
2427 plogf("%s ", pprefix);
2428 plogf(" %-9.9s", m_elementName[j]);
2429 plogf(" %12.3g %12.3g\n", m_elemAbundancesGoal[j], tmp);
2430 }
2431 }
2432 writelogendl();
2433 }
2434
2435 // Make sure all species have positive definite mole numbers Set voltages to
2436 // zero for now, until we figure out what to do
2437 m_deltaMolNumSpecies.assign(m_deltaMolNumSpecies.size(), 0.0);
2438 for (size_t kspec = 0; kspec < m_nsp; ++kspec) {
2440 if (m_molNumSpecies_old[kspec] <= 0.0) {
2441 // HKM Should eventually include logic here for non SS phases
2442 if (!m_SSPhase[kspec]) {
2443 m_molNumSpecies_old[kspec] = 1.0e-30;
2444 }
2445 }
2446 } else {
2447 m_molNumSpecies_old[kspec] = 0.0;
2448 }
2449 }
2450
2451 // Now find the optimized basis that spans the stoichiometric coefficient matrix
2452 vcs_basopt(false, test);
2453
2454 // CALCULATE TOTAL MOLES, CHEMICAL POTENTIALS OF BASIS
2455
2456 // Calculate TMoles and m_tPhaseMoles_old[]
2457 vcs_tmoles();
2458
2459 // m_tPhaseMoles_new[] will consist of just the component moles
2460 m_tPhaseMoles_new.assign(m_numPhases, 1.0e-20);
2461 for (size_t kspec = 0; kspec < m_numComponents; ++kspec) {
2464 }
2465 }
2466 double TMolesMultiphase = 0.0;
2467 for (size_t iph = 0; iph < m_numPhases; iph++) {
2468 if (! m_VolPhaseList[iph]->m_singleSpecies) {
2469 TMolesMultiphase += m_tPhaseMoles_new[iph];
2470 }
2471 }
2473 for (size_t kspec = 0; kspec < m_numComponents; ++kspec) {
2475 m_molNumSpecies_new[kspec] = 0.0;
2476 }
2477 }
2479
2480 for (size_t kspec = 0; kspec < m_numComponents; ++kspec) {
2482 if (! m_SSPhase[kspec]) {
2483 size_t iph = m_phaseID[kspec];
2484 m_feSpecies_new[kspec] += log(m_molNumSpecies_new[kspec] / m_tPhaseMoles_old[iph]);
2485 }
2486 } else {
2487 m_molNumSpecies_new[kspec] = 0.0;
2488 }
2489 }
2490 vcs_deltag(0, true, VCS_STATECALC_NEW);
2491 if (m_debug_print_lvl >= 2) {
2492 for (size_t kspec = 0; kspec < m_nsp; ++kspec) {
2493 plogf("%s", pprefix);
2494 plogf("%-12.12s", m_speciesName[kspec]);
2495 if (kspec < m_numComponents) {
2496 plogf("fe* = %15.5g ff = %15.5g\n", m_feSpecies_new[kspec],
2497 m_SSfeSpecies[kspec]);
2498 } else {
2499 plogf("fe* = %15.5g ff = %15.5g dg* = %15.5g\n",
2501 }
2502 }
2503 }
2504
2505 // ESTIMATE REACTION ADJUSTMENTS
2506 vector<double>& xtphMax = m_TmpPhase;
2507 vector<double>& xtphMin = m_TmpPhase2;
2508 m_deltaPhaseMoles.assign(m_deltaPhaseMoles.size(), 0.0);
2509 for (size_t iph = 0; iph < m_numPhases; iph++) {
2510 xtphMax[iph] = log(m_tPhaseMoles_new[iph] * 1.0E32);
2511 xtphMin[iph] = log(m_tPhaseMoles_new[iph] * 1.0E-32);
2512 }
2513 for (size_t irxn = 0; irxn < nrxn; ++irxn) {
2514 size_t kspec = m_indexRxnToSpecies[irxn];
2515
2516 // For single species phases, we will not estimate the mole numbers. If
2517 // the phase exists, it stays. If it doesn't exist in the estimate, it
2518 // doesn't come into existence here.
2519 if (! m_SSPhase[kspec]) {
2520 size_t iph = m_phaseID[kspec];
2521 if (m_deltaGRxn_new[irxn] > xtphMax[iph]) {
2522 m_deltaGRxn_new[irxn] = 0.8 * xtphMax[iph];
2523 }
2524 if (m_deltaGRxn_new[irxn] < xtphMin[iph]) {
2525 m_deltaGRxn_new[irxn] = 0.8 * xtphMin[iph];
2526 }
2527
2528 // HKM -> The TMolesMultiphase is a change of mine. It more evenly
2529 // distributes the initial moles amongst multiple multispecies
2530 // phases according to the relative values of the standard state
2531 // free energies. There is no change for problems with one
2532 // multispecies phase. It cut diamond4.vin iterations down from 62
2533 // to 14.
2534 m_deltaMolNumSpecies[kspec] = 0.5 * (m_tPhaseMoles_new[iph] + TMolesMultiphase)
2535 * exp(-m_deltaGRxn_new[irxn]);
2536
2537 for (size_t k = 0; k < m_numComponents; ++k) {
2539 }
2540
2541 for (iph = 0; iph < m_numPhases; iph++) {
2543 }
2544 }
2545 }
2546 if (m_debug_print_lvl >= 2) {
2547 for (size_t kspec = 0; kspec < m_nsp; ++kspec) {
2549 plogf("%sdirection (", pprefix);
2550 plogf("%-12.12s", m_speciesName[kspec]);
2551 plogf(") = %g", m_deltaMolNumSpecies[kspec]);
2552 if (m_SSPhase[kspec]) {
2553 if (m_molNumSpecies_old[kspec] > 0.0) {
2554 plogf(" (ssPhase exists at w = %g moles)", m_molNumSpecies_old[kspec]);
2555 } else {
2556 plogf(" (ssPhase doesn't exist -> stability not checked)");
2557 }
2558 }
2559 writelogendl();
2560 }
2561 }
2562 }
2563
2564 // KEEP COMPONENT SPECIES POSITIVE
2565 double par = 0.5;
2566 for (size_t kspec = 0; kspec < m_numComponents; ++kspec) {
2568 par < -m_deltaMolNumSpecies[kspec] / m_molNumSpecies_new[kspec]) {
2569 par = -m_deltaMolNumSpecies[kspec] / m_molNumSpecies_new[kspec];
2570 }
2571 }
2572 par = 1. / par;
2573 if (par <= 1.0 && par > 0.0) {
2574 par *= 0.8;
2575 } else {
2576 par = 1.0;
2577 }
2578
2579 // CALCULATE NEW MOLE NUMBERS
2580 size_t lt = 0;
2581 size_t ikl = 0;
2582 double s1 = 0.0;
2583 while (true) {
2584 for (size_t kspec = 0; kspec < m_numComponents; ++kspec) {
2586 m_molNumSpecies_old[kspec] = m_molNumSpecies_new[kspec] + par * m_deltaMolNumSpecies[kspec];
2587 } else {
2588 m_deltaMolNumSpecies[kspec] = 0.0;
2589 }
2590 }
2591 for (size_t kspec = m_numComponents; kspec < m_nsp; ++kspec) {
2593 m_deltaMolNumSpecies[kspec] != 0.0) {
2594 m_molNumSpecies_old[kspec] = m_deltaMolNumSpecies[kspec] * par;
2595 }
2596 }
2597
2598 // We have a new w[] estimate, go get the TMoles and m_tPhaseMoles_old[]
2599 // values
2600 vcs_tmoles();
2601 if (lt > 0) {
2602 break;
2603 }
2604
2605 // CONVERGENCE FORCING SECTION
2606 vcs_setFlagsVolPhases(false, VCS_STATECALC_OLD);
2608 double s = 0.0;
2609 for (size_t kspec = 0; kspec < m_nsp; ++kspec) {
2610 s += m_deltaMolNumSpecies[kspec] * m_feSpecies_old[kspec];
2611 }
2612 if (s == 0.0) {
2613 break;
2614 }
2615 if (s < 0.0 && ikl == 0) {
2616 break;
2617 }
2618
2619 // TRY HALF STEP SIZE
2620 if (ikl == 0) {
2621 s1 = s;
2622 par *= 0.5;
2623 ikl = 1;
2624 continue;
2625 }
2626
2627 // FIT PARABOLA THROUGH HALF AND FULL STEPS
2628 double xl = (1.0 - s / (s1 - s)) * 0.5;
2629 if (xl < 0.0) {
2630 // POOR DIRECTION, REDUCE STEP SIZE TO 0.2
2631 par *= 0.2;
2632 } else {
2633 if (xl > 1.0) {
2634 // TOO BIG A STEP, TAKE ORIGINAL FULL STEP
2635 par *= 2.0;
2636 } else {
2637 // ACCEPT RESULTS OF FORCER
2638 par = par * 2.0 * xl;
2639 }
2640 }
2641 lt = 1;
2642 }
2643
2644 if (m_debug_print_lvl >= 2) {
2645 plogf("%s Final Mole Numbers produced by inest:\n",
2646 pprefix);
2647 plogf("%s SPECIES MOLE_NUMBER\n", pprefix);
2648 for (size_t kspec = 0; kspec < m_nsp; ++kspec) {
2649 plogf("%s %-12.12s %g\n",
2650 pprefix, m_speciesName[kspec], m_molNumSpecies_old[kspec]);
2651 }
2652 }
2653
2654 // Calculate the elemental abundances
2655 vcs_elab();
2656
2657 // If we still fail to achieve the correct elemental abundances, try to fix
2658 // the problem again by calling the main elemental abundances fixer routine,
2659 // used in the main program. This attempts to tweak the mole numbers of the
2660 // component species to satisfy the element abundance constraints.
2661 //
2662 // Note: We won't do this unless we have to since it involves inverting a
2663 // matrix.
2664 bool rangeCheck = vcs_elabcheck(1);
2665 if (!vcs_elabcheck(0)) {
2666 if (m_debug_print_lvl >= 2) {
2667 plogf("%sInitial guess failed element abundances\n", pprefix);
2668 plogf("%sCall vcs_elcorr to attempt fix\n", pprefix);
2669 }
2670 vcs_elcorr(sm, aw);
2671 rangeCheck = vcs_elabcheck(1);
2672 if (!vcs_elabcheck(0)) {
2673 throw CanteraError("VCS_SOLVE::vcs_inest",
2674 "Initial guess still fails element abundance equations\n"
2675 "Inability to ever satisfy element abundance constraints is probable");
2676 } else {
2677 if (m_debug_print_lvl >= 2) {
2678 if (rangeCheck) {
2679 plogf("%sInitial guess now satisfies element abundances\n", pprefix);
2680 } else {
2681 plogf("%sElement Abundances RANGE ERROR\n", pprefix);
2682 plogf("%s - Initial guess satisfies NC=%d element abundances, "
2683 "BUT not NE=%d element abundances\n", pprefix,
2685 }
2686 }
2687 }
2688 } else {
2689 if (m_debug_print_lvl >= 2) {
2690 if (rangeCheck) {
2691 plogf("%sInitial guess satisfies element abundances\n", pprefix);
2692 } else {
2693 plogf("%sElement Abundances RANGE ERROR\n", pprefix);
2694 plogf("%s - Initial guess satisfies NC=%d element abundances, "
2695 "BUT not NE=%d element abundances\n", pprefix,
2697 }
2698 }
2699 }
2700
2701 if (m_debug_print_lvl >= 2) {
2702 plogf("%sTotal Dimensionless Gibbs Free Energy = %15.7E\n", pprefix,
2704 }
2705
2707}
2708
2710{
2711 vector<int> numPhSpecies(m_numPhases, 0);
2712 for (size_t kspec = 0; kspec < m_nsp; ++kspec) {
2713 numPhSpecies[m_phaseID[kspec]]++;
2714 }
2715
2716 // Handle the special case of a single species in a phase that has been
2717 // earmarked as a multispecies phase. Treat that species as a single-species
2718 // phase
2719 for (size_t iph = 0; iph < m_numPhases; iph++) {
2720 m_VolPhaseList[iph]->m_singleSpecies = (numPhSpecies[iph] <= 1);
2721 }
2722
2723 // Fill in some useful arrays here that have to do with the static
2724 // information concerning the phase ID of species. SSPhase = Boolean
2725 // indicating whether a species is in a single species phase or not.
2726 for (size_t kspec = 0; kspec < m_nsp; kspec++) {
2727 m_SSPhase[kspec] = m_VolPhaseList[m_phaseID[kspec]]->m_singleSpecies;
2728 }
2729}
2730
2732{
2733 m_VCount->Its = 0;
2734 m_VCount->Basis_Opts = 0;
2735 if (ifunc) {
2736 m_VCount->T_Its = 0;
2740 }
2741}
2742
2744{
2745 plogf("\nTCounters: Num_Calls Total_Its\n");
2746 plogf(" vcs_basopt: %5d %5d\n",
2748 plogf(" vcs_TP: %5d %5d\n",
2750 plogf(" vcs_inest: %5d\n", m_VCount->T_Calls_Inest);
2751}
2752
2753void VCS_SOLVE::prob_report(int print_lvl)
2754{
2755 m_printLvl = print_lvl;
2756
2757 // Printout the species information: PhaseID's and mole nums
2758 if (m_printLvl > 0) {
2759 writeline('=', 80, true, true);
2760 writeline('=', 20, false);
2761 plogf(" VCS_PROB: PROBLEM STATEMENT ");
2762 writeline('=', 31);
2763 writeline('=', 80);
2764 plogf("\n");
2765 plogf("\tSolve a constant T, P problem:\n");
2766 plogf("\t\tT = %g K\n", m_temperature);
2767 double pres_atm = m_pressurePA / 1.01325E5;
2768
2769 plogf("\t\tPres = %g atm\n", pres_atm);
2770 plogf("\n");
2771 plogf(" Phase IDs of species\n");
2772 plogf(" species phaseID phaseName ");
2773 plogf(" Initial_Estimated_Moles Species_Type\n");
2774 for (size_t i = 0; i < m_nsp; i++) {
2775 vcs_VolPhase* Vphase = m_VolPhaseList[m_phaseID[i]].get();
2776 plogf("%16s %5d %16s", m_mix->speciesName(i), m_phaseID[i],
2777 Vphase->PhaseName);
2778 if (m_doEstimateEquil >= 0) {
2779 plogf(" %-10.5g", m_molNumSpecies_old[i]);
2780 } else {
2781 plogf(" N/A");
2782 }
2784 plogf(" Mol_Num");
2786 plogf(" Voltage");
2787 } else {
2788 plogf(" ");
2789 }
2790 plogf("\n");
2791 }
2792
2793 // Printout of the Phase structure information
2794 writeline('-', 80, true, true);
2795 plogf(" Information about phases\n");
2796 plogf(" PhaseName PhaseNum SingSpec EqnState NumSpec "
2797 "Tmoles(kmol)\n");
2798
2799 for (size_t iphase = 0; iphase < m_numPhases; iphase++) {
2800 vcs_VolPhase* Vphase = m_VolPhaseList[iphase].get();
2801 plogf("%16s %8d %8d ", Vphase->PhaseName,
2802 Vphase->VP_ID_, Vphase->m_singleSpecies);
2803 plogf("%16s %8d ", Vphase->eos_name(),
2804 Vphase->nSpecies());
2805 if (m_doEstimateEquil >= 0) {
2806 plogf("%16e\n", Vphase->totalMoles());
2807 } else {
2808 plogf(" N/A\n");
2809 }
2810 }
2811
2812 plogf("\nElemental Abundances: ");
2813 plogf(" Target_kmol ElemType ElActive\n");
2814 for (size_t i = 0; i < m_nelem; ++i) {
2815 writeline(' ', 26, false);
2816 plogf("%-2.2s", m_elementName[i]);
2817 plogf("%20.12E ", m_elemAbundancesGoal[i]);
2818 plogf("%3d %3d\n", m_elType[i], m_elementActive[i]);
2819 }
2820
2821 plogf("\nChemical Potentials: (J/kmol)\n");
2822 plogf(" Species (phase) "
2823 " SS0ChemPot StarChemPot\n");
2824 for (size_t iphase = 0; iphase < m_numPhases; iphase++) {
2825 vcs_VolPhase* Vphase = m_VolPhaseList[iphase].get();
2827 for (size_t kindex = 0; kindex < Vphase->nSpecies(); kindex++) {
2828 size_t kglob = Vphase->spGlobalIndexVCS(kindex);
2829 plogf("%16s ", m_mix->speciesName(kglob));
2830 if (kindex == 0) {
2831 plogf("%16s", Vphase->PhaseName);
2832 } else {
2833 plogf(" ");
2834 }
2835
2836 plogf("%16g %16g\n", Vphase->G0_calc_one(kindex),
2837 Vphase->GStar_calc_one(kindex));
2838 }
2839 }
2840 writeline('=', 80, true, true);
2841 writeline('=', 20, false);
2842 plogf(" VCS_PROB: END OF PROBLEM STATEMENT ");
2843 writeline('=', 24);
2844 writeline('=', 80);
2845 plogf("\n");
2846 }
2847}
2848
2849size_t VCS_SOLVE::addOnePhaseSpecies(vcs_VolPhase* volPhase, size_t k, size_t kT)
2850{
2851 if (kT > m_nsp) {
2852 // Need to expand the number of species here
2853 throw CanteraError("VCS_SOLVE::addOnePhaseSpecies", "Shouldn't be here");
2854 }
2855 const Array2D& fm = volPhase->getFormulaMatrix();
2856 for (size_t eVP = 0; eVP < volPhase->nElemConstraints(); eVP++) {
2857 size_t e = volPhase->elemGlobalIndex(eVP);
2858 AssertThrowMsg(e != npos, "VCS_PROB::addOnePhaseSpecies",
2859 "element not found");
2860 m_formulaMatrix(kT,e) = fm(k,eVP);
2861 }
2862
2863 // Tell the phase object about the current position of the species within
2864 // the global species vector
2865 volPhase->setSpGlobalIndexVCS(k, kT);
2866 return kT;
2867}
2868
2869size_t VCS_SOLVE::addElement(const char* elNameNew, int elType, int elactive)
2870{
2871 if (!elNameNew) {
2872 throw CanteraError("VCS_SOLVE::addElement",
2873 "error: element must have a name");
2874 }
2875 m_nelem++;
2877
2880 m_elType.push_back(elType);
2881 m_elementActive.push_back(elactive);
2882 m_elemAbundances.push_back(0.0);
2883 m_elemAbundancesGoal.push_back(0.0);
2884 m_elementMapIndex.push_back(0);
2885 m_elementName.push_back(elNameNew);
2886 return m_nelem - 1;
2887}
2888
2889size_t VCS_SOLVE::elementIndex(const string& elementName) const
2890{
2891 for (size_t e = 0; e < m_nelem; e++) {
2892 if (m_elementName[e] == elementName) {
2893 return e;
2894 }
2895 }
2896 return npos;
2897}
2898
2899}
Headers for the MultiPhase object that is used to set up multiphase equilibrium problems (see Chemica...
Header file for class ThermoPhase, the base class for phases with thermodynamic properties,...
A class for 2D arrays stored in column-major (Fortran-compatible) form.
Definition Array.h:32
void zero()
Set all of the entries to zero.
Definition Array.h:118
size_t nColumns() const
Number of columns.
Definition Array.h:172
span< double > col(size_t j)
Return a writable span over column j.
Definition Array.h:189
virtual void resize(size_t n, size_t m, double v=0.0)
Resize the array, and fill the new entries with 'v'.
Definition Array.cpp:52
Base class for exceptions thrown by Cantera classes.
A class for full (non-sparse) matrices with Fortran-compatible data storage, which adds matrix operat...
Definition DenseMatrix.h:42
A class for multiphase mixtures.
Definition MultiPhase.h:62
double speciesMoles(size_t kGlob) const
Returns the moles of global species k. units = kmol.
double pressure() const
Pressure [Pa].
Definition MultiPhase.h:352
double temperature() const
Temperature [K].
Definition MultiPhase.h:308
string speciesName(size_t kGlob) const
Name of species with global index kGlob.
size_t nElements() const
Number of elements.
Definition MultiPhase.h:86
ThermoPhase & phase(size_t n)
Return a reference to phase n.
double elementMoles(size_t m) const
Total moles of global element m, summed over all phases.
size_t nSpecies() const
Returns the number of species in the phase.
Definition Phase.h:247
size_t nDim() const
Returns the number of spatial dimensions (1, 2, or 3)
Definition Phase.h:570
double molecularWeight(size_t k) const
Molecular weight of species k.
Definition Phase.cpp:398
double charge(size_t k) const
Dimensionless electrical charge of a single molecule of species k The charge is normalized by the the...
Definition Phase.h:562
Base class for a phase with thermodynamic properties.
double electricPotential() const
Returns the electric potential of this phase (V).
Class to keep track of time and iterations.
int Its
Current number of iterations in the main loop of vcs_TP() to solve for thermo equilibrium.
int T_Calls_Inest
Current number of times the initial thermo equilibrium estimator has been called.
int T_Basis_Opts
Total number of optimizations of the components basis set done.
int Basis_Opts
number of optimizations of the components basis set done
int T_Its
Total number of iterations in the main loop of vcs_TP() to solve for thermo equilibrium.
int T_Calls_vcs_TP
Current number of calls to vcs_TP.
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
VCS_SOLVE(MultiPhase *mphase, int printLvl=0)
Initialize the sizes within the VCS_SOLVE object.
Definition vcs_solve.cpp:40
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
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_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_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.
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.
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
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
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.
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
void vcs_elem_rearrange()
Rearrange the constraint equations represented by the Formula Matrix so that the operational ones are...
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
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...
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.
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
vector< double > m_deltaGRxn_old
Last deltag[irxn] from the previous step.
Definition vcs_solve.h:1003
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.
int speciesUnknownType(const size_t k) const
Returns the type of the species unknown.
size_t elemGlobalIndex(const size_t e) const
Returns the global index of the local element index for the phase.
void setMolesFromVCS(const int stateCalc, span< const double > molesSpeciesVCS={})
Set the moles within the phase.
size_t nSpecies() const
Return the number of species in the phase.
string PhaseName
String name for the phase.
size_t VP_ID_
Original ID of the phase in the problem.
const Array2D & getFormulaMatrix() const
Get a constant form of the Species Formula Matrix.
void sendToVCS_ActCoeff(const int stateCalc, span< double > AC)
Fill in an activity coefficients vector within a VCS_SOLVE object.
size_t phiVarIndex() const
Return the index of the species that represents the the voltage of the phase.
bool m_singleSpecies
If true, this phase consists of a single species.
void setElemGlobalIndex(const size_t eLocal, const size_t eGlobal)
sets a local phase element to a global index value
int exists() const
Retrieve the kth Species structure for the species belonging to this phase.
size_t nElemConstraints() const
Returns the number of element constraints.
size_t spGlobalIndexVCS(const size_t spIndex) const
Return the Global VCS index of the kth species in the phase.
double GStar_calc_one(size_t kspec) const
Gibbs free energy calculation for standard state of one species.
int p_activityConvention
Convention for the activity formulation.
void setCreationMoleNumbers(span< const double > n_k, span< const size_t > creationGlobalRxnNumbers)
Sets the creationMoleNum's within the phase object.
span< const double > creationMoleNumbers(span< size_t > creationGlobalRxnNumbers) const
Return a const reference to the creationMoleNumbers stored in the object.
double G0_calc_one(size_t kspec) const
Gibbs free energy calculation at a temperature for the reference state of a species,...
void setState_TP(const double temperature_Kelvin, const double pressure_PA)
Sets the temperature and pressure in this object and underlying ThermoPhase objects.
string eos_name() const
Return the name corresponding to the equation of state.
void setMoleFractionsState(const double molNum, span< const double > moleFracVec, const int vcsStateStatus)
Set the moles and/or mole fractions within the phase.
double sendToVCS_VolPM(span< double > VolPM) const
Fill in the partial molar volume vector for VCS.
void sendToVCS_GStar(span< double > gstar) const
Fill in the standard state Gibbs free energy vector for VCS.
void setSpGlobalIndexVCS(const size_t spIndex, const size_t spGlobalIndex)
set the Global VCS index of the kth species in the phase
void setExistence(const int existence)
Set the existence flag in the object.
double totalMoles() const
Return the total moles in the phase.
Definitions for the classes that are thrown when Cantera experiences an error condition (also contain...
#define AssertThrowMsg(expr, procedure,...)
Assertion must be true or an error is thrown.
void debuglog(const string &msg, int loglevel)
Write a message to the log only if loglevel > 0.
Definition global.h:154
void writelogendl()
Write an end of line character to the screen and flush output.
Definition global.cpp:41
void scale(InputIter begin, InputIter end, OutputIter out, S scale_factor)
Multiply elements of an array by a scale factor.
Definition utilities.h:118
T clip(const T &value, const T &lower, const T &upper)
Clip value such that lower <= value <= upper.
Definition global.h:326
const double Faraday
Faraday constant [C/kmol].
Definition ct_defs.h:134
const double GasConstant
Universal Gas Constant [J/kmol/K].
Definition ct_defs.h:123
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595
const size_t npos
index returned by functions to indicate "no position"
Definition ct_defs.h:183
void solve(DenseMatrix &A, span< double > b, size_t nrhs, size_t ldb)
Solve Ax = b. Array b is overwritten on exit with x.
bool vcs_doubleEqual(double d1, double d2)
Simple routine to check whether two doubles are equal up to roundoff error.
Definition vcs_util.cpp:89
double vcs_l2norm(span< const double > vec)
determine the l2 norm of a vector of doubles
Definition vcs_util.cpp:19
void checkArraySize(const char *procedure, size_t available, size_t required)
Wrapper for throwing ArraySizeError.
Contains const definitions for types of species reference-state thermodynamics managers (see Species ...
Contains declarations for string manipulation functions within Cantera.
Header for the object representing each phase within vcs.
#define VCS_SPECIES_INTERFACIALVOLTAGE
Species refers to an electron in the metal.
Definition vcs_defs.h:101
#define VCS_SPECIES_TYPE_INTERFACIALVOLTAGE
Unknown refers to the voltage level of a phase.
Definition vcs_defs.h:197
#define VCS_STATECALC_OLD
State Calculation based on the old or base mole numbers.
Definition vcs_defs.h:206
#define VCS_DELETE_MINORSPECIES_CUTOFF
Cutoff relative mole number value, below which species are deleted from the equilibrium problem.
Definition vcs_defs.h:32
#define VCS_SPECIES_COMPONENT
These defines are valid values for spStatus()
Definition vcs_defs.h:58
#define VCS_SPECIES_MAJOR
Species is a major species.
Definition vcs_defs.h:65
#define VCS_PHASE_EXIST_NO
Phase doesn't currently exist in the mixture.
Definition vcs_defs.h:147
#define VCS_PHASE_EXIST_ALWAYS
These defines are valid values for the phase existence flag.
Definition vcs_defs.h:141
#define VCS_SPECIES_STOICHZERO
Species lies in a multicomponent phase that is active, but species concentration is zero due to stoic...
Definition vcs_defs.h:131
#define VCS_ELEM_TYPE_ABSPOS
Normal element constraint consisting of positive coefficients for the formula matrix.
Definition vcs_defs.h:170
#define VCS_STATECALC_NEW
State Calculation based on the new or tentative mole numbers.
Definition vcs_defs.h:209
#define VCS_SPECIES_TYPE_MOLNUM
Unknown refers to mole number of a single species.
Definition vcs_defs.h:189
#define VCS_SPECIES_ZEROEDPHASE
Species lies in a multicomponent phase that is zeroed atm.
Definition vcs_defs.h:109
#define VCS_SPECIES_ACTIVEBUTZERO
Species lies in a multicomponent phase that is active, but species concentration is zero.
Definition vcs_defs.h:119
#define VCS_ELEM_TYPE_ELECTRONCHARGE
This refers to conservation of electrons.
Definition vcs_defs.h:176
#define VCS_SPECIES_ZEROEDSS
Species is a SS phase, that is currently zeroed out.
Definition vcs_defs.h:85
#define VCS_DELETE_PHASE_CUTOFF
Cutoff relative moles below which a phase is deleted from the equilibrium problem.
Definition vcs_defs.h:40
#define VCS_SMALL_MULTIPHASE_SPECIES
Relative value of multiphase species mole number for a multiphase species which is small.
Definition vcs_defs.h:36
#define VCS_ELEM_TYPE_CHARGENEUTRALITY
This refers to a charge neutrality of a single phase.
Definition vcs_defs.h:182
#define VCS_SPECIES_MINOR
Species is a major species.
Definition vcs_defs.h:72
#define VCS_DELETE_ELEMENTABS_CUTOFF
Cutoff moles below which a phase or species which comprises the bulk of an element's total concentrat...
Definition vcs_defs.h:44
#define VCS_PHASE_EXIST_YES
Phase is a normal phase that currently exists.
Definition vcs_defs.h:144
#define VCS_SUCCESS
ERROR CODES.
Definition vcs_defs.h:18
#define VCS_STATECALC_PHASESTABILITY
State Calculation based on tentative mole numbers for a phase which is currently zeroed,...
Definition vcs_defs.h:214
#define plogf
define this Cantera function to replace printf
Header file for the internal object that holds the vcs equilibrium problem (see Class VCS_SOLVE and C...