Cantera  3.1.0a1
Loading...
Searching...
No Matches
MultiPhaseEquil.cpp
Go to the documentation of this file.
1/**
2 * @file MultiPhaseEquil.cpp
3 */
4
5// This file is part of Cantera. See License.txt in the top-level directory or
6// at https://cantera.org/license.txt for license and copyright information.
7
11
12#include <cstdio>
13#include <numeric>
14
15namespace Cantera
16{
17
18MultiPhaseEquil::MultiPhaseEquil(MultiPhase* mix, bool start, int loglevel) : m_mix(mix)
19{
20 // store some mixture parameters locally
21 m_nel_mix = mix->nElements();
22 m_nsp_mix = mix->nSpecies();
23 m_press = mix->pressure();
24 m_temp = mix->temperature();
25
26 m_incl_species.resize(m_nsp_mix,1);
27 m_incl_element.resize(m_nel_mix,1);
28 for (size_t m = 0; m < m_nel_mix; m++) {
29 string enm = mix->elementName(m);
30 // element 'E' or 'e' represents an electron; this requires special
31 // handling, so save its index for later use
32 if (enm == "E" || enm == "e") {
33 m_eloc = m;
34 }
35 // if an element other than electrons is not present in the mixture,
36 // then exclude it and all species containing it from the calculation.
37 // Electrons are a special case, since a species can have a negative
38 // number of 'atoms' of electrons (positive ions).
39 if (m_mix->elementMoles(m) <= 0.0 && m != m_eloc) {
40 m_incl_element[m] = 0;
41 for (size_t k = 0; k < m_nsp_mix; k++) {
42 if (m_mix->nAtoms(k,m) != 0.0) {
43 m_incl_species[k] = 0;
44 }
45 }
46 }
47 }
48
49 // Now build the list of elements to be included, starting with
50 // electrons, if they are present.
51 if (m_eloc < m_nel_mix) {
52 m_element.push_back(m_eloc);
53 m_nel++;
54 }
55 // add the included elements other than electrons
56 for (size_t m = 0; m < m_nel_mix; m++) {
57 if (m_incl_element[m] == 1 && m != m_eloc) {
58 m_nel++;
59 m_element.push_back(m);
60 }
61 }
62
63 // include pure single-constituent phases only if their thermo data are
64 // valid for this temperature. This is necessary, since some thermo
65 // polynomial fits are done only for a limited temperature range. For
66 // example, using the NASA polynomial fits for solid ice and liquid water,
67 // if this were not done the calculation would predict solid ice to be
68 // present far above its melting point, since the thermo polynomial fits
69 // only extend to 273.15 K, and give unphysical results above this
70 // temperature, leading (incorrectly) to Gibbs free energies at high
71 // temperature lower than for liquid water.
72 for (size_t k = 0; k < m_nsp_mix; k++) {
73 size_t ip = m_mix->speciesPhaseIndex(k);
74 if (!m_mix->solutionSpecies(k) &&
75 !m_mix->tempOK(ip)) {
76 m_incl_species[k] = 0;
77 if (m_mix->speciesMoles(k) > 0.0) {
78 throw CanteraError("MultiPhaseEquil::MultiPhaseEquil",
79 "condensed-phase species"+ m_mix->speciesName(k)
80 + " is excluded since its thermo properties are \n"
81 "not valid at this temperature, but it has "
82 "non-zero moles in the initial state.");
83 }
84 }
85 }
86
87 // Now build the list of all species to be included in the calculation.
88 for (size_t k = 0; k < m_nsp_mix; k++) {
89 if (m_incl_species[k] ==1) {
90 m_nsp++;
91 m_species.push_back(k);
92 }
93 }
94
95 // some work arrays for internal use
96 m_work.resize(m_nsp);
97 m_work2.resize(m_nsp);
98 m_work3.resize(m_nsp_mix);
99 m_mu.resize(m_nsp_mix);
100
101 // number of moles of each species
102 m_moles.resize(m_nsp);
103 m_lastmoles.resize(m_nsp);
104 m_dxi.resize(nFree());
105
106 // initialize the mole numbers to the mixture composition
107 for (size_t ik = 0; ik < m_nsp; ik++) {
108 m_moles[ik] = m_mix->speciesMoles(m_species[ik]);
109 }
110
111 // Delta G / RT for each reaction
112 m_deltaG_RT.resize(nFree(), 0.0);
113 m_majorsp.resize(m_nsp);
114 m_sortindex.resize(m_nsp,0);
115 m_lastsort.resize(m_nel);
116 m_solnrxn.resize(nFree());
117 m_A.resize(m_nel, m_nsp, 0.0);
118 m_N.resize(m_nsp, nFree());
119 m_order.resize(std::max(m_nsp, m_nel), 0);
120 iota(m_order.begin(), m_order.begin() + m_nsp, 0);
121
122 // if the 'start' flag is set, estimate the initial mole numbers by doing a
123 // linear Gibbs minimization. In this case, only the elemental composition
124 // of the initial mixture state matters.
125 if (start) {
126 setInitialMoles(loglevel-1);
127 }
128 computeN();
129
130 // Take a very small step in composition space, so that no
131 // species has precisely zero moles.
132 vector<double> dxi(nFree(), 1.0e-20);
133 if (!dxi.empty()) {
134 multiply(m_N, dxi.data(), m_work.data());
135 unsort(m_work);
136 }
137
138 for (size_t k = 0; k < m_nsp; k++) {
139 m_moles[k] += m_work[k];
140 m_lastmoles[k] = m_moles[k];
141 if (m_mix->solutionSpecies(m_species[k])) {
142 m_dsoln.push_back(1);
143 } else {
144 m_dsoln.push_back(0);
145 }
146 }
147 m_force = false;
148 updateMixMoles();
149
150 // At this point, the instance has been created, the species to be included
151 // have been determined, and an initial composition has been selected that
152 // has all non-zero mole numbers for the included species.
153}
154
155double MultiPhaseEquil::equilibrate(int XY, double err, int maxsteps, int loglevel)
156{
157 int i;
158 m_iter = 0;
159 for (i = 0; i < maxsteps; i++) {
160 stepComposition(loglevel-1);
161 if (error() < err) {
162 break;
163 }
164 }
165 if (i >= maxsteps) {
166 throw CanteraError("MultiPhaseEquil::equilibrate",
167 "no convergence in {} iterations. Error = {}",
168 maxsteps, error());
169 }
170 finish();
171 return error();
172}
173
174void MultiPhaseEquil::updateMixMoles()
175{
176 fill(m_work3.begin(), m_work3.end(), 0.0);
177 for (size_t k = 0; k < m_nsp; k++) {
178 m_work3[m_species[k]] = m_moles[k];
179 }
180 m_mix->setMoles(m_work3.data());
181}
182
184{
185 fill(m_work3.begin(), m_work3.end(), 0.0);
186 for (size_t k = 0; k < m_nsp; k++) {
187 m_work3[m_species[k]] = (m_moles[k] > 0.0 ? m_moles[k] : 0.0);
188 }
189 m_mix->setMoles(m_work3.data());
190}
191
193{
194 double not_mu = 1.0e12;
195 m_mix->getValidChemPotentials(not_mu, m_mu.data(), true);
196 double dxi_min = 1.0e10;
197 bool redo = true;
198 int iter = 0;
199
200 while (redo) {
201 // choose a set of components based on the current composition
202 computeN();
203 redo = false;
204 iter++;
205 if (iter > 4) {
206 break;
207 }
208
209 // loop over all reactions
210 for (size_t j = 0; j < nFree(); j++) {
211 double dg_rt = 0.0;
212 dxi_min = 1.0e10;
213 for (size_t ik = 0; ik < m_nsp; ik++) {
214 dg_rt += mu(ik) * m_N(ik,j);
215 }
216
217 // fwd or rev direction
218 int idir = (dg_rt < 0.0 ? 1 : -1);
219
220 for (size_t ik = 0; ik < m_nsp; ik++) {
221 double nu = m_N(ik, j);
222
223 // set max change in progress variable by
224 // non-negativity requirement
225 // -> Note, 0.99 factor is so that difference of 2 numbers
226 // isn't zero. This causes differences between
227 // optimized and debug versions of the code
228 if (nu*idir < 0) {
229 double delta_xi = fabs(0.99*moles(ik)/nu);
230 // if a component has nearly zero moles, redo
231 // with a new set of components
232 if (!redo && delta_xi < 1.0e-10 && ik < m_nel) {
233 redo = true;
234 }
235 dxi_min = std::min(dxi_min, delta_xi);
236 }
237 }
238 // step the composition by dxi_min
239 for (size_t ik = 0; ik < m_nsp; ik++) {
240 moles(ik) += m_N(ik, j) * idir*dxi_min;
241 }
242 }
243 // set the moles of the phase objects to match
244 updateMixMoles();
245 }
246 return 0;
247}
248
249void MultiPhaseEquil::getComponents(const vector<size_t>& order)
250{
251 // if the input species array has the wrong size, ignore it
252 // and consider the species for components in declaration order.
253 if (order.size() != m_nsp) {
254 for (size_t k = 0; k < m_nsp; k++) {
255 m_order[k] = k;
256 }
257 } else {
258 for (size_t k = 0; k < m_nsp; k++) {
259 m_order[k] = order[k];
260 }
261 }
262
263 size_t nRows = m_nel;
264 size_t nColumns = m_nsp;
265
266 // set up the atomic composition matrix
267 for (size_t m = 0; m < nRows; m++) {
268 for (size_t k = 0; k < nColumns; k++) {
269 m_A(m, k) = m_mix->nAtoms(m_species[m_order[k]], m_element[m]);
270 }
271 }
272
273 // Do Gaussian elimination
274 for (size_t m = 0; m < nRows; m++) {
275 // Check for rows that are zero
276 bool isZeroRow = true;
277 for (size_t k = m; k < nColumns; k++) {
278 if (fabs(m_A(m,k)) > sqrt(Tiny)) {
279 isZeroRow = false;
280 break;
281 }
282 }
283 if (isZeroRow) {
284 // Find the last non-zero row
285 size_t n = nRows - 1;
286 bool foundSwapCandidate = false;
287 for (; n > m; n--) {
288 for (size_t k = m; k < nColumns; k++) {
289 if (fabs(m_A(n,k)) > sqrt(Tiny)) {
290 foundSwapCandidate = true;
291 break;
292 }
293 }
294 if (foundSwapCandidate) {
295 break;
296 }
297 }
298 if (m != n) {
299 // Swap this row with the last non-zero row
300 for (size_t k = 0; k < nColumns; k++) {
301 std::swap(m_A(n,k), m_A(m,k));
302 }
303 } else {
304 // All remaining rows are zero. Elimination is complete.
305 break;
306 }
307 }
308
309 // If a pivot is zero, exchange columns. This occurs when a species has
310 // an elemental composition that is not linearly independent of the
311 // component species that have already been assigned
312 if (m < nColumns && m_A(m,m) == 0.0) {
313 // First, we need to find a good candidate for a component species
314 // to swap in for the one that has zero pivot. It must contain
315 // element m, be linearly independent of the components processed so
316 // far (m_A(m,k) != 0), and should be a major species if possible.
317 // We'll choose the species with greatest mole fraction that
318 // satisfies these criteria.
319 double maxmoles = -999.0;
320 size_t kmax = 0;
321 for (size_t k = m+1; k < nColumns; k++) {
322 if (m_A(m,k) != 0.0 && fabs(m_moles[m_order[k]]) > maxmoles) {
323 kmax = k;
324 maxmoles = fabs(m_moles[m_order[k]]);
325 }
326 }
327
328 // Now exchange the column with zero pivot with the
329 // column for this major species
330 for (size_t n = 0; n < nRows; n++) {
331 std::swap(m_A(n, m), m_A(n, kmax));
332 }
333
334 // exchange the species labels on the columns
335 std::swap(m_order[m], m_order[kmax]);
336 }
337
338 // scale row m so that the diagonal element is unity
339 double fctr = 1.0/m_A(m,m);
340 for (size_t k = 0; k < nColumns; k++) {
341 m_A(m,k) *= fctr;
342 }
343
344 // For all rows below the diagonal, subtract A(n,m)/A(m,m)
345 // * (row m) from row n, so that A(n,m) = 0.
346 for (size_t n = m+1; n < m_nel; n++) {
347 fctr = m_A(n,m)/m_A(m,m);
348 for (size_t k = 0; k < m_nsp; k++) {
349 m_A(n,k) -= m_A(m,k)*fctr;
350 }
351 }
352 }
353
354 // The left m_nel columns of A are now upper-diagonal. Now
355 // reduce the m_nel columns to diagonal form by back-solving
356 for (size_t m = std::min(nRows,nColumns)-1; m > 0; m--) {
357 for (size_t n = m-1; n != npos; n--) {
358 if (m_A(n,m) != 0.0) {
359 double fctr = m_A(n,m);
360 for (size_t k = m; k < m_nsp; k++) {
361 m_A(n,k) -= fctr*m_A(m,k);
362 }
363 }
364 }
365 }
366
367 // create stoichiometric coefficient matrix.
368 for (size_t n = 0; n < m_nsp; n++) {
369 if (n < m_nel) {
370 for (size_t k = 0; k < nFree(); k++) {
371 m_N(n, k) = -m_A(n, k + m_nel);
372 }
373 } else {
374 for (size_t k = 0; k < nFree(); k++) {
375 m_N(n, k) = 0.0;
376 }
377 m_N(n, n - m_nel) = 1.0;
378 }
379 }
380
381 // find reactions involving solution phase species
382 for (size_t j = 0; j < nFree(); j++) {
383 m_solnrxn[j] = false;
384 for (size_t k = 0; k < m_nsp; k++) {
385 if (m_N(k, j) != 0 && m_mix->solutionSpecies(m_species[m_order[k]])) {
386 m_solnrxn[j] = true;
387 }
388 }
389 }
390}
391
392void MultiPhaseEquil::unsort(vector<double>& x)
393{
394 m_work2 = x;
395 for (size_t k = 0; k < m_nsp; k++) {
396 x[m_order[k]] = m_work2[k];
397 }
398}
399
400void MultiPhaseEquil::step(double omega, vector<double>& deltaN, int loglevel)
401{
402 if (omega < 0.0) {
403 throw CanteraError("MultiPhaseEquil::step","negative omega");
404 }
405
406 for (size_t ik = 0; ik < m_nel; ik++) {
407 size_t k = m_order[ik];
408 m_lastmoles[k] = m_moles[k];
409 m_moles[k] += omega * deltaN[k];
410 }
411
412 for (size_t ik = m_nel; ik < m_nsp; ik++) {
413 size_t k = m_order[ik];
414 m_lastmoles[k] = m_moles[k];
415 if (m_majorsp[k]) {
416 m_moles[k] += omega * deltaN[k];
417 } else {
418 m_moles[k] = fabs(m_moles[k])*std::min(10.0,
419 exp(-m_deltaG_RT[ik - m_nel]));
420 }
421 }
422 updateMixMoles();
423}
424
426{
427 m_iter++;
428 double grad0 = computeReactionSteps(m_dxi);
429
430 // compute the mole fraction changes.
431 if (nFree()) {
432 multiply(m_N, m_dxi.data(), m_work.data());
433 }
434
435 // change to sequential form
436 unsort(m_work);
437
438 // scale omega to keep the major species non-negative
439 double FCTR = 0.99;
440 const double MAJOR_THRESHOLD = 1.0e-12;
441 double omegamax = 1.0;
442 for (size_t ik = 0; ik < m_nsp; ik++) {
443 size_t k = m_order[ik];
444 if (ik < m_nel) {
445 FCTR = 0.99;
446 if (m_moles[k] < MAJOR_THRESHOLD) {
447 m_force = true;
448 }
449 } else {
450 FCTR = 0.9;
451 }
452 // if species k is in a multi-species solution phase, then its mole
453 // number must remain positive, unless the entire phase goes away. First
454 // we'll determine an upper bound on omega, such that all
455 if (m_dsoln[k] == 1) {
456 if ((m_moles[k] > MAJOR_THRESHOLD) || (ik < m_nel)) {
457 if (m_moles[k] < MAJOR_THRESHOLD) {
458 m_force = true;
459 }
460 double omax = m_moles[k]*FCTR/(fabs(m_work[k]) + Tiny);
461 if (m_work[k] < 0.0 && omax < omegamax) {
462 omegamax = omax;
463 if (omegamax < 1.0e-5) {
464 m_force = true;
465 }
466 }
467 m_majorsp[k] = true;
468 } else {
469 m_majorsp[k] = false;
470 }
471 } else {
472 if (m_work[k] < 0.0 && m_moles[k] > 0.0) {
473 double omax = -m_moles[k]/m_work[k];
474 if (omax < omegamax) {
475 omegamax = omax;
476 if (omegamax < 1.0e-5) {
477 m_force = true;
478 }
479 }
480 }
481 m_majorsp[k] = true;
482 }
483 }
484
485 // now take a step with this scaled omega
486 step(omegamax, m_work);
487 // compute the gradient of G at this new position in the current direction.
488 // If it is positive, then we have overshot the minimum. In this case,
489 // interpolate back.
490 double not_mu = 1.0e12;
491 m_mix->getValidChemPotentials(not_mu, m_mu.data());
492 double grad1 = 0.0;
493 for (size_t k = 0; k < m_nsp; k++) {
494 grad1 += m_work[k] * m_mu[m_species[k]];
495 }
496
497 double omega = omegamax;
498 if (grad1 > 0.0) {
499 omega *= fabs(grad0) / (grad1 + fabs(grad0));
500 for (size_t k = 0; k < m_nsp; k++) {
501 m_moles[k] = m_lastmoles[k];
502 }
503 step(omega, m_work);
504 }
505 return omega;
506}
507
508double MultiPhaseEquil::computeReactionSteps(vector<double>& dxi)
509{
510 vector<double> nu;
511 double grad = 0.0;
512 dxi.resize(nFree());
513 computeN();
514 double not_mu = 1.0e12;
515 m_mix->getValidChemPotentials(not_mu, m_mu.data());
516
517 for (size_t j = 0; j < nFree(); j++) {
518 // get stoichiometric vector
519 getStoichVector(j, nu);
520
521 // compute Delta G
522 double dg_rt = 0.0;
523 for (size_t k = 0; k < m_nsp; k++) {
524 dg_rt += m_mu[m_species[k]] * nu[k];
525 }
526 dg_rt /= (m_temp * GasConstant);
527
528 m_deltaG_RT[j] = dg_rt;
529 double fctr = 1.0;
530
531 // if this is a formation reaction for a single-component phase,
532 // check whether reaction should be included
533 size_t k = m_order[j + m_nel];
534 if (!m_dsoln[k]) {
535 if (m_moles[k] <= 0.0 && dg_rt > 0.0) {
536 fctr = 0.0;
537 } else {
538 fctr = 0.5;
539 }
540 } else if (!m_solnrxn[j]) {
541 fctr = 1.0;
542 } else {
543 // component sum
544 double csum = 0.0;
545 for (k = 0; k < m_nel; k++) {
546 size_t kc = m_order[k];
547 double nmoles = fabs(m_mix->speciesMoles(m_species[kc])) + Tiny;
548 csum += pow(nu[kc], 2)*m_dsoln[kc]/nmoles;
549 }
550
551 // noncomponent term
552 size_t kc = m_order[j + m_nel];
553 double nmoles = fabs(m_mix->speciesMoles(m_species[kc])) + Tiny;
554 double term1 = m_dsoln[kc]/nmoles;
555
556 // sum over solution phases
557 double sum = 0.0;
558 for (size_t ip = 0; ip < m_mix->nPhases(); ip++) {
559 ThermoPhase& p = m_mix->phase(ip);
560 if (p.nSpecies() > 1) {
561 double psum = 0.0;
562 for (k = 0; k < m_nsp; k++) {
563 kc = m_species[k];
564 if (m_mix->speciesPhaseIndex(kc) == ip) {
565 psum += pow(nu[k], 2);
566 }
567 }
568 sum -= psum / (fabs(m_mix->phaseMoles(ip)) + Tiny);
569 }
570 }
571 double rfctr = term1 + csum + sum;
572 if (fabs(rfctr) < Tiny) {
573 fctr = 1.0;
574 } else {
575 fctr = 1.0/(term1 + csum + sum);
576 }
577 }
578 dxi[j] = -fctr*dg_rt;
579
580 for (size_t m = 0; m < m_nel; m++) {
581 if (m_moles[m_order[m]] <= 0.0 && (m_N(m, j)*dxi[j] < 0.0)) {
582 dxi[j] = 0.0;
583 }
584 }
585 grad += dxi[j]*dg_rt;
586
587 }
588 return grad*GasConstant*m_temp;
589}
590
591void MultiPhaseEquil::computeN()
592{
593 // Sort the list of species by mole fraction (decreasing order)
594 vector<pair<double, size_t>> moleFractions(m_nsp);
595 for (size_t k = 0; k < m_nsp; k++) {
596 // use -Xk to generate reversed sort order
597 moleFractions[k] = {-m_mix->speciesMoles(m_species[k]), k};
598 }
599 std::sort(moleFractions.begin(), moleFractions.end());
600 for (size_t k = 0; k < m_nsp; k++) {
601 m_sortindex[k] = moleFractions[k].second;
602 }
603
604 for (size_t m = 0; m < m_nel; m++) {
605 size_t k = 0;
606 for (size_t ik = 0; ik < m_nsp; ik++) {
607 k = m_sortindex[ik];
608 if (m_mix->nAtoms(m_species[k],m_element[m]) != 0) {
609 break;
610 }
611 }
612 bool ok = false;
613 for (size_t ij = 0; ij < m_nel; ij++) {
614 if (k == m_order[ij]) {
615 ok = true;
616 }
617 }
618 if (!ok || m_force) {
619 getComponents(m_sortindex);
620 m_force = true;
621 break;
622 }
623 }
624}
625
626double MultiPhaseEquil::error()
627{
628 double err, maxerr = 0.0;
629
630 // examine every reaction
631 for (size_t j = 0; j < nFree(); j++) {
632 size_t ik = j + m_nel;
633
634 // don't require formation reactions for solution species
635 // present in trace amounts to be equilibrated
636 if (!isStoichPhase(ik) && fabs(moles(ik)) <= SmallNumber) {
637 err = 0.0;
638 }
639
640 // for stoichiometric phase species, no error if not present and
641 // delta G for the formation reaction is positive
642 if (isStoichPhase(ik) && moles(ik) <= 0.0 &&
643 m_deltaG_RT[j] >= 0.0) {
644 err = 0.0;
645 } else {
646 err = fabs(m_deltaG_RT[j]);
647 }
648 maxerr = std::max(maxerr, err);
649 }
650 return maxerr;
651}
652
653double MultiPhaseEquil::phaseMoles(size_t iph) const
654{
655 return m_mix->phaseMoles(iph);
656}
657
658void MultiPhaseEquil::reportCSV(const string& reportFile)
659{
660 FILE* FP = fopen(reportFile.c_str(), "w");
661 if (!FP) {
662 throw CanteraError("MultiPhaseEquil::reportCSV", "Failure to open file");
663 }
664 vector<double> mf(m_nsp_mix, 1.0);
665 vector<double> fe(m_nsp_mix, 0.0);
666 vector<double> VolPM;
667 vector<double> activity;
668 vector<double> ac;
669 vector<double> mu;
670 vector<double> mu0;
671 vector<double> molalities;
672
673 double vol = 0.0;
674 for (size_t iphase = 0; iphase < m_mix->nPhases(); iphase++) {
675 size_t istart = m_mix->speciesIndex(0, iphase);
676 ThermoPhase& tref = m_mix->phase(iphase);
677 size_t nSpecies = tref.nSpecies();
678 VolPM.resize(nSpecies, 0.0);
679 tref.getMoleFractions(&mf[istart]);
680 tref.getPartialMolarVolumes(VolPM.data());
681
682 double TMolesPhase = phaseMoles(iphase);
683 double VolPhaseVolumes = 0.0;
684 for (size_t k = 0; k < nSpecies; k++) {
685 VolPhaseVolumes += VolPM[k] * mf[istart + k];
686 }
687 VolPhaseVolumes *= TMolesPhase;
688 vol += VolPhaseVolumes;
689 }
690 fprintf(FP,"--------------------- VCS_MULTIPHASE_EQUIL FINAL REPORT"
691 " -----------------------------\n");
692 fprintf(FP,"Temperature = %11.5g kelvin\n", m_mix->temperature());
693 fprintf(FP,"Pressure = %11.5g Pascal\n", m_mix->pressure());
694 fprintf(FP,"Total Volume = %11.5g m**3\n", vol);
695
696 for (size_t iphase = 0; iphase < m_mix->nPhases(); iphase++) {
697 size_t istart = m_mix->speciesIndex(0, iphase);
698 ThermoPhase& tref = m_mix->phase(iphase);
699 ThermoPhase* tp = &tref;
700 tp->getMoleFractions(&mf[istart]);
701 string phaseName = tref.name();
702 double TMolesPhase = phaseMoles(iphase);
703 size_t nSpecies = tref.nSpecies();
704 activity.resize(nSpecies, 0.0);
705 ac.resize(nSpecies, 0.0);
706 mu0.resize(nSpecies, 0.0);
707 mu.resize(nSpecies, 0.0);
708 VolPM.resize(nSpecies, 0.0);
709 molalities.resize(nSpecies, 0.0);
710 int actConvention = tp->activityConvention();
711 tp->getActivities(activity.data());
712 tp->getActivityCoefficients(ac.data());
713 tp->getStandardChemPotentials(mu0.data());
714 tp->getPartialMolarVolumes(VolPM.data());
715 tp->getChemPotentials(mu.data());
716 double VolPhaseVolumes = 0.0;
717 for (size_t k = 0; k < nSpecies; k++) {
718 VolPhaseVolumes += VolPM[k] * mf[istart + k];
719 }
720 VolPhaseVolumes *= TMolesPhase;
721 vol += VolPhaseVolumes;
722 if (actConvention == 1) {
723 MolalityVPSSTP* mTP = static_cast<MolalityVPSSTP*>(tp);
724 mTP->getMolalities(molalities.data());
725 tp->getChemPotentials(mu.data());
726
727 if (iphase == 0) {
728 fprintf(FP," Name, Phase, PhaseMoles, Mole_Fract, "
729 "Molalities, ActCoeff, Activity,"
730 "ChemPot_SS0, ChemPot, mole_num, PMVol, Phase_Volume\n");
731
732 fprintf(FP," , , (kmol), , "
733 ", , ,"
734 " (kJ/gmol), (kJ/gmol), (kmol), (m**3/kmol), (m**3)\n");
735 }
736 for (size_t k = 0; k < nSpecies; k++) {
737 string sName = tp->speciesName(k);
738 fprintf(FP,"%12s, %11s, %11.3e, %11.3e, %11.3e, %11.3e, %11.3e,"
739 "%11.3e, %11.3e, %11.3e, %11.3e, %11.3e\n",
740 sName.c_str(),
741 phaseName.c_str(), TMolesPhase,
742 mf[istart + k], molalities[k], ac[k], activity[k],
743 mu0[k]*1.0E-6, mu[k]*1.0E-6,
744 mf[istart + k] * TMolesPhase,
745 VolPM[k], VolPhaseVolumes);
746 }
747 } else {
748 if (iphase == 0) {
749 fprintf(FP," Name, Phase, PhaseMoles, Mole_Fract, "
750 "Molalities, ActCoeff, Activity,"
751 " ChemPotSS0, ChemPot, mole_num, PMVol, Phase_Volume\n");
752
753 fprintf(FP," , , (kmol), , "
754 ", , ,"
755 " (kJ/gmol), (kJ/gmol), (kmol), (m**3/kmol), (m**3)\n");
756 }
757 for (size_t k = 0; k < nSpecies; k++) {
758 molalities[k] = 0.0;
759 }
760 for (size_t k = 0; k < nSpecies; k++) {
761 string sName = tp->speciesName(k);
762 fprintf(FP,"%12s, %11s, %11.3e, %11.3e, %11.3e, %11.3e, %11.3e, "
763 "%11.3e, %11.3e,% 11.3e, %11.3e, %11.3e\n",
764 sName.c_str(),
765 phaseName.c_str(), TMolesPhase,
766 mf[istart + k], molalities[k], ac[k],
767 activity[k], mu0[k]*1.0E-6, mu[k]*1.0E-6,
768 mf[istart + k] * TMolesPhase,
769 VolPM[k], VolPhaseVolumes);
770 }
771 }
772 }
773 fclose(FP);
774}
775
776}
Header for intermediate ThermoPhase object for phases which employ molality based activity coefficien...
vector< double > & data()
Return a reference to the data vector.
Definition Array.h:186
Base class for exceptions thrown by Cantera classes.
void resize(size_t n, size_t m, double v=0.0) override
Resize the matrix.
double computeReactionSteps(vector< double > &dxi)
Compute the change in extent of reaction for each reaction.
double stepComposition(int loglevel=0)
Take one step in composition, given the gradient of G at the starting point, and a vector of reaction...
void finish()
Clean up the composition.
void unsort(vector< double > &x)
Re-arrange a vector of species properties in sorted form (components first) into unsorted,...
size_t nFree() const
Number of degrees of freedom.
MultiPhaseEquil(MultiPhase *mix, bool start=true, int loglevel=0)
Construct a multiphase equilibrium manager for a multiphase mixture.
int setInitialMoles(int loglevel=0)
Estimate the initial mole numbers.
void getComponents(const vector< size_t > &order)
This method finds a set of component species and a complete set of formation reactions for the non-co...
A class for multiphase mixtures.
Definition MultiPhase.h:57
size_t speciesIndex(size_t k, size_t p) const
Return the global index of the species belonging to phase number p with local index k within the phas...
Definition MultiPhase.h:226
bool solutionSpecies(size_t kGlob) const
Return true is species kGlob is a species in a multicomponent solution phase.
double nAtoms(const size_t kGlob, const size_t mGlob) const
Returns the Number of atoms of global element mGlob in global species kGlob.
void setMoles(const double *n)
Sets all of the global species mole numbers.
double speciesMoles(size_t kGlob) const
Returns the moles of global species k. units = kmol.
void getValidChemPotentials(double not_mu, double *mu, bool standard=false) const
Returns a vector of Valid chemical potentials.
size_t nSpecies() const
Number of species, summed over all phases.
Definition MultiPhase.h:123
double pressure() const
Pressure [Pa].
Definition MultiPhase.h:388
size_t nPhases() const
Number of phases.
Definition MultiPhase.h:425
double temperature() const
Temperature [K].
Definition MultiPhase.h:328
bool tempOK(size_t p) const
Return true if the phase p has valid thermo data for the current temperature.
size_t speciesPhaseIndex(const size_t kGlob) const
Returns the phase index of the Kth "global" species.
string speciesName(const size_t kGlob) const
Name of species with global index kGlob.
double phaseMoles(const size_t n) const
Return the number of moles in phase n.
size_t nElements() const
Number of elements.
Definition MultiPhase.h:97
ThermoPhase & phase(size_t n)
Return a reference to phase n.
string elementName(size_t m) const
Returns the name of the global element m.
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:231
void getMoleFractions(double *const x) const
Get the species mole fraction vector.
Definition Phase.cpp:434
Base class for a phase with thermodynamic properties.
const double GasConstant
Universal Gas Constant [J/kmol/K].
Definition ct_defs.h:120
Namespace for the Cantera kernel.
Definition AnyMap.cpp:564
const size_t npos
index returned by functions to indicate "no position"
Definition ct_defs.h:180
const double Tiny
Small number to compare differences of mole fractions against.
Definition ct_defs.h:173
void multiply(const DenseMatrix &A, const double *const b, double *const prod)
Multiply A*b and return the result in prod. Uses BLAS routine DGEMV.
const double SmallNumber
smallest number to compare to zero.
Definition ct_defs.h:158
Contains declarations for string manipulation functions within Cantera.