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++) {
32 if (enm ==
"E" || enm ==
"e") {
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;
51 if (m_eloc < m_nel_mix) {
52 m_element.push_back(m_eloc);
56 for (
size_t m = 0; m < m_nel_mix; m++) {
57 if (m_incl_element[m] == 1 && m != m_eloc) {
59 m_element.push_back(m);
72 for (
size_t k = 0; k < m_nsp_mix; k++) {
76 m_incl_species[k] = 0;
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.");
88 for (
size_t k = 0; k < m_nsp_mix; k++) {
89 if (m_incl_species[k] ==1) {
91 m_species.push_back(k);
97 m_work2.resize(m_nsp);
98 m_work3.resize(m_nsp_mix);
99 m_mu.resize(m_nsp_mix);
102 m_moles.resize(m_nsp);
103 m_lastmoles.resize(m_nsp);
104 m_dxi.resize(
nFree());
107 for (
size_t ik = 0; ik < m_nsp; ik++) {
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);
119 m_order.resize(std::max(m_nsp, m_nel), 0);
120 iota(m_order.begin(), m_order.begin() + m_nsp, 0);
132 vector<double> dxi(
nFree(), 1.0e-20);
138 for (
size_t k = 0; k < m_nsp; k++) {
139 m_moles[k] += m_work[k];
140 m_lastmoles[k] = m_moles[k];
142 m_dsoln.push_back(1);
144 m_dsoln.push_back(0);
155double MultiPhaseEquil::equilibrate(
int XY,
double err,
int maxsteps,
int loglevel)
159 for (i = 0; i < maxsteps; i++) {
166 throw CanteraError(
"MultiPhaseEquil::equilibrate",
167 "no convergence in {} iterations. Error = {}",
174void MultiPhaseEquil::updateMixMoles()
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];
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);
194 double not_mu = 1.0e12;
196 double dxi_min = 1.0e10;
210 for (
size_t j = 0; j <
nFree(); j++) {
213 for (
size_t ik = 0; ik < m_nsp; ik++) {
214 dg_rt += mu(ik) * m_N(ik,j);
218 int idir = (dg_rt < 0.0 ? 1 : -1);
220 for (
size_t ik = 0; ik < m_nsp; ik++) {
221 double nu = m_N(ik, j);
229 double delta_xi = fabs(0.99*moles(ik)/nu);
232 if (!redo && delta_xi < 1.0e-10 && ik < m_nel) {
235 dxi_min = std::min(dxi_min, delta_xi);
239 for (
size_t ik = 0; ik < m_nsp; ik++) {
240 moles(ik) += m_N(ik, j) * idir*dxi_min;
253 if (order.size() != m_nsp) {
254 for (
size_t k = 0; k < m_nsp; k++) {
258 for (
size_t k = 0; k < m_nsp; k++) {
259 m_order[k] = order[k];
263 size_t nRows = m_nel;
264 size_t nColumns = m_nsp;
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]);
274 for (
size_t m = 0; m < nRows; m++) {
276 bool isZeroRow =
true;
277 for (
size_t k = m; k < nColumns; k++) {
278 if (fabs(m_A(m,k)) > sqrt(
Tiny)) {
285 size_t n = nRows - 1;
286 bool foundSwapCandidate =
false;
288 for (
size_t k = m; k < nColumns; k++) {
289 if (fabs(m_A(n,k)) > sqrt(
Tiny)) {
290 foundSwapCandidate =
true;
294 if (foundSwapCandidate) {
300 for (
size_t k = 0; k < nColumns; k++) {
301 std::swap(m_A(n,k), m_A(m,k));
312 if (m < nColumns && m_A(m,m) == 0.0) {
319 double maxmoles = -999.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) {
324 maxmoles = fabs(m_moles[m_order[k]]);
330 for (
size_t n = 0; n < nRows; n++) {
331 std::swap(m_A(n, m), m_A(n, kmax));
335 std::swap(m_order[m], m_order[kmax]);
339 double fctr = 1.0/m_A(m,m);
340 for (
size_t k = 0; k < nColumns; k++) {
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;
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);
368 for (
size_t n = 0; n < m_nsp; n++) {
370 for (
size_t k = 0; k <
nFree(); k++) {
371 m_N(n, k) = -m_A(n, k + m_nel);
374 for (
size_t k = 0; k <
nFree(); k++) {
377 m_N(n, n - m_nel) = 1.0;
382 for (
size_t j = 0; j <
nFree(); j++) {
383 m_solnrxn[j] =
false;
384 for (
size_t k = 0; k < m_nsp; k++) {
395 m_work2.assign(x.begin(), x.end());
396 for (
size_t k = 0; k < m_nsp; k++) {
397 x[m_order[k]] = m_work2[k];
401void MultiPhaseEquil::step(
double omega, span<double> deltaN,
int loglevel)
404 throw CanteraError(
"MultiPhaseEquil::step",
"negative omega");
406 if (deltaN.size() != m_nsp) {
407 throw CanteraError(
"MultiPhaseEquil::step",
408 "Expected deltaN size {}, got {}", m_nsp, deltaN.size());
411 for (
size_t ik = 0; ik < m_nel; ik++) {
412 size_t k = m_order[ik];
413 m_lastmoles[k] = m_moles[k];
414 m_moles[k] += omega * deltaN[k];
417 for (
size_t ik = m_nel; ik < m_nsp; ik++) {
418 size_t k = m_order[ik];
419 m_lastmoles[k] = m_moles[k];
421 m_moles[k] += omega * deltaN[k];
423 m_moles[k] = fabs(m_moles[k])*std::min(10.0,
424 exp(-m_deltaG_RT[ik - m_nel]));
445 const double MAJOR_THRESHOLD = 1.0e-12;
446 double omegamax = 1.0;
447 for (
size_t ik = 0; ik < m_nsp; ik++) {
448 size_t k = m_order[ik];
451 if (m_moles[k] < MAJOR_THRESHOLD) {
460 if (m_dsoln[k] == 1) {
461 if ((m_moles[k] > MAJOR_THRESHOLD) || (ik < m_nel)) {
462 if (m_moles[k] < MAJOR_THRESHOLD) {
465 double omax = m_moles[k]*FCTR/(fabs(m_work[k]) +
Tiny);
466 if (m_work[k] < 0.0 && omax < omegamax) {
468 if (omegamax < 1.0e-5) {
474 m_majorsp[k] =
false;
477 if (m_work[k] < 0.0 && m_moles[k] > 0.0) {
478 double omax = -m_moles[k]/m_work[k];
479 if (omax < omegamax) {
481 if (omegamax < 1.0e-5) {
491 step(omegamax, m_work);
495 double not_mu = 1.0e12;
498 for (
size_t k = 0; k < m_nsp; k++) {
499 grad1 += m_work[k] * m_mu[m_species[k]];
502 double omega = omegamax;
504 omega *= fabs(grad0) / (grad1 + fabs(grad0));
505 for (
size_t k = 0; k < m_nsp; k++) {
506 m_moles[k] = m_lastmoles[k];
516 vector<double> nu(m_nsp, 0.0);
519 double not_mu = 1.0e12;
522 for (
size_t j = 0; j <
nFree(); j++) {
524 getStoichVector(j, nu);
528 for (
size_t k = 0; k < m_nsp; k++) {
529 dg_rt += m_mu[m_species[k]] * nu[k];
533 m_deltaG_RT[j] = dg_rt;
538 size_t k = m_order[j + m_nel];
540 if (m_moles[k] <= 0.0 && dg_rt > 0.0) {
545 }
else if (!m_solnrxn[j]) {
550 for (k = 0; k < m_nel; k++) {
551 size_t kc = m_order[k];
553 csum += pow(nu[kc], 2)*m_dsoln[kc]/nmoles;
557 size_t kc = m_order[j + m_nel];
559 double term1 = m_dsoln[kc]/nmoles;
563 for (
size_t ip = 0; ip < m_mix->
nPhases(); ip++) {
567 for (k = 0; k < m_nsp; k++) {
570 psum += pow(nu[k], 2);
576 double rfctr = term1 + csum + sum;
577 if (fabs(rfctr) <
Tiny) {
580 fctr = 1.0/(term1 + csum + sum);
583 dxi[j] = -fctr*dg_rt;
585 for (
size_t m = 0; m < m_nel; m++) {
586 if (m_moles[m_order[m]] <= 0.0 && (m_N(m, j)*dxi[j] < 0.0)) {
590 grad += dxi[j]*dg_rt;
596void MultiPhaseEquil::computeN()
599 vector<pair<double, size_t>> moleFractions(m_nsp);
600 for (
size_t k = 0; k < m_nsp; k++) {
602 moleFractions[k] = {-m_mix->
speciesMoles(m_species[k]), k};
604 std::sort(moleFractions.begin(), moleFractions.end());
605 for (
size_t k = 0; k < m_nsp; k++) {
606 m_sortindex[k] = moleFractions[k].second;
609 for (
size_t m = 0; m < m_nel; m++) {
611 for (
size_t ik = 0; ik < m_nsp; ik++) {
613 if (m_mix->
nAtoms(m_species[k],m_element[m]) != 0) {
618 for (
size_t ij = 0; ij < m_nel; ij++) {
619 if (k == m_order[ij]) {
623 if (!ok || m_force) {
631double MultiPhaseEquil::error()
633 double err, maxerr = 0.0;
636 for (
size_t j = 0; j <
nFree(); j++) {
637 size_t ik = j + m_nel;
641 if (!isStoichPhase(ik) && fabs(moles(ik)) <=
SmallNumber) {
647 if (isStoichPhase(ik) && moles(ik) <= 0.0 &&
648 m_deltaG_RT[j] >= 0.0) {
651 err = fabs(m_deltaG_RT[j]);
653 maxerr = std::max(maxerr, err);
658double MultiPhaseEquil::phaseMoles(
size_t iph)
const
663void MultiPhaseEquil::reportCSV(
const string& reportFile)
665 FILE* FP = fopen(reportFile.c_str(),
"w");
667 throw CanteraError(
"MultiPhaseEquil::reportCSV",
"Failure to open file");
669 vector<double> mf(m_nsp_mix, 1.0);
670 vector<double> fe(m_nsp_mix, 0.0);
671 vector<double> VolPM;
672 vector<double> activity;
676 vector<double> molalities;
679 for (
size_t iphase = 0; iphase < m_mix->
nPhases(); iphase++) {
681 ThermoPhase& tref = m_mix->
phase(iphase);
683 VolPM.resize(nSpecies, 0.0);
684 tref.getMoleFractions(span<double>(&mf[istart], tref.nSpecies()));
685 tref.getPartialMolarVolumes(VolPM);
687 double TMolesPhase = phaseMoles(iphase);
688 double VolPhaseVolumes = 0.0;
689 for (
size_t k = 0; k < nSpecies; k++) {
690 VolPhaseVolumes += VolPM[k] * mf[istart + k];
692 VolPhaseVolumes *= TMolesPhase;
693 vol += VolPhaseVolumes;
695 fprintf(FP,
"--------------------- VCS_MULTIPHASE_EQUIL FINAL REPORT"
696 " -----------------------------\n");
697 fprintf(FP,
"Temperature = %11.5g kelvin\n", m_mix->
temperature());
698 fprintf(FP,
"Pressure = %11.5g Pascal\n", m_mix->
pressure());
699 fprintf(FP,
"Total Volume = %11.5g m**3\n", vol);
701 for (
size_t iphase = 0; iphase < m_mix->
nPhases(); iphase++) {
703 ThermoPhase& tref = m_mix->
phase(iphase);
704 ThermoPhase* tp = &tref;
706 string phaseName = tref.name();
707 double TMolesPhase = phaseMoles(iphase);
708 size_t nSpecies = tref.nSpecies();
709 activity.resize(nSpecies, 0.0);
710 ac.resize(nSpecies, 0.0);
711 mu0.resize(nSpecies, 0.0);
712 mu.resize(nSpecies, 0.0);
713 VolPM.resize(nSpecies, 0.0);
714 molalities.resize(nSpecies, 0.0);
715 int actConvention = tp->activityConvention();
716 tp->getActivities(activity);
717 tp->getActivityCoefficients(ac);
718 tp->getStandardChemPotentials(mu0);
719 tp->getPartialMolarVolumes(VolPM);
720 tp->getChemPotentials(mu);
721 double VolPhaseVolumes = 0.0;
722 for (
size_t k = 0; k < nSpecies; k++) {
723 VolPhaseVolumes += VolPM[k] * mf[istart + k];
725 VolPhaseVolumes *= TMolesPhase;
726 vol += VolPhaseVolumes;
727 if (actConvention == 1) {
728 MolalityVPSSTP* mTP =
static_cast<MolalityVPSSTP*
>(tp);
729 mTP->getMolalities(molalities);
730 tp->getChemPotentials(mu);
733 fprintf(FP,
" Name, Phase, PhaseMoles, Mole_Fract, "
734 "Molalities, ActCoeff, Activity,"
735 "ChemPot_SS0, ChemPot, mole_num, PMVol, Phase_Volume\n");
737 fprintf(FP,
" , , (kmol), , "
739 " (kJ/gmol), (kJ/gmol), (kmol), (m**3/kmol), (m**3)\n");
741 for (
size_t k = 0; k < nSpecies; k++) {
742 string sName = tp->speciesName(k);
743 fprintf(FP,
"%12s, %11s, %11.3e, %11.3e, %11.3e, %11.3e, %11.3e,"
744 "%11.3e, %11.3e, %11.3e, %11.3e, %11.3e\n",
746 phaseName.c_str(), TMolesPhase,
747 mf[istart + k], molalities[k], ac[k], activity[k],
748 mu0[k]*1.0E-6, mu[k]*1.0E-6,
749 mf[istart + k] * TMolesPhase,
750 VolPM[k], VolPhaseVolumes);
754 fprintf(FP,
" Name, Phase, PhaseMoles, Mole_Fract, "
755 "Molalities, ActCoeff, Activity,"
756 " ChemPotSS0, ChemPot, mole_num, PMVol, Phase_Volume\n");
758 fprintf(FP,
" , , (kmol), , "
760 " (kJ/gmol), (kJ/gmol), (kmol), (m**3/kmol), (m**3)\n");
762 for (
size_t k = 0; k < nSpecies; k++) {
765 for (
size_t k = 0; k < nSpecies; k++) {
766 string sName = tp->speciesName(k);
767 fprintf(FP,
"%12s, %11s, %11.3e, %11.3e, %11.3e, %11.3e, %11.3e, "
768 "%11.3e, %11.3e,% 11.3e, %11.3e, %11.3e\n",
770 phaseName.c_str(), TMolesPhase,
771 mf[istart + k], molalities[k], ac[k],
772 activity[k], mu0[k]*1.0E-6, mu[k]*1.0E-6,
773 mf[istart + k] * TMolesPhase,
774 VolPM[k], VolPhaseVolumes);
Header for intermediate ThermoPhase object for phases which employ molality based activity coefficien...
Base class for exceptions thrown by Cantera classes.
void resize(size_t n, size_t m, double v=0.0) override
Resize the matrix.
void getComponents(span< const size_t > order)
This method finds a set of component species and a complete set of formation reactions for the non-co...
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 unsort(span< double > x)
Re-arrange a vector of species properties in sorted form (components first) into unsorted,...
void finish()
Clean up the composition.
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.
double computeReactionSteps(span< double > dxi)
Compute the change in extent of reaction for each reaction.
int setInitialMoles(int loglevel=0)
Estimate the initial mole numbers.
A class for multiphase mixtures.
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...
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(span< 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.
size_t nSpecies() const
Number of species, summed over all phases.
double pressure() const
Pressure [Pa].
size_t nPhases() const
Number of phases.
double temperature() const
Temperature [K].
string speciesName(size_t kGlob) const
Name of species with global index kGlob.
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.
void getValidChemPotentials(double not_mu, span< double > mu, bool standard=false) const
Returns a vector of Valid chemical potentials.
double phaseMoles(const size_t n) const
Return the number of moles in phase n.
size_t nElements() const
Number of elements.
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.
void getMoleFractions(span< double > x) const
Get the species mole fraction vector.
size_t nSpecies() const
Returns the number of species in the phase.
Base class for a phase with thermodynamic properties.
const double GasConstant
Universal Gas Constant [J/kmol/K].
Namespace for the Cantera kernel.
const size_t npos
index returned by functions to indicate "no position"
const double Tiny
Small number to compare differences of mole fractions against.
void multiply(const DenseMatrix &A, span< const double > b, span< double > prod)
Multiply A*b and return the result in prod. Uses BLAS routine DGEMV.
const double SmallNumber
smallest number to compare to zero.
void checkArraySize(const char *procedure, size_t available, size_t required)
Wrapper for throwing ArraySizeError.
Contains declarations for string manipulation functions within Cantera.