Cantera  3.3.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
420size_t VCS_SOLVE::vcs_popPhaseID(vector<size_t> & phasePopPhaseIDs)
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 size_t irxn = kspec - m_numComponents;
446 if (irxn > m_deltaGRxn_old.size()) {
447 throw CanteraError("VCS_SOLVE::vcs_popPhaseID",
448 "Index out of bounds due to logic error.");
449 }
450 double deltaGRxn = m_deltaGRxn_old[irxn];
451 Fephase = exp(-deltaGRxn) - 1.0;
452 if (Fephase > 0.0) {
453 note = "(ready to be birthed)";
454 if (Fephase > FephaseMax) {
455 iphasePop = iph;
456 FephaseMax = Fephase;
457 note = "(chosen to be birthed)";
458 }
459 }
460 if (Fephase < 0.0) {
461 note = "(not stable)";
463 "VCS_SOLVE::vcs_popPhaseID", "shouldn't be here");
464 }
465
466 if (m_debug_print_lvl >= 2) {
467 plogf(" --- %18s %5d %10.3g %10.3g %s\n",
468 Vphase->PhaseName, existence, Fephase,
469 m_tPhaseMoles_old[iph], note.c_str());
470 }
471 } else if (vcs_popPhasePossible(iph)) {
472 // MultiSpecies Phase Stability Resolution
473 Fephase = vcs_phaseStabilityTest(iph);
474 if (Fephase > 0.0) {
475 if (Fephase > FephaseMax) {
476 iphasePop = iph;
477 FephaseMax = Fephase;
478 }
479 } else {
480 FephaseMax = std::max(FephaseMax, Fephase);
481 }
482 if (m_debug_print_lvl >= 2) {
483 plogf(" --- %18s %5d %11.3g %11.3g\n",
484 Vphase->PhaseName, existence, Fephase,
485 m_tPhaseMoles_old[iph]);
486 }
487 } else if (m_debug_print_lvl >= 2) {
488 plogf(" --- %18s %5d blocked %11.3g\n",
489 Vphase->PhaseName,
490 existence, m_tPhaseMoles_old[iph]);
491 }
492 }
493 }
494 phasePopPhaseIDs.resize(0);
495 if (iphasePop != npos) {
496 phasePopPhaseIDs.push_back(iphasePop);
497 }
498
499 // Insert logic here to figure out if phase pops are linked together. Only
500 // do one linked pop at a time.
501 if (m_debug_print_lvl >= 2) {
502 plogf(" ---------------------------------------------------------------------\n");
503 }
504 return iphasePop;
505}
506
507int VCS_SOLVE::vcs_popPhaseRxnStepSizes(const size_t iphasePop)
508{
509 vcs_VolPhase* Vphase = m_VolPhaseList[iphasePop].get();
510 // Identify the first species in the phase
511 size_t kspec = Vphase->spGlobalIndexVCS(0);
512 // Identify the formation reaction for that species
513 size_t irxn = kspec - m_numComponents;
514 vector<size_t> creationGlobalRxnNumbers;
515
516 // Calculate the initial moles of the phase being born.
517 // Here we set it to 10x of the value which would cause the phase to be
518 // zeroed out within the algorithm. We may later adjust the value.
519 double tPhaseMoles = 10. * m_totalMolNum * VCS_DELETE_PHASE_CUTOFF;
520
521 AssertThrowMsg(!Vphase->exists(), "VCS_SOLVE::vcs_popPhaseRxnStepSizes",
522 "called for a phase that exists!");
523 if (m_debug_print_lvl >= 2) {
524 plogf(" --- vcs_popPhaseRxnStepSizes() called to pop phase %s %d into existence\n",
525 Vphase->PhaseName, iphasePop);
526 }
527 // Section for a single-species phase
528 if (Vphase->m_singleSpecies) {
529 double s = 0.0;
530 for (size_t j = 0; j < m_numComponents; ++j) {
531 if (!m_SSPhase[j] && m_molNumSpecies_old[j] > 0.0) {
532 s += pow(m_stoichCoeffRxnMatrix(j,irxn), 2) / m_molNumSpecies_old[j];
533 }
534 }
535 for (size_t j = 0; j < m_numPhases; j++) {
536 Vphase = m_VolPhaseList[j].get();
537 if (! Vphase->m_singleSpecies && m_tPhaseMoles_old[j] > 0.0) {
538 s -= pow(m_deltaMolNumPhase(j,irxn), 2) / m_tPhaseMoles_old[j];
539 }
540 }
541 if (s != 0.0) {
542 s = vcs_Hessian_diag_adj(irxn, s);
543 m_deltaMolNumSpecies[kspec] = -m_deltaGRxn_new[irxn] / s;
544 } else {
545 // Ok, s is equal to zero. We can not apply a sophisticated theory
546 // to birth the phase. Just pick a small delta and go with it.
547 m_deltaMolNumSpecies[kspec] = tPhaseMoles;
548 }
549
550 // section to do damping of the m_deltaMolNumSpecies[]
551 for (size_t j = 0; j < m_numComponents; ++j) {
552 double stoicC = m_stoichCoeffRxnMatrix(j,irxn);
553 if (stoicC != 0.0 && m_elType[j] == VCS_ELEM_TYPE_ABSPOS) {
554 double negChangeComp = - stoicC * m_deltaMolNumSpecies[kspec];
555 if (negChangeComp > m_molNumSpecies_old[j]) {
556 if (m_molNumSpecies_old[j] > 0.0) {
557 m_deltaMolNumSpecies[kspec] = - 0.5 * m_molNumSpecies_old[j] / stoicC;
558 } else {
559 m_deltaMolNumSpecies[kspec] = 0.0;
560 }
561 }
562 }
563 }
564 // Implement a damping term that limits m_deltaMolNumSpecies to the size
565 // of the mole number
566 if (-m_deltaMolNumSpecies[kspec] > m_molNumSpecies_old[kspec]) {
568 }
569 } else {
570 vector<double> fracDelta(Vphase->nSpecies());
571 vector<double> X_est(Vphase->nSpecies());
572 fracDelta = Vphase->creationMoleNumbers(creationGlobalRxnNumbers);
573
574 double sumFrac = 0.0;
575 for (size_t k = 0; k < Vphase->nSpecies(); k++) {
576 sumFrac += fracDelta[k];
577 }
578 for (size_t k = 0; k < Vphase->nSpecies(); k++) {
579 X_est[k] = fracDelta[k] / sumFrac;
580 }
581
582 double deltaMolNumPhase = tPhaseMoles;
583 double damp = 1.0;
585 double* molNumSpecies_tmp = m_deltaGRxn_tmp.data();
586
587 for (size_t k = 0; k < Vphase->nSpecies(); k++) {
588 kspec = Vphase->spGlobalIndexVCS(k);
589 double delmol = deltaMolNumPhase * X_est[k];
590 if (kspec >= m_numComponents) {
591 irxn = kspec - m_numComponents;
592 if (irxn > m_stoichCoeffRxnMatrix.nColumns()) {
593 throw CanteraError("VCS_SOLVE::vcs_popPhaseRxnStepSizes",
594 "Index out of bounds due to logic error.");
595 }
596 for (size_t j = 0; j < m_numComponents; ++j) {
597 double stoicC = m_stoichCoeffRxnMatrix(j,irxn);
598 if (stoicC != 0.0 && m_elType[j] == VCS_ELEM_TYPE_ABSPOS) {
599 molNumSpecies_tmp[j] += stoicC * delmol;
600 }
601 }
602 }
603 }
604
605 double ratioComp = 0.0;
606 for (size_t j = 0; j < m_numComponents; ++j) {
607 double deltaJ = m_molNumSpecies_old[j] - molNumSpecies_tmp[j];
608 if (molNumSpecies_tmp[j] < 0.0) {
609 ratioComp = 1.0;
610 if (deltaJ > 0.0) {
611 double delta0 = m_molNumSpecies_old[j];
612 damp = std::min(damp, delta0 / deltaJ * 0.9);
613 }
614 } else {
615 if (m_elType[j] == VCS_ELEM_TYPE_ABSPOS) {
616 size_t jph = m_phaseID[j];
617 if ((jph != iphasePop) && (!m_SSPhase[j])) {
618 double fdeltaJ = fabs(deltaJ);
619 if (m_molNumSpecies_old[j] > 0.0) {
620 ratioComp = std::max(ratioComp, fdeltaJ/ m_molNumSpecies_old[j]);
621 }
622 }
623 }
624 }
625 }
626
627 // We may have greatly underestimated the deltaMoles for the phase pop
628 // Here we create a damp > 1 to account for this possibility. We adjust
629 // upwards to make sure that a component in an existing multispecies
630 // phase is modified by a factor of 1/1000.
631 if (ratioComp > 1.0E-30 && ratioComp < 0.001) {
632 damp = 0.001 / ratioComp;
633 }
634 if (damp <= 1.0E-6) {
635 return 3;
636 }
637
638 for (size_t k = 0; k < Vphase->nSpecies(); k++) {
639 kspec = Vphase->spGlobalIndexVCS(k);
640 if (kspec < m_numComponents) {
642 } else {
643 m_deltaMolNumSpecies[kspec] = deltaMolNumPhase * X_est[k] * damp;
644 m_speciesStatus[kspec] = (X_est[k] > 1.0E-3) ? VCS_SPECIES_MAJOR
646 }
647 }
648 }
649 return 0;
650}
651
652size_t VCS_SOLVE::vcs_RxnStepSizes(int& forceComponentCalc, size_t& kSpecial)
653{
654 size_t iphDel = npos;
655 size_t k = 0;
656 if (m_debug_print_lvl >= 2) {
657 plogf(" ");
658 for (int j = 0; j < 82; j++) {
659 plogf("-");
660 }
661 plogf("\n");
662 plogf(" --- Subroutine vcs_RxnStepSizes called - Details:\n");
663 plogf(" ");
664 for (int j = 0; j < 82; j++) {
665 plogf("-");
666 }
667 plogf("\n");
668 plogf(" --- Species KMoles Rxn_Adjustment DeltaG\n");
669 }
670
671 // We update the matrix dlnActCoeffdmolNumber[][] at the top of the loop,
672 // when necessary
673 if (m_useActCoeffJac) {
675 }
676
677 // LOOP OVER THE FORMATION REACTIONS
678 for (size_t irxn = 0; irxn < m_numRxnRdc; ++irxn) {
679
680 size_t kspec = m_indexRxnToSpecies[irxn];
682 m_deltaMolNumSpecies[kspec] = 0.0;
684 if (m_molNumSpecies_old[kspec] == 0.0 && (!m_SSPhase[kspec])) {
685 // MULTISPECIES PHASE WITH total moles equal to zero
686 //
687 // If dg[irxn] is negative, then the multispecies phase should
688 // come alive again. Add a small positive step size to make it
689 // come alive.
690 if (m_deltaGRxn_new[irxn] < -1.0e-4) {
691 // First decide if this species is part of a multiphase that
692 // is nontrivial in size.
693 size_t iph = m_phaseID[kspec];
694 double tphmoles = m_tPhaseMoles_old[iph];
695 double trphmoles = tphmoles / m_totalMolNum;
696 vcs_VolPhase* Vphase = m_VolPhaseList[iph].get();
697 if (Vphase->exists() && (trphmoles > VCS_DELETE_PHASE_CUTOFF)) {
700 m_deltaMolNumSpecies[kspec] = 0.0;
701 } else {
703 }
704 } else {
705 m_deltaMolNumSpecies[kspec] = 0.0;
706 if (Vphase->exists() > 0 && trphmoles > 0.0) {
708 }
709 }
710 } else {
711 m_deltaMolNumSpecies[kspec] = 0.0;
712 }
713 } else {
714 // REGULAR PROCESSING
715 //
716 // First take care of cases where we want to bail out. Don't
717 // bother if superconvergence has already been achieved in this
718 // mode.
719 if (fabs(m_deltaGRxn_new[irxn]) <= m_tolmaj2) {
720 if (m_debug_print_lvl >= 2) {
721 plogf(" --- %-12.12s", m_speciesName[kspec]);
722 plogf(" %12.4E %12.4E %12.4E\n",
724 m_deltaGRxn_new[irxn]);
725 }
726 continue;
727 }
728
729 // Don't calculate for minor or nonexistent species if their
730 // values are to be decreasing anyway.
731 if ((m_speciesStatus[kspec] != VCS_SPECIES_MAJOR) && (m_deltaGRxn_new[irxn] >= 0.0)) {
732 if (m_debug_print_lvl >= 2) {
733 plogf(" --- %-12.12s", m_speciesName[kspec]);
734 plogf(" %12.4E %12.4E %12.4E\n",
736 m_deltaGRxn_new[irxn]);
737 }
738 continue;
739 }
740
741 // Start of the regular processing
742 double s = (m_SSPhase[kspec]) ? 0.0 : 1.0 / m_molNumSpecies_old[kspec];
743 for (size_t j = 0; j < m_numComponents; ++j) {
744 if (!m_SSPhase[j] && m_molNumSpecies_old[j] > 0.0) {
745 s += pow(m_stoichCoeffRxnMatrix(j,irxn), 2) / m_molNumSpecies_old[j];
746 }
747 }
748 for (size_t j = 0; j < m_numPhases; j++) {
749 vcs_VolPhase* Vphase = m_VolPhaseList[j].get();
750 if (!Vphase->m_singleSpecies && m_tPhaseMoles_old[j] > 0.0) {
751 s -= pow(m_deltaMolNumPhase(j,irxn), 2) / m_tPhaseMoles_old[j];
752 }
753 }
754 if (s != 0.0) {
755 // Take into account of the derivatives of the activity
756 // coefficients with respect to the mole numbers, even in
757 // our diagonal approximation.
758 if (m_useActCoeffJac) {
759 double s_old = s;
760 s = vcs_Hessian_diag_adj(irxn, s_old);
761 }
762
763 m_deltaMolNumSpecies[kspec] = -m_deltaGRxn_new[irxn] / s;
764 // New section to do damping of the m_deltaMolNumSpecies[]
765 for (size_t j = 0; j < m_numComponents; ++j) {
766 double stoicC = m_stoichCoeffRxnMatrix(j,irxn);
767 if (stoicC != 0.0) {
768 double negChangeComp = -stoicC * m_deltaMolNumSpecies[kspec];
769 if (negChangeComp > m_molNumSpecies_old[j]) {
770 if (m_molNumSpecies_old[j] > 0.0) {
771 m_deltaMolNumSpecies[kspec] = -m_molNumSpecies_old[j] / stoicC;
772 } else {
773 m_deltaMolNumSpecies[kspec] = 0.0;
774 }
775 }
776 }
777 }
778 // Implement a damping term that limits m_deltaMolNumSpecies
779 // to the size of the mole number
780 if (-m_deltaMolNumSpecies[kspec] > m_molNumSpecies_old[kspec]) {
782 }
783 } else {
784 // REACTION IS ENTIRELY AMONGST SINGLE SPECIES PHASES.
785 // DELETE ONE OF THE PHASES AND RECOMPUTE BASIS.
786 //
787 // Either the species L will disappear or one of the
788 // component single species phases will disappear. The sign
789 // of DG(I) will indicate which way the reaction will go.
790 // Then, we need to follow the reaction to see which species
791 // will zero out first. The species to be zeroed out will be
792 // "k".
793 double dss;
794 if (m_deltaGRxn_new[irxn] > 0.0) {
795 dss = m_molNumSpecies_old[kspec];
796 k = kspec;
797 for (size_t j = 0; j < m_numComponents; ++j) {
798 if (m_stoichCoeffRxnMatrix(j,irxn) > 0.0) {
799 double xx = m_molNumSpecies_old[j] / m_stoichCoeffRxnMatrix(j,irxn);
800 if (xx < dss) {
801 dss = xx;
802 k = j;
803 }
804 }
805 }
806 dss = -dss;
807 } else {
808 dss = 1.0e10;
809 for (size_t j = 0; j < m_numComponents; ++j) {
810 if (m_stoichCoeffRxnMatrix(j,irxn) < 0.0) {
811 double xx = -m_molNumSpecies_old[j] / m_stoichCoeffRxnMatrix(j,irxn);
812 if (xx < dss) {
813 dss = xx;
814 k = j;
815 }
816 }
817 }
818 }
819
820 // Here we adjust the mole fractions according to DSS and
821 // the stoichiometric array to take into account that we are
822 // eliminating the kth species. DSS contains the amount of
823 // moles of the kth species that needs to be added back into
824 // the component species.
825 if (dss != 0.0) {
826 if ((k == kspec) && (m_SSPhase[kspec] != 1)) {
827 // Found out that we can be in this spot, when
828 // components of multispecies phases are zeroed,
829 // leaving noncomponent species of the same phase
830 // having all of the mole numbers of that phases. it
831 // seems that we can suggest a zero of the species
832 // and the code will recover.
834 if (m_debug_print_lvl >= 2) {
835 plogf(" --- %-12.12s", m_speciesName[kspec]);
836 plogf(" %12.4E %12.4E %12.4E\n",
838 m_deltaGRxn_new[irxn]);
839 }
840 continue;
841 }
842
843 // Delete the single species phase
844 for (size_t j = 0; j < m_nsp; j++) {
845 m_deltaMolNumSpecies[j] = 0.0;
846 }
847 m_deltaMolNumSpecies[kspec] = dss;
848 for (size_t j = 0; j < m_numComponents; ++j) {
850 }
851
852 iphDel = m_phaseID[k];
853 kSpecial = k;
854
855 if (k != kspec) {
856 } else {
857 }
858 if (m_debug_print_lvl >= 2) {
859 plogf(" --- %-12.12s", m_speciesName[kspec]);
860 plogf(" %12.4E %12.4E %12.4E\n",
862 m_deltaGRxn_new[irxn]);
863 plogf(" --- vcs_RxnStepSizes Special section to set up to delete %s\n",
864 m_speciesName[k]);
865 }
866 if (k != kspec) {
867 forceComponentCalc = 1;
868 debuglog(" --- Force a component recalculation\n\n", m_debug_print_lvl >= 2);
869 }
870 if (m_debug_print_lvl >= 2) {
871 plogf(" ");
872 writeline('-', 82);
873 }
874 return iphDel;
875 }
876 }
877 } // End of regular processing
878 if (m_debug_print_lvl >= 2) {
879 plogf(" --- %-12.12s", m_speciesName[kspec]);
880 plogf(" %12.4E %12.4E %12.4E\n",
882 m_deltaGRxn_new[irxn]);
883 }
884 } // End of loop over m_speciesUnknownType
885 } // End of loop over non-component stoichiometric formation reactions
886 if (m_debug_print_lvl >= 2) {
887 plogf(" ");
888 writeline('-', 82);
889 }
890 return iphDel;
891}
892
893double VCS_SOLVE::vcs_phaseStabilityTest(const size_t iph)
894{
895 // We will use the _new state calc here
896 vcs_VolPhase* Vphase = m_VolPhaseList[iph].get();
897 const size_t nsp = Vphase->nSpecies();
898 int minNumberIterations = 3;
899 if (nsp <= 1) {
900 minNumberIterations = 1;
901 }
902
903 // We will do a full Newton calculation later, but for now, ...
904 bool doSuccessiveSubstitution = true;
905 double funcPhaseStability;
906 vector<double> X_est(nsp, 0.0);
907 vector<double> delFrac(nsp, 0.0);
908 vector<double> E_phi(nsp, 0.0);
909 vector<double> fracDelta_old(nsp, 0.0);
910 vector<double> fracDelta_raw(nsp, 0.0);
911 vector<size_t> creationGlobalRxnNumbers(nsp, npos);
913 vector<double> feSpecies_Deficient = m_feSpecies_old;
914
915 // get the activity coefficients
917
918 // Get the stored estimate for the composition of the phase if
919 // it gets created
920 vector<double> fracDelta_new = Vphase->creationMoleNumbers(creationGlobalRxnNumbers);
921
922 vector<size_t> componentList;
923 for (size_t k = 0; k < nsp; k++) {
924 size_t kspec = Vphase->spGlobalIndexVCS(k);
925 if (kspec < m_numComponents) {
926 componentList.push_back(k);
927 }
928 }
929
930 double normUpdate = 0.1 * vcs_l2norm(fracDelta_new);
931 double damp = 1.0E-2;
932
933 if (doSuccessiveSubstitution) {
934 int KP = 0;
935 if (m_debug_print_lvl >= 2) {
936 plogf(" --- vcs_phaseStabilityTest() called\n");
937 plogf(" --- Its X_old[%2d] FracDel_old[%2d] deltaF[%2d] FracDel_new[%2d]"
938 " normUpdate damp FuncPhaseStability\n", KP, KP, KP, KP);
939 plogf(" --------------------------------------------------------------"
940 "--------------------------------------------------------\n");
941 } else if (m_debug_print_lvl == 1) {
942 plogf(" --- vcs_phaseStabilityTest() called for phase %d\n", iph);
943 }
944
945 for (size_t k = 0; k < nsp; k++) {
946 if (fracDelta_new[k] < 1.0E-13) {
947 fracDelta_new[k] =1.0E-13;
948 }
949 }
950 bool converged = false;
951 double dirProd = 0.0;
952 for (int its = 0; its < 200 && (!converged); its++) {
953 double dampOld = damp;
954 double normUpdateOld = normUpdate;
955 fracDelta_old = fracDelta_new;
956 double dirProdOld = dirProd;
957
958 // Given a set of fracDelta's, we calculate the fracDelta's
959 // for the component species, if any
960 for (size_t i = 0; i < componentList.size(); i++) {
961 size_t kc = componentList[i];
962 size_t kc_spec = Vphase->spGlobalIndexVCS(kc);
963 fracDelta_old[kc] = 0.0;
964 for (size_t k = 0; k < nsp; k++) {
965 size_t kspec = Vphase->spGlobalIndexVCS(k);
966 if (kspec >= m_numComponents) {
967 size_t irxn = kspec - m_numComponents;
968 fracDelta_old[kc] += m_stoichCoeffRxnMatrix(kc_spec,irxn) * fracDelta_old[k];
969 }
970 }
971 }
972
973 // Now, calculate the predicted mole fractions, X_est[k]
974 double sumFrac = 0.0;
975 for (size_t k = 0; k < nsp; k++) {
976 sumFrac += fracDelta_old[k];
977 }
978 // Necessary because this can be identically zero. -> we need to fix
979 // this algorithm!
980 if (sumFrac <= 0.0) {
981 sumFrac = 1.0;
982 }
983 double sum_Xcomp = 0.0;
984 for (size_t k = 0; k < nsp; k++) {
985 X_est[k] = fracDelta_old[k] / sumFrac;
986 if (Vphase->spGlobalIndexVCS(k) < m_numComponents) {
987 sum_Xcomp += X_est[k];
988 }
989 }
990
991 // Feed the newly formed estimate of the mole fractions back into the
992 // ThermoPhase object
994
995 // get the activity coefficients
997
998 // First calculate altered chemical potentials for component species
999 // belonging to this phase.
1000 for (size_t i = 0; i < componentList.size(); i++) {
1001 size_t kc = componentList[i];
1002 size_t kc_spec = Vphase->spGlobalIndexVCS(kc);
1003 double moleFrac = std::max(X_est[kc], VCS_DELETE_MINORSPECIES_CUTOFF);
1004 feSpecies_Deficient[kc_spec] = m_feSpecies_old[kc_spec]
1005 + log(m_actCoeffSpecies_new[kc_spec] * moleFrac);
1006 }
1007
1008 for (size_t i = 0; i < componentList.size(); i++) {
1009 size_t kc_spec = Vphase->spGlobalIndexVCS(componentList[i]);
1010 for (size_t k = 0; k < Vphase->nSpecies(); k++) {
1011 size_t kspec = Vphase->spGlobalIndexVCS(k);
1012 if (kspec >= m_numComponents) {
1013 size_t irxn = kspec - m_numComponents;
1014 if (i == 0) {
1016 }
1017 if (m_stoichCoeffRxnMatrix(kc_spec,irxn) != 0.0) {
1018 m_deltaGRxn_Deficient[irxn] +=
1019 m_stoichCoeffRxnMatrix(kc_spec,irxn) * (feSpecies_Deficient[kc_spec]- m_feSpecies_old[kc_spec]);
1020 }
1021 }
1022 }
1023 }
1024
1025 // Calculate the E_phi's
1026 double sum = 0.0;
1027 funcPhaseStability = sum_Xcomp - 1.0;
1028 for (size_t k = 0; k < nsp; k++) {
1029 size_t kspec = Vphase->spGlobalIndexVCS(k);
1030 if (kspec >= m_numComponents) {
1031 size_t irxn = kspec - m_numComponents;
1032 double deltaGRxn = clip(m_deltaGRxn_Deficient[irxn], -50.0, 50.0);
1033 E_phi[k] = std::exp(-deltaGRxn) / m_actCoeffSpecies_new[kspec];
1034 sum += E_phi[k];
1035 funcPhaseStability += E_phi[k];
1036 } else {
1037 E_phi[k] = 0.0;
1038 }
1039 }
1040
1041 // Calculate the raw estimate of the new fracs
1042 for (size_t k = 0; k < nsp; k++) {
1043 size_t kspec = Vphase->spGlobalIndexVCS(k);
1044 double b = E_phi[k] / sum * (1.0 - sum_Xcomp);
1045 if (kspec >= m_numComponents) {
1046 fracDelta_raw[k] = b;
1047 }
1048 }
1049
1050 // Given a set of fracDelta's, we calculate the fracDelta's
1051 // for the component species, if any
1052 for (size_t i = 0; i < componentList.size(); i++) {
1053 size_t kc = componentList[i];
1054 size_t kc_spec = Vphase->spGlobalIndexVCS(kc);
1055 fracDelta_raw[kc] = 0.0;
1056 for (size_t k = 0; k < nsp; k++) {
1057 size_t kspec = Vphase->spGlobalIndexVCS(k);
1058 if (kspec >= m_numComponents) {
1059 size_t irxn = kspec - m_numComponents;
1060 fracDelta_raw[kc] += m_stoichCoeffRxnMatrix(kc_spec,irxn) * fracDelta_raw[k];
1061 }
1062 }
1063 }
1064
1065 // Now possibly dampen the estimate.
1066 for (size_t k = 0; k < nsp; k++) {
1067 delFrac[k] = fracDelta_raw[k] - fracDelta_old[k];
1068 }
1069 normUpdate = vcs_l2norm(delFrac);
1070
1071 dirProd = 0.0;
1072 for (size_t k = 0; k < nsp; k++) {
1073 dirProd += fracDelta_old[k] * delFrac[k];
1074 }
1075 bool crossedSign = false;
1076 if (dirProd * dirProdOld < 0.0) {
1077 crossedSign = true;
1078 }
1079
1080 damp = 0.5;
1081 if (dampOld < 0.25) {
1082 damp = 2.0 * dampOld;
1083 }
1084 if (crossedSign) {
1085 if (normUpdate *1.5 > normUpdateOld) {
1086 damp = 0.5 * dampOld;
1087 } else if (normUpdate *2.0 > normUpdateOld) {
1088 damp = 0.8 * dampOld;
1089 }
1090 } else {
1091 if (normUpdate > normUpdateOld * 2.0) {
1092 damp = 0.6 * dampOld;
1093 } else if (normUpdate > normUpdateOld * 1.2) {
1094 damp = 0.9 * dampOld;
1095 }
1096 }
1097
1098 for (size_t k = 0; k < nsp; k++) {
1099 if (fabs(damp * delFrac[k]) > 0.3*fabs(fracDelta_old[k])) {
1100 damp = std::max(0.3*fabs(fracDelta_old[k]) / fabs(delFrac[k]),
1101 1.0E-8/fabs(delFrac[k]));
1102 }
1103 if (delFrac[k] < 0.0 && 2.0 * damp * (-delFrac[k]) > fracDelta_old[k]) {
1104 damp = fracDelta_old[k] / (2.0 * -delFrac[k]);
1105 }
1106 if (delFrac[k] > 0.0 && 2.0 * damp * delFrac[k] > fracDelta_old[k]) {
1107 damp = fracDelta_old[k] / (2.0 * delFrac[k]);
1108 }
1109 }
1110 damp = std::max(damp, 0.000001);
1111 for (size_t k = 0; k < nsp; k++) {
1112 fracDelta_new[k] = fracDelta_old[k] + damp * delFrac[k];
1113 }
1114
1115 if (m_debug_print_lvl >= 2) {
1116 plogf(" --- %3d %12g %12g %12g %12g %12g %12g %12g\n", its, X_est[KP], fracDelta_old[KP],
1117 delFrac[KP], fracDelta_new[KP], normUpdate, damp, funcPhaseStability);
1118 }
1119
1120 if (normUpdate < 1.0E-5 * damp) {
1121 converged = true;
1122 if (its < minNumberIterations) {
1123 converged = false;
1124 }
1125 }
1126 }
1127
1128 if (converged) {
1129 // Save the final optimized stated back into the VolPhase object for later use
1130 Vphase->setMoleFractionsState(0.0, &X_est[0], VCS_STATECALC_PHASESTABILITY);
1131
1132 // Save fracDelta for later use to initialize the problem better
1133 // @todo creationGlobalRxnNumbers needs to be calculated here and stored.
1134 Vphase->setCreationMoleNumbers(&fracDelta_new[0], creationGlobalRxnNumbers);
1135 }
1136 } else {
1137 throw CanteraError("VCS_SOLVE::vcs_phaseStabilityTest", "not done yet");
1138 }
1139 if (m_debug_print_lvl >= 2) {
1140 plogf(" ------------------------------------------------------------"
1141 "-------------------------------------------------------------\n");
1142 } else if (m_debug_print_lvl == 1) {
1143 if (funcPhaseStability > 0.0) {
1144 plogf(" --- phase %d with func = %g is to be born\n", iph, funcPhaseStability);
1145 } else {
1146 plogf(" --- phase %d with func = %g stays dead\n", iph, funcPhaseStability);
1147 }
1148 }
1149 return funcPhaseStability;
1150}
1151
1152int VCS_SOLVE::vcs_evalSS_TP(int ipr, int ip1, double Temp, double pres)
1153{
1154 for (size_t iph = 0; iph < m_numPhases; iph++) {
1155 vcs_VolPhase* vph = m_VolPhaseList[iph].get();
1158 }
1159 for (size_t k = 0; k < m_nsp; k++) {
1161 }
1162
1163 return VCS_SUCCESS;
1164}
1165
1166void VCS_SOLVE::vcs_prep(int printLvl)
1167{
1168 m_debug_print_lvl = printLvl;
1169
1170 // Calculate the Single Species status of phases
1171 // Also calculate the number of species per phase
1172 vcs_SSPhase();
1173
1174 // Set an initial estimate for the number of noncomponent species equal to
1175 // nspecies - nelements. This may be changed below
1176 m_numRxnTot = (m_nelem > m_nsp) ? 0 : m_nsp - m_nelem;
1179 for (size_t i = 0; i < m_numRxnRdc; ++i) {
1181 }
1182
1183 for (size_t kspec = 0; kspec < m_nsp; ++kspec) {
1184 double sz = 0.0;
1185 for (size_t e = 0; e < m_nelem; e++) {
1186 sz += fabs(m_formulaMatrix(kspec, e));
1187 }
1188 m_spSize[kspec] = (sz > 0.0) ? sz : 1.0;
1189 }
1190
1191 // DETERMINE THE NUMBER OF COMPONENTS
1192 //
1193 // Obtain a valid estimate of the mole fraction. This will be used as an
1194 // initial ordering vector for prioritizing which species are defined as
1195 // components.
1196 //
1197 // If a mole number estimate was supplied from the input file, use that mole
1198 // number estimate.
1199 //
1200 // If a solution estimate wasn't supplied from the input file, supply an
1201 // initial estimate for the mole fractions based on the relative reverse
1202 // ordering of the chemical potentials.
1203 //
1204 // For voltage unknowns, set these to zero for the moment.
1205 double test = -1.0e-10;
1206 bool modifiedSoln = false;
1207 if (m_doEstimateEquil < 0) {
1208 double sum = 0.0;
1209 for (size_t kspec = 0; kspec < m_nsp; ++kspec) {
1211 sum += fabs(m_molNumSpecies_old[kspec]);
1212 }
1213 }
1214 if (fabs(sum) < 1.0E-6) {
1215 modifiedSoln = true;
1216 double pres = (m_pressurePA <= 0.0) ? 1.01325E5 : m_pressurePA;
1217 vcs_evalSS_TP(0, 0, m_temperature, pres);
1218 for (size_t kspec = 0; kspec < m_nsp; ++kspec) {
1220 m_molNumSpecies_old[kspec] = - m_SSfeSpecies[kspec];
1221 } else {
1222 m_molNumSpecies_old[kspec] = 0.0;
1223 }
1224 }
1225 }
1226 test = -1.0e20;
1227 }
1228
1229 // NC = number of components is in the vcs.h common block. This call to
1230 // BASOPT doesn't calculate the stoichiometric reaction matrix.
1231 vcs_basopt(true, test);
1232
1233 if (m_nsp >= m_numComponents) {
1235 for (size_t i = 0; i < m_numRxnRdc; ++i) {
1237 }
1238 } else {
1240 }
1241
1242 // The elements might need to be rearranged.
1244
1245 // If we mucked up the solution unknowns because they were all
1246 // zero to start with, set them back to zero here
1247 if (modifiedSoln) {
1248 for (size_t kspec = 0; kspec < m_nsp; ++kspec) {
1249 m_molNumSpecies_old[kspec] = 0.0;
1250 }
1251 }
1252
1253 // Initialize various arrays in the data to zero
1254 m_feSpecies_old.assign(m_feSpecies_old.size(), 0.0);
1255 m_feSpecies_new.assign(m_feSpecies_new.size(), 0.0);
1256 m_molNumSpecies_new.assign(m_molNumSpecies_new.size(), 0.0);
1259 m_deltaPhaseMoles.assign(m_deltaPhaseMoles.size(), 0.0);
1260 m_tPhaseMoles_new.assign(m_tPhaseMoles_new.size(), 0.0);
1261
1262 // Calculate the total number of moles in all phases.
1263 vcs_tmoles();
1264
1265 // Check to see if the current problem is well posed.
1266 double sum = 0.0;
1267 for (size_t e = 0; e < m_mix->nElements(); e++) {
1268 sum += m_mix->elementMoles(e);
1269 }
1270 if (sum < 1.0E-20) {
1271 throw CanteraError("VCS_SOLVE::vcs_prep", "The problem is not well posed"
1272 " because the total number of element moles is zero.");
1273 }
1274}
1275
1277{
1278 vector<double> awSpace(m_nsp + (m_nelem + 2)*(m_nelem), 0.0);
1279 vector<double> aw(m_nelem), sa(m_nelem), ss(m_nelem);
1280 vector<double> sm(m_nelem*m_nelem);
1281
1282 size_t ncomponents = m_numComponents;
1283 if (m_debug_print_lvl >= 2) {
1284 plogf(" ");
1285 writeline('-', 77);
1286 plogf(" --- Subroutine elem_rearrange() called to ");
1287 plogf("check stoich. coefficient matrix\n");
1288 plogf(" --- and to rearrange the element ordering once\n");
1289 }
1290
1291 // Use a temporary work array for the element numbers
1292 // Also make sure the value of test is unique.
1293 bool lindep = true;
1294 double test = -1.0E10;
1295 while (lindep) {
1296 lindep = false;
1297 for (size_t i = 0; i < m_nelem; ++i) {
1298 test -= 1.0;
1299 aw[i] = m_elemAbundancesGoal[i];
1300 if (test == aw[i]) {
1301 lindep = true;
1302 }
1303 }
1304 }
1305
1306 // Top of a loop of some sort based on the index JR. JR is the current
1307 // number independent elements found.
1308 size_t jr = 0;
1309 while (jr < ncomponents) {
1310 size_t k;
1311
1312 // Top of another loop point based on finding a linearly independent
1313 // species
1314 while (true) {
1315 // Search the remaining part of the mole fraction vector, AW, for
1316 // the largest remaining species. Return its identity in K.
1317 k = m_nelem;
1318 for (size_t ielem = jr; ielem < m_nelem; ielem++) {
1319 if (m_elementActive[ielem] && aw[ielem] != test) {
1320 k = ielem;
1321 break;
1322 }
1323 }
1324 if (k == m_nelem) {
1325 throw CanteraError("VCS_SOLVE::vcs_elem_rearrange",
1326 "Shouldn't be here. Algorithm misfired.");
1327 }
1328
1329 // Assign a large negative number to the element that we have just
1330 // found, in order to take it out of further consideration.
1331 aw[k] = test;
1332
1333 // CHECK LINEAR INDEPENDENCE OF CURRENT FORMULA MATRIX LINE WITH
1334 // PREVIOUS LINES OF THE FORMULA MATRIX
1335 //
1336 // Modified Gram-Schmidt Method, p. 202 Dalquist QR factorization of
1337 // a matrix without row pivoting.
1338 size_t jl = jr;
1339
1340 // Fill in the row for the current element, k, under consideration
1341 // The row will contain the Formula matrix value for that element
1342 // from the current component.
1343 for (size_t j = 0; j < ncomponents; ++j) {
1344 sm[j + jr*ncomponents] = m_formulaMatrix(j,k);
1345 }
1346 if (jl > 0) {
1347 // Compute the coefficients of JA column of the the upper
1348 // triangular R matrix, SS(J) = R_J_JR (this is slightly
1349 // different than Dalquist) R_JA_JA = 1
1350 for (size_t j = 0; j < jl; ++j) {
1351 ss[j] = 0.0;
1352 for (size_t i = 0; i < ncomponents; ++i) {
1353 ss[j] += sm[i + jr*ncomponents] * sm[i + j*ncomponents];
1354 }
1355 ss[j] /= sa[j];
1356 }
1357
1358 // Now make the new column, (*,JR), orthogonal to the previous
1359 // columns
1360 for (size_t j = 0; j < jl; ++j) {
1361 for (size_t i = 0; i < ncomponents; ++i) {
1362 sm[i + jr*ncomponents] -= ss[j] * sm[i + j*ncomponents];
1363 }
1364 }
1365 }
1366
1367 // Find the new length of the new column in Q. It will be used in
1368 // the denominator in future row calcs.
1369 sa[jr] = 0.0;
1370 for (size_t ml = 0; ml < ncomponents; ++ml) {
1371 sa[jr] += pow(sm[ml + jr*ncomponents], 2);
1372 }
1373 // IF NORM OF NEW ROW .LT. 1E-6 REJECT
1374 if (sa[jr] > 1.0e-6) {
1375 break;
1376 }
1377 }
1378 // REARRANGE THE DATA
1379 if (jr != k) {
1380 if (m_debug_print_lvl >= 2) {
1381 plogf(" --- %-2.2s(%9.2g) replaces %-2.2s(%9.2g) as element %3d\n",
1383 m_elementName[jr], m_elemAbundancesGoal[jr], jr);
1384 }
1385 vcs_switch_elem_pos(jr, k);
1386 std::swap(aw[jr], aw[k]);
1387 }
1388
1389 // If we haven't found enough components, go back and find some more.
1390 jr++;
1391 }
1392}
1393
1394void VCS_SOLVE::vcs_switch_elem_pos(size_t ipos, size_t jpos)
1395{
1396 if (ipos == jpos) {
1397 return;
1398 }
1399 AssertThrowMsg(ipos < m_nelem && jpos < m_nelem,
1400 "vcs_switch_elem_pos",
1401 "inappropriate args: {} {}", ipos, jpos);
1402
1403 // Change the element Global Index list in each vcs_VolPhase object
1404 // to reflect the switch in the element positions.
1405 for (size_t iph = 0; iph < m_numPhases; iph++) {
1406 vcs_VolPhase* volPhase = m_VolPhaseList[iph].get();
1407 for (size_t e = 0; e < volPhase->nElemConstraints(); e++) {
1408 if (volPhase->elemGlobalIndex(e) == ipos) {
1409 volPhase->setElemGlobalIndex(e, jpos);
1410 }
1411 if (volPhase->elemGlobalIndex(e) == jpos) {
1412 volPhase->setElemGlobalIndex(e, ipos);
1413 }
1414 }
1415 }
1416 std::swap(m_elemAbundancesGoal[ipos], m_elemAbundancesGoal[jpos]);
1417 std::swap(m_elemAbundances[ipos], m_elemAbundances[jpos]);
1418 std::swap(m_elementMapIndex[ipos], m_elementMapIndex[jpos]);
1419 std::swap(m_elType[ipos], m_elType[jpos]);
1420 std::swap(m_elementActive[ipos], m_elementActive[jpos]);
1421 for (size_t j = 0; j < m_nsp; ++j) {
1422 std::swap(m_formulaMatrix(j,ipos), m_formulaMatrix(j,jpos));
1423 }
1424 std::swap(m_elementName[ipos], m_elementName[jpos]);
1425}
1426
1427double VCS_SOLVE::vcs_Hessian_diag_adj(size_t irxn, double hessianDiag_Ideal)
1428{
1429 double diag = hessianDiag_Ideal;
1430 double hessActCoef = vcs_Hessian_actCoeff_diag(irxn);
1431 if (hessianDiag_Ideal <= 0.0) {
1432 throw CanteraError("VCS_SOLVE::vcs_Hessian_diag_adj",
1433 "We shouldn't be here");
1434 }
1435 if (hessActCoef >= 0.0) {
1436 diag += hessActCoef;
1437 } else if (fabs(hessActCoef) < 0.6666 * hessianDiag_Ideal) {
1438 diag += hessActCoef;
1439 } else {
1440 diag -= 0.6666 * hessianDiag_Ideal;
1441 }
1442 return diag;
1443}
1444
1446{
1447 size_t kspec = m_indexRxnToSpecies[irxn];
1448 size_t kph = m_phaseID[kspec];
1449 double np_kspec = std::max(m_tPhaseMoles_old[kph], 1e-13);
1450 double* sc_irxn = m_stoichCoeffRxnMatrix.ptrColumn(irxn);
1451
1452 // First the diagonal term of the Jacobian
1453 double s = m_np_dLnActCoeffdMolNum(kspec,kspec) / np_kspec;
1454
1455 // Next, the other terms. Note this only a loop over the components So, it's
1456 // not too expensive to calculate.
1457 for (size_t j = 0; j < m_numComponents; j++) {
1458 if (!m_SSPhase[j]) {
1459 for (size_t k = 0; k < m_numComponents; ++k) {
1460 if (m_phaseID[k] == m_phaseID[j]) {
1461 double np = m_tPhaseMoles_old[m_phaseID[k]];
1462 if (np > 0.0) {
1463 s += sc_irxn[k] * sc_irxn[j] * m_np_dLnActCoeffdMolNum(j,k) / np;
1464 }
1465 }
1466 }
1467 if (kph == m_phaseID[j]) {
1468 s += sc_irxn[j] * (m_np_dLnActCoeffdMolNum(j,kspec) + m_np_dLnActCoeffdMolNum(kspec,j)) / np_kspec;
1469 }
1470 }
1471 }
1472 return s;
1473}
1474
1475void VCS_SOLVE::vcs_CalcLnActCoeffJac(const double* const moleSpeciesVCS)
1476{
1477 for (auto& Vphase : m_VolPhaseList) {
1478 // We don't need to call single species phases;
1479 if (!Vphase->m_singleSpecies && !Vphase->isIdealSoln()) {
1480 // update the mole numbers
1481 Vphase->setMolesFromVCS(VCS_STATECALC_OLD, moleSpeciesVCS);
1482
1483 // Download the resulting calculation into the full vector. This
1484 // scatter calculation is carried out in the vcs_VolPhase object.
1485 Vphase->sendToVCS_LnActCoeffJac(m_np_dLnActCoeffdMolNum);
1486 }
1487 }
1488}
1489
1491{
1492 // SORT DEPENDENT SPECIES IN DECREASING ORDER OF MOLES
1493 vector<pair<double, size_t>> x_order;
1494 for (size_t i = 0; i < m_nsp; i++) {
1495 x_order.push_back({-m_molNumSpecies_old[i], i});
1496 }
1497 std::sort(x_order.begin() + m_numComponents,
1498 x_order.begin() + m_numSpeciesRdc);
1499
1500 vcs_setFlagsVolPhases(false, VCS_STATECALC_OLD);
1502
1503 // PRINT OUT RESULTS
1504 plogf("\n\n\n\n");
1505 writeline('-', 80);
1506 writeline('-', 80);
1507 plogf("\t\t VCS_TP REPORT\n");
1508 writeline('-', 80);
1509 writeline('-', 80);
1510 if (iconv < 0) {
1511 plogf(" ERROR: CONVERGENCE CRITERION NOT SATISFIED.\n");
1512 } else if (iconv == 1) {
1513 plogf(" RANGE SPACE ERROR: Equilibrium Found but not all Element Abundances are Satisfied\n");
1514 }
1515
1516 // Calculate some quantities that may need updating
1517 vcs_tmoles();
1518 double totalVolume = 0.0;
1519 for (size_t iphase = 0; iphase < m_numPhases; iphase++) {
1520 vcs_VolPhase* Vphase = m_VolPhaseList[iphase].get();
1523 double Volp = Vphase->sendToVCS_VolPM(m_PMVolumeSpecies.data());
1524 totalVolume += Volp;
1525 }
1526
1527 plogf("\t\tTemperature = %15.2g Kelvin\n", m_temperature);
1528 plogf("\t\tPressure = %15.5g Pa \n", m_pressurePA);
1529 plogf("\t\ttotal Volume = %15.5g m**3\n", totalVolume);
1530
1531 // TABLE OF SPECIES IN DECREASING MOLE NUMBERS
1532 plogf("\n\n");
1533 writeline('-', 80);
1534 plogf(" Species Equilibrium kmoles ");
1535 plogf("Mole Fraction ChemPot/RT SpecUnkType\n");
1536 writeline('-', 80);
1537 for (size_t i = 0; i < m_numComponents; ++i) {
1538 plogf(" %-12.12s", m_speciesName[i]);
1539 writeline(' ', 13, false);
1540 plogf("%14.7E %14.7E %12.4E", m_molNumSpecies_old[i],
1542 plogf(" %3d", m_speciesUnknownType[i]);
1543 plogf("\n");
1544 }
1545 for (size_t i = m_numComponents; i < m_numSpeciesRdc; ++i) {
1546 size_t j = x_order[i].second;
1547 plogf(" %-12.12s", m_speciesName[j]);
1548 writeline(' ', 13, false);
1549
1551 plogf("%14.7E %14.7E %12.4E", m_molNumSpecies_old[j],
1553 plogf(" KMolNum ");
1555 plogf(" NA %14.7E %12.4E", 1.0, m_feSpecies_old[j]);
1556 plogf(" Voltage = %14.7E", m_molNumSpecies_old[j]);
1557 } else {
1558 throw CanteraError("VCS_SOLVE::vcs_report", "we have a problem");
1559 }
1560 plogf("\n");
1561 }
1562 if (m_numSpeciesRdc != m_nsp) {
1563 plogf("\n SPECIES WITH LESS THAN 1.0E-32 KMOLES:\n\n");
1564 for (size_t kspec = m_numSpeciesRdc; kspec < m_nsp; ++kspec) {
1565 plogf(" %-12.12s", m_speciesName[kspec]);
1566 // Note m_deltaGRxn_new[] stores in kspec slot not irxn slot, after solve
1567 plogf(" %14.7E %14.7E %12.4E",
1568 m_molNumSpecies_old[kspec],
1569 m_molNumSpecies_new[kspec], m_deltaGRxn_new[kspec]);
1571 plogf(" KMol_Num");
1573 plogf(" Voltage");
1574 } else {
1575 plogf(" Unknown");
1576 }
1577 plogf("\n");
1578 }
1579 }
1580 writeline('-', 80);
1581 plogf("\n");
1582
1583 // TABLE OF SPECIES FORMATION REACTIONS
1584 writeline('-', m_numComponents*10 + 45, true, true);
1585 plogf(" |ComponentID|");
1586 for (size_t j = 0; j < m_numComponents; j++) {
1587 plogf(" %3d", j);
1588 }
1589 plogf(" | |\n");
1590 plogf(" | Components|");
1591 for (size_t j = 0; j < m_numComponents; j++) {
1592 plogf(" %10.10s", m_speciesName[j]);
1593 }
1594 plogf(" | |\n");
1595 plogf(" NonComponent | Moles |");
1596 for (size_t j = 0; j < m_numComponents; j++) {
1597 plogf(" %10.3g", m_molNumSpecies_old[j]);
1598 }
1599 plogf(" | DG/RT Rxn |\n");
1600 writeline('-', m_numComponents*10 + 45);
1601 for (size_t irxn = 0; irxn < m_numRxnTot; irxn++) {
1602 size_t kspec = m_indexRxnToSpecies[irxn];
1603 plogf(" %3d ", kspec);
1604 plogf("%-10.10s", m_speciesName[kspec]);
1605 plogf("|%10.3g |", m_molNumSpecies_old[kspec]);
1606 for (size_t j = 0; j < m_numComponents; j++) {
1607 plogf(" %6.2f", m_stoichCoeffRxnMatrix(j,irxn));
1608 }
1609 plogf(" |%10.3g |", m_deltaGRxn_new[irxn]);
1610 plogf("\n");
1611 }
1612 writeline('-', m_numComponents*10 + 45);
1613 plogf("\n");
1614
1615 // TABLE OF PHASE INFORMATION
1616 vector<double> gaTPhase(m_nelem, 0.0);
1617 double totalMoles = 0.0;
1618 double gibbsPhase = 0.0;
1619 double gibbsTotal = 0.0;
1620 plogf("\n\n");
1621 plogf("\n");
1622 writeline('-', m_nelem*10 + 58);
1623 plogf(" | ElementID |");
1624 for (size_t j = 0; j < m_nelem; j++) {
1625 plogf(" %3d", j);
1626 }
1627 plogf(" | |\n");
1628 plogf(" | Element |");
1629 for (size_t j = 0; j < m_nelem; j++) {
1630 plogf(" %10.10s", m_elementName[j]);
1631 }
1632 plogf(" | |\n");
1633 plogf(" PhaseName |KMolTarget |");
1634 for (size_t j = 0; j < m_nelem; j++) {
1635 plogf(" %10.3g", m_elemAbundancesGoal[j]);
1636 }
1637 plogf(" | Gibbs Total |\n");
1638 writeline('-', m_nelem*10 + 58);
1639 for (size_t iphase = 0; iphase < m_numPhases; iphase++) {
1640 plogf(" %3d ", iphase);
1641 vcs_VolPhase* VPhase = m_VolPhaseList[iphase].get();
1642 plogf("%-12.12s |",VPhase->PhaseName);
1643 plogf("%10.3e |", m_tPhaseMoles_old[iphase]);
1644 totalMoles += m_tPhaseMoles_old[iphase];
1645 if (m_tPhaseMoles_old[iphase] != VPhase->totalMoles() &&
1646 !vcs_doubleEqual(m_tPhaseMoles_old[iphase], VPhase->totalMoles())) {
1647 throw CanteraError("VCS_SOLVE::vcs_report", "we have a problem");
1648 }
1649 // Compute the elemental abundances for each element in this phase
1650 for (size_t j = 0; j < m_nelem; ++j) {
1651 double abundance_j = 0.0;
1652 for (size_t i = 0; i < m_nsp; ++i) {
1654 && m_phaseID[i] == iphase)
1655 {
1656 abundance_j += m_formulaMatrix(i,j) * m_molNumSpecies_old[i];
1657 }
1658 }
1659 plogf(" %10.3g", abundance_j);
1660 gaTPhase[j] += abundance_j;
1661 }
1662 gibbsPhase = 0.0;
1663 for (size_t kspec = 0; kspec < m_numSpeciesRdc; ++kspec) {
1664 if (m_phaseID[kspec] == iphase && m_speciesUnknownType[kspec] != VCS_SPECIES_TYPE_INTERFACIALVOLTAGE) {
1665 gibbsPhase += m_molNumSpecies_old[kspec] * m_feSpecies_old[kspec];
1666 }
1667 }
1668 gibbsTotal += gibbsPhase;
1669 plogf(" | %18.11E |\n", gibbsPhase);
1670 }
1671 writeline('-', m_nelem*10 + 58);
1672 plogf(" TOTAL |%10.3e |", totalMoles);
1673 for (size_t j = 0; j < m_nelem; j++) {
1674 plogf(" %10.3g", gaTPhase[j]);
1675 }
1676 plogf(" | %18.11E |\n", gibbsTotal);
1677
1678 writeline('-', m_nelem*10 + 58);
1679 plogf("\n");
1680
1681 // GLOBAL SATISFACTION INFORMATION
1682
1683 // Calculate the total dimensionless Gibbs Free Energy.
1685 &m_tPhaseMoles_old[0]);
1686 plogf("\n\tTotal Dimensionless Gibbs Free Energy = G/RT = %15.7E\n", g);
1687 plogf("\nElemental Abundances (kmol): ");
1688 plogf(" Actual Target Type ElActive\n");
1689 for (size_t i = 0; i < m_nelem; ++i) {
1690 writeline(' ', 26, false);
1691 plogf("%-2.2s", m_elementName[i]);
1692 plogf("%20.12E %20.12E", m_elemAbundances[i], m_elemAbundancesGoal[i]);
1693 plogf(" %3d %3d\n", m_elType[i], m_elementActive[i]);
1694 }
1695 plogf("\n");
1696
1697 // TABLE OF SPECIES CHEM POTS
1698 writeline('-', 93, true, true);
1699 plogf("Chemical Potentials of the Species: (dimensionless)\n");
1700
1701 plogf("\t\t(RT = %g J/kmol)\n", GasConstant * m_temperature);
1702 plogf(" Name TKMoles StandStateChemPot "
1703 " ln(AC) ln(X_i) | F z_i phi | ChemPot | (-lnMnaught)");
1704 plogf("| (MolNum ChemPot)|");
1705 writeline('-', 147, true, true);
1706 for (size_t i = 0; i < m_nsp; ++i) {
1707 size_t j = x_order[i].second;
1708 size_t pid = m_phaseID[j];
1709 plogf(" %-12.12s", m_speciesName[j]);
1710 plogf(" %14.7E ", m_molNumSpecies_old[j]);
1711 plogf("%14.7E ", m_SSfeSpecies[j]);
1712 plogf("%14.7E ", log(m_actCoeffSpecies_old[j]));
1713 double tpmoles = m_tPhaseMoles_old[pid];
1714 double phi = m_phasePhi[pid];
1715 double eContrib = phi * m_chargeSpecies[j] * m_Faraday_dim;
1716 double lx = 0.0;
1718 lx = 0.0;
1719 } else {
1720 if (tpmoles > 0.0 && m_molNumSpecies_old[j] > 0.0) {
1721 double tmp = std::max(VCS_DELETE_MINORSPECIES_CUTOFF, m_molNumSpecies_old[j]);
1722 lx = log(tmp) - log(tpmoles);
1723 } else {
1724 lx = m_feSpecies_old[j] - m_SSfeSpecies[j]
1726 }
1727 }
1728 plogf("%14.7E |", lx);
1729 plogf("%14.7E | ", eContrib);
1730 double tmp = m_SSfeSpecies[j] + log(m_actCoeffSpecies_old[j])
1731 + lx - m_lnMnaughtSpecies[j] + eContrib;
1732 if (fabs(m_feSpecies_old[j] - tmp) > 1.0E-7) {
1733 throw CanteraError("VCS_SOLVE::vcs_report",
1734 "we have a problem - doesn't add up");
1735 }
1736 plogf(" %12.4E |", m_feSpecies_old[j]);
1737 if (m_lnMnaughtSpecies[j] != 0.0) {
1738 plogf("(%11.5E)", - m_lnMnaughtSpecies[j]);
1739 } else {
1740 plogf(" ");
1741 }
1742
1743 plogf("| %20.9E |", m_feSpecies_old[j] * m_molNumSpecies_old[j]);
1744 plogf("\n");
1745 }
1746 for (size_t i = 0; i < 125; i++) {
1747 plogf(" ");
1748 }
1749 plogf(" %20.9E\n", g);
1750 writeline('-', 147);
1751
1752 // TABLE OF SOLUTION COUNTERS
1753 plogf("\n");
1754 plogf("\nCounters: Iterations\n");
1755 plogf(" vcs_basopt: %5d\n", m_VCount->Basis_Opts);
1756 plogf(" vcs_TP: %5d\n", m_VCount->Its);
1757 writeline('-', 80);
1758 writeline('-', 80);
1759}
1760
1762{
1763 for (size_t j = 0; j < m_nelem; ++j) {
1764 m_elemAbundances[j] = 0.0;
1765 for (size_t i = 0; i < m_nsp; ++i) {
1768 }
1769 }
1770 }
1771}
1772
1774{
1775 size_t top = (ibound) ? m_nelem : m_numComponents;
1776 for (size_t i = 0; i < top; ++i) {
1777 // Require 12 digits of accuracy on non-zero constraints.
1778 if (m_elementActive[i] && fabs(m_elemAbundances[i] - m_elemAbundancesGoal[i]) > fabs(m_elemAbundancesGoal[i]) * 1.0e-12) {
1779 // This logic is for charge neutrality condition
1781 m_elemAbundancesGoal[i] != 0.0) {
1782 throw CanteraError("VCS_SOLVE::vcs_elabcheck",
1783 "Problem with charge neutrality condition");
1784 }
1787
1788 // Find out if the constraint is a multisign constraint. If it
1789 // is, then we have to worry about roundoff error in the
1790 // addition of terms. We are limited to 13 digits of finite
1791 // arithmetic accuracy.
1792 bool multisign = false;
1793 for (size_t kspec = 0; kspec < m_nsp; kspec++) {
1794 double eval = m_formulaMatrix(kspec,i);
1795 if (eval < 0.0) {
1796 multisign = true;
1797 }
1798 if (eval != 0.0) {
1799 scale = std::max(scale, fabs(eval * m_molNumSpecies_old[kspec]));
1800 }
1801 }
1802 if (multisign) {
1803 if (fabs(m_elemAbundances[i] - m_elemAbundancesGoal[i]) > 1e-11 * scale) {
1804 return false;
1805 }
1806 } else {
1808 return false;
1809 }
1810 }
1811 } else {
1812 return false;
1813 }
1814 }
1815 }
1816 return true;
1817}
1818
1819int VCS_SOLVE::vcs_elcorr(double aa[], double x[])
1820{
1821 int retn = 0;
1822
1823 vector<double> ga_save(m_elemAbundances);
1824 if (m_debug_print_lvl >= 2) {
1825 plogf(" --- vcsc_elcorr: Element abundances correction routine");
1826 if (m_nelem != m_numComponents) {
1827 plogf(" (m_numComponents != m_nelem)");
1828 }
1829 plogf("\n");
1830 }
1831
1832 for (size_t i = 0; i < m_nelem; ++i) {
1834 }
1835 double l2before = 0.0;
1836 for (size_t i = 0; i < m_nelem; ++i) {
1837 l2before += x[i] * x[i];
1838 }
1839 l2before = sqrt(l2before/m_nelem);
1840
1841 // Special section to take out single species, single component,
1842 // moles. These are species which have non-zero entries in the
1843 // formula matrix, and no other species have zero values either.
1844 bool changed = false;
1845 for (size_t i = 0; i < m_nelem; ++i) {
1846 int numNonZero = 0;
1847 bool multisign = false;
1848 for (size_t kspec = 0; kspec < m_nsp; kspec++) {
1850 double eval = m_formulaMatrix(kspec,i);
1851 if (eval < 0.0) {
1852 multisign = true;
1853 }
1854 if (eval != 0.0) {
1855 numNonZero++;
1856 }
1857 }
1858 }
1859 if (!multisign) {
1860 if (numNonZero < 2) {
1861 for (size_t kspec = 0; kspec < m_nsp; kspec++) {
1863 double eval = m_formulaMatrix(kspec,i);
1864 if (eval > 0.0) {
1865 m_molNumSpecies_old[kspec] = m_elemAbundancesGoal[i] / eval;
1866 changed = true;
1867 }
1868 }
1869 }
1870 } else {
1871 int numCompNonZero = 0;
1872 size_t compID = npos;
1873 for (size_t kspec = 0; kspec < m_numComponents; kspec++) {
1875 double eval = m_formulaMatrix(kspec,i);
1876 if (eval > 0.0) {
1877 compID = kspec;
1878 numCompNonZero++;
1879 }
1880 }
1881 }
1882 if (numCompNonZero == 1) {
1883 double diff = m_elemAbundancesGoal[i];
1884 for (size_t kspec = m_numComponents; kspec < m_nsp; kspec++) {
1886 double eval = m_formulaMatrix(kspec,i);
1887 diff -= eval * m_molNumSpecies_old[kspec];
1888 }
1889 m_molNumSpecies_old[compID] = std::max(0.0,diff/m_formulaMatrix(compID,i));
1890 changed = true;
1891 }
1892 }
1893 }
1894 }
1895 }
1896 if (changed) {
1897 vcs_elab();
1898 }
1899
1900 // Section to check for maximum bounds errors on all species due to
1901 // elements. This may only be tried on element types which are
1902 // VCS_ELEM_TYPE_ABSPOS. This is because no other species may have a
1903 // negative number of these.
1904 //
1905 // Note, also we can do this over ne, the number of elements, not just the
1906 // number of components.
1907 changed = false;
1908 for (size_t i = 0; i < m_nelem; ++i) {
1909 int elType = m_elType[i];
1910 if (elType == VCS_ELEM_TYPE_ABSPOS) {
1911 for (size_t kspec = 0; kspec < m_nsp; kspec++) {
1913 double atomComp = m_formulaMatrix(kspec,i);
1914 if (atomComp > 0.0) {
1915 double maxPermissible = m_elemAbundancesGoal[i] / atomComp;
1916 if (m_molNumSpecies_old[kspec] > maxPermissible) {
1917 if (m_debug_print_lvl >= 3) {
1918 plogf(" --- vcs_elcorr: Reduced species %s from %g to %g "
1919 "due to %s max bounds constraint\n",
1920 m_speciesName[kspec], m_molNumSpecies_old[kspec],
1921 maxPermissible, m_elementName[i]);
1922 }
1923 m_molNumSpecies_old[kspec] = maxPermissible;
1924 changed = true;
1926 m_molNumSpecies_old[kspec] = 0.0;
1927 if (m_SSPhase[kspec]) {
1929 } else {
1931 }
1932 if (m_debug_print_lvl >= 2) {
1933 plogf(" --- vcs_elcorr: Zeroed species %s and changed "
1934 "status to %d due to max bounds constraint\n",
1935 m_speciesName[kspec], m_speciesStatus[kspec]);
1936 }
1937 }
1938 }
1939 }
1940 }
1941 }
1942 }
1943 }
1944
1945 // Recalculate the element abundances if something has changed.
1946 if (changed) {
1947 vcs_elab();
1948 }
1949
1950 // Ok, do the general case. Linear algebra problem is of length nc, not ne,
1951 // as there may be degenerate rows when nc .ne. ne.
1953 for (size_t i = 0; i < m_numComponents; ++i) {
1955 if (fabs(x[i]) > 1.0E-13) {
1956 retn = 1;
1957 }
1958 for (size_t j = 0; j < m_numComponents; ++j) {
1959 A(j, i) = - m_formulaMatrix(i,j);
1960 }
1961 }
1962
1963 solve(A, x, 1, m_nelem);
1964
1965 // Now apply the new direction without creating negative species.
1966 double par = 0.5;
1967 for (size_t i = 0; i < m_numComponents; ++i) {
1968 if (m_molNumSpecies_old[i] > 0.0) {
1969 par = std::max(par, -x[i] / m_molNumSpecies_old[i]);
1970 }
1971 }
1972 par = std::min(par, 100.0);
1973 par = 1.0 / par;
1974 if (par < 1.0 && par > 0.0) {
1975 retn = 2;
1976 par *= 0.9999;
1977 } else {
1978 par = 1.0;
1979 }
1980 for (size_t i = 0; i < m_numComponents; ++i) {
1981 double tmp = m_molNumSpecies_old[i] + par * x[i];
1982 if (tmp > 0.0) {
1983 m_molNumSpecies_old[i] = tmp;
1984 } else if (m_SSPhase[i]) {
1985 m_molNumSpecies_old[i] = 0.0;
1986 } else {
1987 m_molNumSpecies_old[i] *= 0.0001;
1988 }
1989 }
1990
1991 // We have changed the element abundances. Calculate them again
1992 vcs_elab();
1993
1994 // We have changed the total moles in each phase. Calculate them again
1995 vcs_tmoles();
1996
1997 // cleanup/logging function to be called before returning
1998 auto finalize = [&]() {
1999 vcs_tmoles();
2000 double l2after = 0.0;
2001 for (size_t i = 0; i < m_nelem; ++i) {
2002 l2after += pow(m_elemAbundances[i] - m_elemAbundancesGoal[i], 2);
2003 }
2004 l2after = sqrt(l2after/m_nelem);
2005 if (m_debug_print_lvl >= 2) {
2006 plogf(" --- Elem_Abund: Correct Initial "
2007 " Final\n");
2008 for (size_t i = 0; i < m_nelem; ++i) {
2009 plogf(" --- ");
2010 plogf("%-2.2s", m_elementName[i]);
2011 plogf(" %20.12E %20.12E %20.12E\n", m_elemAbundancesGoal[i], ga_save[i], m_elemAbundances[i]);
2012 }
2013 plogf(" --- Diff_Norm: %20.12E %20.12E\n",
2014 l2before, l2after);
2015 }
2016 };
2017
2018 // Try some ad hoc procedures for fixing the problem
2019 if (retn >= 2) {
2020 // First find a species whose adjustment is a win-win situation.
2021 for (size_t kspec = 0; kspec < m_nsp; kspec++) {
2023 continue;
2024 }
2025 double saveDir = 0.0;
2026 bool goodSpec = true;
2027 for (size_t i = 0; i < m_numComponents; ++i) {
2028 double dir = m_formulaMatrix(kspec,i) * (m_elemAbundancesGoal[i] - m_elemAbundances[i]);
2029 if (fabs(dir) > 1.0E-10) {
2030 if (dir > 0.0) {
2031 if (saveDir < 0.0) {
2032 goodSpec = false;
2033 break;
2034 }
2035 } else if (saveDir > 0.0) {
2036 goodSpec = false;
2037 break;
2038 }
2039 saveDir = dir;
2040 } else if (m_formulaMatrix(kspec,i) != 0.0) {
2041 goodSpec = false;
2042 break;
2043 }
2044 }
2045 if (goodSpec) {
2046 int its = 0;
2047 double xx = 0.0;
2048 for (size_t i = 0; i < m_numComponents; ++i) {
2049 if (m_formulaMatrix(kspec,i) != 0.0) {
2050 xx += (m_elemAbundancesGoal[i] - m_elemAbundances[i]) / m_formulaMatrix(kspec,i);
2051 its++;
2052 }
2053 }
2054 if (its > 0) {
2055 xx /= its;
2056 }
2057 m_molNumSpecies_old[kspec] += xx;
2058 m_molNumSpecies_old[kspec] = std::max(m_molNumSpecies_old[kspec], 1.0E-10);
2059
2060 // If we are dealing with a deleted species, then we need to
2061 // reinsert it into the active list.
2062 if (kspec >= m_numSpeciesRdc) {
2063 vcs_reinsert_deleted(kspec);
2065 vcs_elab();
2066 finalize();
2067 return retn;
2068 }
2069 vcs_elab();
2070 }
2071 }
2072 }
2073 if (vcs_elabcheck(0)) {
2074 finalize();
2075 return 1;
2076 }
2077
2078 for (size_t i = 0; i < m_nelem; ++i) {
2080 (m_elType[i] == VCS_ELEM_TYPE_ABSPOS && m_elemAbundancesGoal[i] == 0.0)) {
2081 for (size_t kspec = 0; kspec < m_numSpeciesRdc; kspec++) {
2082 if (m_elemAbundances[i] > 0.0 && m_formulaMatrix(kspec,i) < 0.0) {
2084 m_molNumSpecies_old[kspec] = std::max(m_molNumSpecies_old[kspec], 0.0);
2085 vcs_elab();
2086 break;
2087 }
2088 if (m_elemAbundances[i] < 0.0 && m_formulaMatrix(kspec,i) > 0.0) {
2090 m_molNumSpecies_old[kspec] = std::max(m_molNumSpecies_old[kspec], 0.0);
2091 vcs_elab();
2092 break;
2093 }
2094 }
2095 }
2096 }
2097 if (vcs_elabcheck(1)) {
2098 finalize();
2099 return 1;
2100 }
2101
2102 // For electron charges element types, we try positive deltas in the species
2103 // concentrations to match the desired electron charge exactly.
2104 for (size_t i = 0; i < m_nelem; ++i) {
2105 double dev = m_elemAbundancesGoal[i] - m_elemAbundances[i];
2106 if (m_elType[i] == VCS_ELEM_TYPE_ELECTRONCHARGE && (fabs(dev) > 1.0E-300)) {
2107 bool useZeroed = true;
2108 for (size_t kspec = 0; kspec < m_numSpeciesRdc; kspec++) {
2109 if (dev < 0.0) {
2110 if (m_formulaMatrix(kspec,i) < 0.0 && m_molNumSpecies_old[kspec] > 0.0) {
2111 useZeroed = false;
2112 }
2113 } else {
2114 if (m_formulaMatrix(kspec,i) > 0.0 && m_molNumSpecies_old[kspec] > 0.0) {
2115 useZeroed = false;
2116 }
2117 }
2118 }
2119 for (size_t kspec = 0; kspec < m_numSpeciesRdc; kspec++) {
2120 if (m_molNumSpecies_old[kspec] > 0.0 || useZeroed) {
2121 if (dev < 0.0 && m_formulaMatrix(kspec,i) < 0.0) {
2122 double delta = dev / m_formulaMatrix(kspec,i);
2123 m_molNumSpecies_old[kspec] += delta;
2124 m_molNumSpecies_old[kspec] = std::max(m_molNumSpecies_old[kspec], 0.0);
2125 vcs_elab();
2126 break;
2127 }
2128 if (dev > 0.0 && m_formulaMatrix(kspec,i) > 0.0) {
2129 double delta = dev / m_formulaMatrix(kspec,i);
2130 m_molNumSpecies_old[kspec] += delta;
2131 m_molNumSpecies_old[kspec] = std::max(m_molNumSpecies_old[kspec], 0.0);
2132 vcs_elab();
2133 break;
2134 }
2135 }
2136 }
2137 }
2138 }
2139 if (vcs_elabcheck(1)) {
2140 finalize();
2141 return 1;
2142 }
2143 return retn;
2144}
2145
2147{
2148 double test = -1.0E-10;
2149
2150 if (m_debug_print_lvl >= 2) {
2151 plogf(" --- call setInitialMoles\n");
2152 }
2153
2154 int iter = 0;
2155 bool abundancesOK = true;
2156 vector<double> sm(m_nelem * m_nelem, 0.0);
2157 vector<double> ss(m_nelem, 0.0);
2158 vector<double> sa(m_nelem, 0.0);
2159 vector<double> wx(m_nelem, 0.0);
2160 vector<double> aw(m_nsp, 0.0);
2161
2162 for (size_t ik = 0; ik < m_nsp; ik++) {
2164 m_molNumSpecies_old[ik] = max(0.0, m_molNumSpecies_old[ik]);
2165 }
2166 }
2167
2168 if (m_debug_print_lvl >= 2) {
2170 }
2171
2172 bool redo = true;
2173 while (redo) {
2174 if (!vcs_elabcheck(0)) {
2175 if (m_debug_print_lvl >= 2) {
2176 plogf(" --- seMolesLinProg Mole numbers failing element abundances\n");
2177 plogf(" --- seMolesLinProg Call vcs_elcorr to attempt fix\n");
2178 }
2179 int retn = vcs_elcorr(&sm[0], &wx[0]);
2180 abundancesOK = (retn < 2);
2181 } else {
2182 abundancesOK = true;
2183 }
2184
2185 // Now find the optimized basis that spans the stoichiometric
2186 // coefficient matrix, based on the current composition,
2187 // m_molNumSpecies_old[] We also calculate sc[][], the reaction matrix.
2188 vcs_basopt(false, test);
2189
2190 if (m_debug_print_lvl >= 2) {
2191 plogf("iteration %d\n", iter);
2192 }
2193 redo = false;
2194 iter++;
2195 if (iter > 15) {
2196 break;
2197 }
2198
2199 // loop over all reactions
2200 for (size_t irxn = 0; irxn < m_numRxnTot; irxn++) {
2201 // dg_rt is the Delta_G / RT value for the reaction
2202 size_t ik = m_numComponents + irxn;
2203 double dg_rt = m_SSfeSpecies[ik];
2204 double dxi_min = 1.0e10;
2205 const double* sc_irxn = m_stoichCoeffRxnMatrix.ptrColumn(irxn);
2206 for (size_t jcomp = 0; jcomp < m_nelem; jcomp++) {
2207 dg_rt += m_SSfeSpecies[jcomp] * sc_irxn[jcomp];
2208 }
2209 // fwd or rev direction.
2210 // idir > 0 implies increasing the current species
2211 // idir < 0 implies decreasing the current species
2212 int idir = (dg_rt < 0.0 ? 1 : -1);
2213 if (idir < 0) {
2214 dxi_min = m_molNumSpecies_old[ik];
2215 }
2216
2217 for (size_t jcomp = 0; jcomp < m_numComponents; jcomp++) {
2218 double nu = sc_irxn[jcomp];
2219 // set max change in progress variable by
2220 // non-negativity requirement
2221 if (nu*idir < 0) {
2222 double delta_xi = fabs(m_molNumSpecies_old[jcomp]/nu);
2223 // if a component has nearly zero moles, redo
2224 // with a new set of components
2225 if (!redo && delta_xi < 1.0e-10 && (m_molNumSpecies_old[ik] >= 1.0E-10)) {
2226 if (m_debug_print_lvl >= 2) {
2227 plogf(" --- Component too small: %s\n", m_speciesName[jcomp]);
2228 }
2229 redo = true;
2230 }
2231 dxi_min = std::min(dxi_min, delta_xi);
2232 }
2233 }
2234
2235 // step the composition by dxi_min, check against zero, since
2236 // we are zeroing components and species on every step.
2237 // Redo the iteration, if a component went from positive to zero on this step.
2238 double dsLocal = idir*dxi_min;
2239 m_molNumSpecies_old[ik] += dsLocal;
2240 m_molNumSpecies_old[ik] = max(0.0, m_molNumSpecies_old[ik]);
2241 for (size_t jcomp = 0; jcomp < m_numComponents; jcomp++) {
2242 bool full = false;
2243 if (m_molNumSpecies_old[jcomp] > 1.0E-15) {
2244 full = true;
2245 }
2246 m_molNumSpecies_old[jcomp] += sc_irxn[jcomp] * dsLocal;
2247 m_molNumSpecies_old[jcomp] = max(0.0, m_molNumSpecies_old[jcomp]);
2248 if (full && m_molNumSpecies_old[jcomp] < 1.0E-60) {
2249 redo = true;
2250 }
2251 }
2252 }
2253
2254 if (m_debug_print_lvl >= 2) {
2256 }
2257 }
2258
2259 if (m_debug_print_lvl == 1) {
2261 plogf(" --- setInitialMoles end\n");
2262 }
2263 if (!abundancesOK) {
2264 return -1;
2265 } else if (iter > 15) {
2266 return 1;
2267 }
2268 return 0;
2269}
2270
2271double VCS_SOLVE::vcs_Total_Gibbs(double* molesSp, double* chemPot,
2272 double* tPhMoles)
2273{
2274 double g = 0.0;
2275
2276 for (size_t kspec = 0; kspec < m_numSpeciesRdc; ++kspec) {
2278 g += molesSp[kspec] * chemPot[kspec];
2279 }
2280 }
2281
2282 return g;
2283}
2284
2286{
2287 // Whether we have an estimate or not gets overwritten on
2288 // the call to the equilibrium solver.
2289 m_temperature = m_mix->temperature();
2290 m_pressurePA = m_mix->pressure();
2292
2293 for (size_t iphase = 0; iphase < m_numPhases; iphase++) {
2294 vcs_VolPhase* volPhase = m_VolPhaseList[iphase].get();
2295
2297
2298 // Loop through each species in the current phase
2299 size_t nSpPhase = m_mix->phase(iphase).nSpecies();
2300 if ((nSpPhase == 1) && (volPhase->phiVarIndex() == 0)) {
2302 } else if (volPhase->totalMoles() > 0.0) {
2304 } else {
2306 }
2307 }
2308
2309 // Printout the species information: PhaseID's and mole nums
2310 if (m_printLvl > 1) {
2311 writeline('=', 80, true, true);
2312 writeline('=', 20, false);
2313 plogf(" Cantera_to_vprob: START OF PROBLEM STATEMENT ");
2314 writeline('=', 20);
2315 writeline('=', 80);
2316 plogf("\n");
2317 plogf(" Phase IDs of species\n");
2318 plogf(" species phaseID phaseName ");
2319 plogf(" Initial_Estimated_kMols\n");
2320 for (size_t i = 0; i < m_nsp; i++) {
2321 size_t iphase = m_phaseID[i];
2322 vcs_VolPhase* VolPhase = m_VolPhaseList[iphase].get();
2323 plogf("%16s %5d %16s", m_speciesName[i].c_str(), iphase,
2324 VolPhase->PhaseName.c_str());
2326 plogf(" Volts = %-10.5g\n", m_molNumSpecies_old[i]);
2327 } else {
2328 plogf(" %-10.5g\n", m_molNumSpecies_old[i]);
2329 }
2330 }
2331
2332 // Printout of the Phase structure information
2333 writeline('-', 80, true, true);
2334 plogf(" Information about phases\n");
2335 plogf(" PhaseName PhaseNum SingSpec EqnState NumSpec "
2336 "Tmoles(kmol)\n");
2337
2338 for (size_t iphase = 0; iphase < m_numPhases; iphase++) {
2339 vcs_VolPhase* VolPhase = m_VolPhaseList[iphase].get();
2340 plogf("%16s %8d %8d %16s %8d ", VolPhase->PhaseName.c_str(),
2341 VolPhase->VP_ID_, VolPhase->m_singleSpecies,
2342 VolPhase->eos_name(), VolPhase->nSpecies());
2343 plogf("%16e\n", VolPhase->totalMoles());
2344 }
2345
2346 writeline('=', 80, true, true);
2347 writeline('=', 20, false);
2348 plogf(" Cantera_to_vprob: END OF PROBLEM STATEMENT ");
2349 writeline('=', 20);
2350 writeline('=', 80);
2351 plogf("\n");
2352 }
2353
2354 // m_numRxnTot = number of noncomponents, also equal to the number of
2355 // reactions. Note, it's possible that the number of elements is greater
2356 // than the number of species. In that case set the number of reactions to
2357 // zero.
2358 m_numRxnTot = (m_nelem > m_nsp) ? 0 : m_nsp - m_nelem;
2360}
2361
2363{
2364 const char* pprefix = " --- vcs_inest: ";
2365 if (m_doEstimateEquil > 0) {
2366 // Calculate the elemental abundances
2367 vcs_elab();
2368 if (vcs_elabcheck(0)) {
2369 if (m_debug_print_lvl >= 2) {
2370 plogf("%s Initial guess passed element abundances on input\n", pprefix);
2371 plogf("%s m_doEstimateEquil = 1 so will use the input mole "
2372 "numbers as estimates\n", pprefix);
2373 }
2374 return;
2375 } else if (m_debug_print_lvl >= 2) {
2376 plogf("%s Initial guess failed element abundances on input\n", pprefix);
2377 plogf("%s m_doEstimateEquil = 1 so will discard input "
2378 "mole numbers and find our own estimate\n", pprefix);
2379 }
2380 }
2381
2382 // temporary space for usage in this routine and in subroutines
2383 vector<double> sm(m_nelem*m_nelem, 0.0);
2384 vector<double> ss(m_nelem, 0.0);
2385 vector<double> sa(m_nelem, 0.0);
2386 vector<double> aw(m_nsp + m_nelem, 0.0);
2387
2388 // Go get the estimate of the solution
2389 if (m_debug_print_lvl >= 2) {
2390 plogf("%sGo find an initial estimate for the equilibrium problem\n",
2391 pprefix);
2392 }
2393 double test = -1.0E20;
2394 size_t nrxn = m_numRxnTot;
2395
2396 // CALL ROUTINE TO SOLVE MAX(CC*molNum) SUCH THAT AX*molNum = BB AND
2397 // molNum(I) .GE. 0.0. Note, both of these programs do this.
2399
2400 if (m_debug_print_lvl >= 2) {
2401 plogf("%s Mole Numbers returned from linear programming (vcs_inest initial guess):\n",
2402 pprefix);
2403 plogf("%s SPECIES MOLE_NUMBER -SS_ChemPotential\n", pprefix);
2404 for (size_t kspec = 0; kspec < m_nsp; ++kspec) {
2405 plogf("%s ", pprefix);
2406 plogf("%-12.12s", m_speciesName[kspec]);
2407 plogf(" %15.5g %12.3g\n", m_molNumSpecies_old[kspec], -m_SSfeSpecies[kspec]);
2408 }
2409 plogf("%s Element Abundance Agreement returned from linear "
2410 "programming (vcs_inest initial guess):\n", pprefix);
2411 plogf("%s Element Goal Actual\n", pprefix);
2412 for (size_t j = 0; j < m_nelem; j++) {
2413 if (m_elementActive[j]) {
2414 double tmp = 0.0;
2415 for (size_t kspec = 0; kspec < m_nsp; ++kspec) {
2416 tmp += m_formulaMatrix(kspec,j) * m_molNumSpecies_old[kspec];
2417 }
2418 plogf("%s ", pprefix);
2419 plogf(" %-9.9s", m_elementName[j]);
2420 plogf(" %12.3g %12.3g\n", m_elemAbundancesGoal[j], tmp);
2421 }
2422 }
2423 writelogendl();
2424 }
2425
2426 // Make sure all species have positive definite mole numbers Set voltages to
2427 // zero for now, until we figure out what to do
2428 m_deltaMolNumSpecies.assign(m_deltaMolNumSpecies.size(), 0.0);
2429 for (size_t kspec = 0; kspec < m_nsp; ++kspec) {
2431 if (m_molNumSpecies_old[kspec] <= 0.0) {
2432 // HKM Should eventually include logic here for non SS phases
2433 if (!m_SSPhase[kspec]) {
2434 m_molNumSpecies_old[kspec] = 1.0e-30;
2435 }
2436 }
2437 } else {
2438 m_molNumSpecies_old[kspec] = 0.0;
2439 }
2440 }
2441
2442 // Now find the optimized basis that spans the stoichiometric coefficient matrix
2443 vcs_basopt(false, test);
2444
2445 // CALCULATE TOTAL MOLES, CHEMICAL POTENTIALS OF BASIS
2446
2447 // Calculate TMoles and m_tPhaseMoles_old[]
2448 vcs_tmoles();
2449
2450 // m_tPhaseMoles_new[] will consist of just the component moles
2451 m_tPhaseMoles_new.assign(m_numPhases, 1.0e-20);
2452 for (size_t kspec = 0; kspec < m_numComponents; ++kspec) {
2455 }
2456 }
2457 double TMolesMultiphase = 0.0;
2458 for (size_t iph = 0; iph < m_numPhases; iph++) {
2459 if (! m_VolPhaseList[iph]->m_singleSpecies) {
2460 TMolesMultiphase += m_tPhaseMoles_new[iph];
2461 }
2462 }
2464 for (size_t kspec = 0; kspec < m_numComponents; ++kspec) {
2466 m_molNumSpecies_new[kspec] = 0.0;
2467 }
2468 }
2470
2471 for (size_t kspec = 0; kspec < m_numComponents; ++kspec) {
2473 if (! m_SSPhase[kspec]) {
2474 size_t iph = m_phaseID[kspec];
2475 m_feSpecies_new[kspec] += log(m_molNumSpecies_new[kspec] / m_tPhaseMoles_old[iph]);
2476 }
2477 } else {
2478 m_molNumSpecies_new[kspec] = 0.0;
2479 }
2480 }
2481 vcs_deltag(0, true, VCS_STATECALC_NEW);
2482 if (m_debug_print_lvl >= 2) {
2483 for (size_t kspec = 0; kspec < m_nsp; ++kspec) {
2484 plogf("%s", pprefix);
2485 plogf("%-12.12s", m_speciesName[kspec]);
2486 if (kspec < m_numComponents) {
2487 plogf("fe* = %15.5g ff = %15.5g\n", m_feSpecies_new[kspec],
2488 m_SSfeSpecies[kspec]);
2489 } else {
2490 plogf("fe* = %15.5g ff = %15.5g dg* = %15.5g\n",
2492 }
2493 }
2494 }
2495
2496 // ESTIMATE REACTION ADJUSTMENTS
2497 vector<double>& xtphMax = m_TmpPhase;
2498 vector<double>& xtphMin = m_TmpPhase2;
2499 m_deltaPhaseMoles.assign(m_deltaPhaseMoles.size(), 0.0);
2500 for (size_t iph = 0; iph < m_numPhases; iph++) {
2501 xtphMax[iph] = log(m_tPhaseMoles_new[iph] * 1.0E32);
2502 xtphMin[iph] = log(m_tPhaseMoles_new[iph] * 1.0E-32);
2503 }
2504 for (size_t irxn = 0; irxn < nrxn; ++irxn) {
2505 size_t kspec = m_indexRxnToSpecies[irxn];
2506
2507 // For single species phases, we will not estimate the mole numbers. If
2508 // the phase exists, it stays. If it doesn't exist in the estimate, it
2509 // doesn't come into existence here.
2510 if (! m_SSPhase[kspec]) {
2511 size_t iph = m_phaseID[kspec];
2512 if (m_deltaGRxn_new[irxn] > xtphMax[iph]) {
2513 m_deltaGRxn_new[irxn] = 0.8 * xtphMax[iph];
2514 }
2515 if (m_deltaGRxn_new[irxn] < xtphMin[iph]) {
2516 m_deltaGRxn_new[irxn] = 0.8 * xtphMin[iph];
2517 }
2518
2519 // HKM -> The TMolesMultiphase is a change of mine. It more evenly
2520 // distributes the initial moles amongst multiple multispecies
2521 // phases according to the relative values of the standard state
2522 // free energies. There is no change for problems with one
2523 // multispecies phase. It cut diamond4.vin iterations down from 62
2524 // to 14.
2525 m_deltaMolNumSpecies[kspec] = 0.5 * (m_tPhaseMoles_new[iph] + TMolesMultiphase)
2526 * exp(-m_deltaGRxn_new[irxn]);
2527
2528 for (size_t k = 0; k < m_numComponents; ++k) {
2530 }
2531
2532 for (iph = 0; iph < m_numPhases; iph++) {
2534 }
2535 }
2536 }
2537 if (m_debug_print_lvl >= 2) {
2538 for (size_t kspec = 0; kspec < m_nsp; ++kspec) {
2540 plogf("%sdirection (", pprefix);
2541 plogf("%-12.12s", m_speciesName[kspec]);
2542 plogf(") = %g", m_deltaMolNumSpecies[kspec]);
2543 if (m_SSPhase[kspec]) {
2544 if (m_molNumSpecies_old[kspec] > 0.0) {
2545 plogf(" (ssPhase exists at w = %g moles)", m_molNumSpecies_old[kspec]);
2546 } else {
2547 plogf(" (ssPhase doesn't exist -> stability not checked)");
2548 }
2549 }
2550 writelogendl();
2551 }
2552 }
2553 }
2554
2555 // KEEP COMPONENT SPECIES POSITIVE
2556 double par = 0.5;
2557 for (size_t kspec = 0; kspec < m_numComponents; ++kspec) {
2559 par < -m_deltaMolNumSpecies[kspec] / m_molNumSpecies_new[kspec]) {
2560 par = -m_deltaMolNumSpecies[kspec] / m_molNumSpecies_new[kspec];
2561 }
2562 }
2563 par = 1. / par;
2564 if (par <= 1.0 && par > 0.0) {
2565 par *= 0.8;
2566 } else {
2567 par = 1.0;
2568 }
2569
2570 // CALCULATE NEW MOLE NUMBERS
2571 size_t lt = 0;
2572 size_t ikl = 0;
2573 double s1 = 0.0;
2574 while (true) {
2575 for (size_t kspec = 0; kspec < m_numComponents; ++kspec) {
2577 m_molNumSpecies_old[kspec] = m_molNumSpecies_new[kspec] + par * m_deltaMolNumSpecies[kspec];
2578 } else {
2579 m_deltaMolNumSpecies[kspec] = 0.0;
2580 }
2581 }
2582 for (size_t kspec = m_numComponents; kspec < m_nsp; ++kspec) {
2584 m_deltaMolNumSpecies[kspec] != 0.0) {
2585 m_molNumSpecies_old[kspec] = m_deltaMolNumSpecies[kspec] * par;
2586 }
2587 }
2588
2589 // We have a new w[] estimate, go get the TMoles and m_tPhaseMoles_old[]
2590 // values
2591 vcs_tmoles();
2592 if (lt > 0) {
2593 break;
2594 }
2595
2596 // CONVERGENCE FORCING SECTION
2597 vcs_setFlagsVolPhases(false, VCS_STATECALC_OLD);
2599 double s = 0.0;
2600 for (size_t kspec = 0; kspec < m_nsp; ++kspec) {
2601 s += m_deltaMolNumSpecies[kspec] * m_feSpecies_old[kspec];
2602 }
2603 if (s == 0.0) {
2604 break;
2605 }
2606 if (s < 0.0 && ikl == 0) {
2607 break;
2608 }
2609
2610 // TRY HALF STEP SIZE
2611 if (ikl == 0) {
2612 s1 = s;
2613 par *= 0.5;
2614 ikl = 1;
2615 continue;
2616 }
2617
2618 // FIT PARABOLA THROUGH HALF AND FULL STEPS
2619 double xl = (1.0 - s / (s1 - s)) * 0.5;
2620 if (xl < 0.0) {
2621 // POOR DIRECTION, REDUCE STEP SIZE TO 0.2
2622 par *= 0.2;
2623 } else {
2624 if (xl > 1.0) {
2625 // TOO BIG A STEP, TAKE ORIGINAL FULL STEP
2626 par *= 2.0;
2627 } else {
2628 // ACCEPT RESULTS OF FORCER
2629 par = par * 2.0 * xl;
2630 }
2631 }
2632 lt = 1;
2633 }
2634
2635 if (m_debug_print_lvl >= 2) {
2636 plogf("%s Final Mole Numbers produced by inest:\n",
2637 pprefix);
2638 plogf("%s SPECIES MOLE_NUMBER\n", pprefix);
2639 for (size_t kspec = 0; kspec < m_nsp; ++kspec) {
2640 plogf("%s %-12.12s %g\n",
2641 pprefix, m_speciesName[kspec], m_molNumSpecies_old[kspec]);
2642 }
2643 }
2644
2645 // Calculate the elemental abundances
2646 vcs_elab();
2647
2648 // If we still fail to achieve the correct elemental abundances, try to fix
2649 // the problem again by calling the main elemental abundances fixer routine,
2650 // used in the main program. This attempts to tweak the mole numbers of the
2651 // component species to satisfy the element abundance constraints.
2652 //
2653 // Note: We won't do this unless we have to since it involves inverting a
2654 // matrix.
2655 bool rangeCheck = vcs_elabcheck(1);
2656 if (!vcs_elabcheck(0)) {
2657 if (m_debug_print_lvl >= 2) {
2658 plogf("%sInitial guess failed element abundances\n", pprefix);
2659 plogf("%sCall vcs_elcorr to attempt fix\n", pprefix);
2660 }
2661 vcs_elcorr(&sm[0], &aw[0]);
2662 rangeCheck = vcs_elabcheck(1);
2663 if (!vcs_elabcheck(0)) {
2664 throw CanteraError("VCS_SOLVE::vcs_inest",
2665 "Initial guess still fails element abundance equations\n"
2666 "Inability to ever satisfy element abundance constraints is probable");
2667 } else {
2668 if (m_debug_print_lvl >= 2) {
2669 if (rangeCheck) {
2670 plogf("%sInitial guess now satisfies element abundances\n", pprefix);
2671 } else {
2672 plogf("%sElement Abundances RANGE ERROR\n", pprefix);
2673 plogf("%s - Initial guess satisfies NC=%d element abundances, "
2674 "BUT not NE=%d element abundances\n", pprefix,
2676 }
2677 }
2678 }
2679 } else {
2680 if (m_debug_print_lvl >= 2) {
2681 if (rangeCheck) {
2682 plogf("%sInitial guess satisfies element abundances\n", pprefix);
2683 } else {
2684 plogf("%sElement Abundances RANGE ERROR\n", pprefix);
2685 plogf("%s - Initial guess satisfies NC=%d element abundances, "
2686 "BUT not NE=%d element abundances\n", pprefix,
2688 }
2689 }
2690 }
2691
2692 if (m_debug_print_lvl >= 2) {
2693 plogf("%sTotal Dimensionless Gibbs Free Energy = %15.7E\n", pprefix,
2695 &m_tPhaseMoles_old[0]));
2696 }
2697
2699}
2700
2702{
2703 vector<int> numPhSpecies(m_numPhases, 0);
2704 for (size_t kspec = 0; kspec < m_nsp; ++kspec) {
2705 numPhSpecies[m_phaseID[kspec]]++;
2706 }
2707
2708 // Handle the special case of a single species in a phase that has been
2709 // earmarked as a multispecies phase. Treat that species as a single-species
2710 // phase
2711 for (size_t iph = 0; iph < m_numPhases; iph++) {
2712 m_VolPhaseList[iph]->m_singleSpecies = (numPhSpecies[iph] <= 1);
2713 }
2714
2715 // Fill in some useful arrays here that have to do with the static
2716 // information concerning the phase ID of species. SSPhase = Boolean
2717 // indicating whether a species is in a single species phase or not.
2718 for (size_t kspec = 0; kspec < m_nsp; kspec++) {
2719 m_SSPhase[kspec] = m_VolPhaseList[m_phaseID[kspec]]->m_singleSpecies;
2720 }
2721}
2722
2724{
2725 m_VCount->Its = 0;
2726 m_VCount->Basis_Opts = 0;
2727 if (ifunc) {
2728 m_VCount->T_Its = 0;
2732 }
2733}
2734
2736{
2737 plogf("\nTCounters: Num_Calls Total_Its\n");
2738 plogf(" vcs_basopt: %5d %5d\n",
2740 plogf(" vcs_TP: %5d %5d\n",
2742 plogf(" vcs_inest: %5d\n", m_VCount->T_Calls_Inest);
2743}
2744
2745void VCS_SOLVE::prob_report(int print_lvl)
2746{
2747 m_printLvl = print_lvl;
2748
2749 // Printout the species information: PhaseID's and mole nums
2750 if (m_printLvl > 0) {
2751 writeline('=', 80, true, true);
2752 writeline('=', 20, false);
2753 plogf(" VCS_PROB: PROBLEM STATEMENT ");
2754 writeline('=', 31);
2755 writeline('=', 80);
2756 plogf("\n");
2757 plogf("\tSolve a constant T, P problem:\n");
2758 plogf("\t\tT = %g K\n", m_temperature);
2759 double pres_atm = m_pressurePA / 1.01325E5;
2760
2761 plogf("\t\tPres = %g atm\n", pres_atm);
2762 plogf("\n");
2763 plogf(" Phase IDs of species\n");
2764 plogf(" species phaseID phaseName ");
2765 plogf(" Initial_Estimated_Moles Species_Type\n");
2766 for (size_t i = 0; i < m_nsp; i++) {
2767 vcs_VolPhase* Vphase = m_VolPhaseList[m_phaseID[i]].get();
2768 plogf("%16s %5d %16s", m_mix->speciesName(i), m_phaseID[i],
2769 Vphase->PhaseName);
2770 if (m_doEstimateEquil >= 0) {
2771 plogf(" %-10.5g", m_molNumSpecies_old[i]);
2772 } else {
2773 plogf(" N/A");
2774 }
2776 plogf(" Mol_Num");
2778 plogf(" Voltage");
2779 } else {
2780 plogf(" ");
2781 }
2782 plogf("\n");
2783 }
2784
2785 // Printout of the Phase structure information
2786 writeline('-', 80, true, true);
2787 plogf(" Information about phases\n");
2788 plogf(" PhaseName PhaseNum SingSpec EqnState NumSpec "
2789 "Tmoles(kmol)\n");
2790
2791 for (size_t iphase = 0; iphase < m_numPhases; iphase++) {
2792 vcs_VolPhase* Vphase = m_VolPhaseList[iphase].get();
2793 plogf("%16s %8d %8d ", Vphase->PhaseName,
2794 Vphase->VP_ID_, Vphase->m_singleSpecies);
2795 plogf("%16s %8d ", Vphase->eos_name(),
2796 Vphase->nSpecies());
2797 if (m_doEstimateEquil >= 0) {
2798 plogf("%16e\n", Vphase->totalMoles());
2799 } else {
2800 plogf(" N/A\n");
2801 }
2802 }
2803
2804 plogf("\nElemental Abundances: ");
2805 plogf(" Target_kmol ElemType ElActive\n");
2806 for (size_t i = 0; i < m_nelem; ++i) {
2807 writeline(' ', 26, false);
2808 plogf("%-2.2s", m_elementName[i]);
2809 plogf("%20.12E ", m_elemAbundancesGoal[i]);
2810 plogf("%3d %3d\n", m_elType[i], m_elementActive[i]);
2811 }
2812
2813 plogf("\nChemical Potentials: (J/kmol)\n");
2814 plogf(" Species (phase) "
2815 " SS0ChemPot StarChemPot\n");
2816 for (size_t iphase = 0; iphase < m_numPhases; iphase++) {
2817 vcs_VolPhase* Vphase = m_VolPhaseList[iphase].get();
2819 for (size_t kindex = 0; kindex < Vphase->nSpecies(); kindex++) {
2820 size_t kglob = Vphase->spGlobalIndexVCS(kindex);
2821 plogf("%16s ", m_mix->speciesName(kglob));
2822 if (kindex == 0) {
2823 plogf("%16s", Vphase->PhaseName);
2824 } else {
2825 plogf(" ");
2826 }
2827
2828 plogf("%16g %16g\n", Vphase->G0_calc_one(kindex),
2829 Vphase->GStar_calc_one(kindex));
2830 }
2831 }
2832 writeline('=', 80, true, true);
2833 writeline('=', 20, false);
2834 plogf(" VCS_PROB: END OF PROBLEM STATEMENT ");
2835 writeline('=', 24);
2836 writeline('=', 80);
2837 plogf("\n");
2838 }
2839}
2840
2841size_t VCS_SOLVE::addOnePhaseSpecies(vcs_VolPhase* volPhase, size_t k, size_t kT)
2842{
2843 if (kT > m_nsp) {
2844 // Need to expand the number of species here
2845 throw CanteraError("VCS_SOLVE::addOnePhaseSpecies", "Shouldn't be here");
2846 }
2847 const Array2D& fm = volPhase->getFormulaMatrix();
2848 for (size_t eVP = 0; eVP < volPhase->nElemConstraints(); eVP++) {
2849 size_t e = volPhase->elemGlobalIndex(eVP);
2850 AssertThrowMsg(e != npos, "VCS_PROB::addOnePhaseSpecies",
2851 "element not found");
2852 m_formulaMatrix(kT,e) = fm(k,eVP);
2853 }
2854
2855 // Tell the phase object about the current position of the species within
2856 // the global species vector
2857 volPhase->setSpGlobalIndexVCS(k, kT);
2858 return kT;
2859}
2860
2861size_t VCS_SOLVE::addElement(const char* elNameNew, int elType, int elactive)
2862{
2863 if (!elNameNew) {
2864 throw CanteraError("VCS_SOLVE::addElement",
2865 "error: element must have a name");
2866 }
2867 m_nelem++;
2869
2872 m_elType.push_back(elType);
2873 m_elementActive.push_back(elactive);
2874 m_elemAbundances.push_back(0.0);
2875 m_elemAbundancesGoal.push_back(0.0);
2876 m_elementMapIndex.push_back(0);
2877 m_elementName.push_back(elNameNew);
2878 return m_nelem - 1;
2879}
2880
2881size_t VCS_SOLVE::elementIndex(const string& elementName) const
2882{
2883 for (size_t e = 0; e < m_nelem; e++) {
2884 if (m_elementName[e] == elementName) {
2885 return e;
2886 }
2887 }
2888 return npos;
2889}
2890
2891}
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:127
size_t nColumns() const
Number of columns.
Definition Array.h:181
double * ptrColumn(size_t j)
Return a pointer to the top of column j, columns are contiguous in memory.
Definition Array.h:203
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:47
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:409
double temperature() const
Temperature [K].
Definition MultiPhase.h:349
string speciesName(size_t kGlob) const
Name of species with global index kGlob.
size_t nElements() const
Number of elements.
Definition MultiPhase.h:113
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:245
size_t nDim() const
Returns the number of spatial dimensions (1, 2, or 3)
Definition Phase.h:582
double molecularWeight(size_t k) const
Molecular weight of species k.
Definition Phase.cpp:408
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:574
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.
size_t vcs_popPhaseID(vector< size_t > &phasePopPhaseIDs)
Decision as to whether a phase pops back into existence.
vector< double > m_deltaPhaseMoles
Change in the total moles in each phase.
Definition vcs_solve.h:1074
size_t m_numPhases
Number of Phases in the problem.
Definition vcs_solve.h:876
int m_doEstimateEquil
Setting for whether to do an initial estimate.
Definition vcs_solve.h:953
vector< double > m_deltaMolNumSpecies
Reaction Adjustments for each species during the current step.
Definition vcs_solve.h:1021
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:919
vector< double > m_scSize
Absolute size of the stoichiometric coefficients.
Definition vcs_solve.h:912
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:798
vector< double > m_wtSpecies
Molecular weight of each species.
Definition vcs_solve.h:1245
vector< size_t > m_indexRxnToSpecies
Mapping between the species index for noncomponent species and the full species index.
Definition vcs_solve.h:1146
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:800
vector< size_t > m_phaseID
Mapping from the species number to the phase number.
Definition vcs_solve.h:1157
vector< string > m_speciesName
Species string name for the kth species.
Definition vcs_solve.h:1164
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:1089
vector< double > m_deltaGRxn_new
Delta G(irxn) for the noncomponent species in the mechanism.
Definition vcs_solve.h:1001
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:1237
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:855
vector< double > m_PMVolumeSpecies
Partial molar volumes of the species.
Definition vcs_solve.h:1264
vector< int > m_elementActive
Specifies whether an element constraint is active.
Definition vcs_solve.h:1188
int m_useActCoeffJac
Choice of Hessians.
Definition vcs_solve.h:1257
vector< size_t > m_speciesLocalPhaseIndex
Index that keeps track of the index of the species within the local phase.
Definition vcs_solve.h:1120
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:1170
vector< int > m_speciesUnknownType
Specifies the species unknown type.
Definition vcs_solve.h:973
int vcs_elcorr(double aa[], double x[])
This subroutine corrects for element abundances.
vector< double > m_TmpPhase
Temporary vector of length NPhase.
Definition vcs_solve.h:1065
vector< size_t > m_speciesMapIndex
Index vector that keeps track of the species vector rearrangement.
Definition vcs_solve.h:1106
Array2D m_phaseParticipation
This is 1 if the phase, iphase, participates in the formation reaction irxn, and zero otherwise.
Definition vcs_solve.h:984
vector< int > m_actConventionSpecies
specifies the activity convention of the phase containing the species
Definition vcs_solve.h:1200
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:1030
vector< double > m_tPhaseMoles_new
total kmols of species in each phase in the tentative soln vector
Definition vcs_solve.h:1062
vector< double > m_actCoeffSpecies_old
Molar-based Activity Coefficients for Species based on old mole numbers.
Definition vcs_solve.h:1227
vector< double > m_SSfeSpecies
Standard state chemical potentials for species K at the current temperature and pressure.
Definition vcs_solve.h:927
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:1216
vector< int > m_elType
Type of the element constraint.
Definition vcs_solve.h:1182
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:960
double vcs_Total_Gibbs(double *w, double *fe, double *tPhMoles)
Calculate the total dimensionless Gibbs free energy.
void vcs_CalcLnActCoeffJac(const double *const moleSpeciesVCS)
Recalculate all of the activity coefficients in all of the phases based on input mole numbers.
vector< size_t > m_elementMapIndex
Index vector that keeps track of the rearrangement of the elements.
Definition vcs_solve.h:1133
vector< double > m_aw
work array of mole fractions used in vcs_basopt()
Definition vcs_solve.h:801
size_t m_numSpeciesRdc
Current number of species in the problems.
Definition vcs_solve.h:865
vector< double > m_feSpecies_old
Free energy vector from the start of the current iteration.
Definition vcs_solve.h:934
int m_debug_print_lvl
Debug printing lvl.
Definition vcs_solve.h:1284
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:799
vector< double > m_chargeSpecies
Charge of each species. Length = number of species.
Definition vcs_solve.h:1248
Array2D m_formulaMatrix
Formula matrix for the problem.
Definition vcs_solve.h:884
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:903
size_t m_numRxnRdc
Current number of non-component species in the problem.
Definition vcs_solve.h:869
VCS_COUNTERS * m_VCount
Timing and iteration counters for the vcs object.
Definition vcs_solve.h:1270
vector< double > m_deltaGRxn_Deficient
Last deltag[irxn] from the previous step with additions for possible births of zeroed phases.
Definition vcs_solve.h:1008
vector< unique_ptr< vcs_VolPhase > > m_VolPhaseList
Array of Phase Structures. Length = number of phases.
Definition vcs_solve.h:1191
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:858
void vcs_prob_specifyFully()
Fully specify the problem to be solved.
double m_pressurePA
Pressure.
Definition vcs_solve.h:1080
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:943
double m_Faraday_dim
dimensionless value of Faraday's constant, F / RT (1/volt)
Definition vcs_solve.h:1267
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:852
vector< double > m_actCoeffSpecies_new
Molar-based Activity Coefficients for Species.
Definition vcs_solve.h:1220
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:1161
vector< double > m_phasePhi
electric potential of the iph phase
Definition vcs_solve.h:987
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:1039
size_t vcs_RxnStepSizes(int &forceComponentCalc, size_t &kSpecial)
Calculates formation reaction step sizes.
Array2D m_deltaMolNumPhase
Change in the number of moles of phase, iphase, due to the noncomponent formation reaction,...
Definition vcs_solve.h:980
vector< int > m_phaseActConvention
specifies the activity convention of the phase.
Definition vcs_solve.h:1209
vector< int > m_speciesStatus
Major -Minor status vector for the species in the problem.
Definition vcs_solve.h:1154
vector< double > m_deltaGRxn_tmp
Temporary vector of Rxn DeltaG's.
Definition vcs_solve.h:1014
vector< double > m_TmpPhase2
Temporary vector of length NPhase.
Definition vcs_solve.h:1068
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:806
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:1053
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:861
vector< double > m_molNumSpecies_new
Tentative value of the mole number vector.
Definition vcs_solve.h:991
double m_temperature
Temperature (Kelvin)
Definition vcs_solve.h:1077
vector< double > m_deltaGRxn_old
Last deltag[irxn] from the previous step.
Definition vcs_solve.h:1004
void vcs_TCounters_report()
Create a report on the plog file containing iteration counters.
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:1045
Phase information and Phase calculations for vcs.
void setCreationMoleNumbers(const double *const n_k, const vector< size_t > &creationGlobalRxnNumbers)
Sets the creationMoleNum's within the phase object.
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.
size_t nSpecies() const
Return the number of species in the phase.
string PhaseName
String name for the phase.
void setMoleFractionsState(const double molNum, const double *const moleFracVec, const int vcsStateStatus)
Set the moles and/or mole fractions within 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.
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 setMolesFromVCS(const int stateCalc, const double *molesSpeciesVCS=0)
Set the moles within the phase.
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 sendToVCS_GStar(double *const gstar) const
Fill in the standard state Gibbs free energy vector for VCS.
void sendToVCS_ActCoeff(const int stateCalc, double *const AC)
Fill in an activity coefficients vector within a VCS_SOLVE object.
double sendToVCS_VolPM(double *const VolPM) const
Fill in the partial molar volume vector for VCS.
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 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.
const vector< double > & creationMoleNumbers(vector< size_t > &creationGlobalRxnNumbers) const
Return a const reference to the creationMoleNumbers stored 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:104
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:131
const double GasConstant
Universal Gas Constant [J/kmol/K].
Definition ct_defs.h:120
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:180
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
void solve(DenseMatrix &A, double *b, size_t nrhs, size_t ldb)
Solve Ax = b. Array b is overwritten on exit with x.
double vcs_l2norm(const vector< double > &vec)
determine the l2 norm of a vector of doubles
Definition vcs_util.cpp:19
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...