4 #include "cantera/equil/MultiPhaseEquil.h"
19 const doublereal TINY = 1.0e-20;
21 #if defined(WITH_HTML_LOGS)
31 static string coeffString(
bool first, doublereal nu,
string sym)
40 if (nu == 1.0 || nu == -1.0) {
43 string s =
fp2str(fabs(nu));
44 return strt + s +
" " + sym;
54 MultiPhaseEquil::MultiPhaseEquil(
MultiPhase* mix,
bool start,
int loglevel) : m_mix(mix)
71 m_incl_species.resize(m_nsp_mix,1);
72 m_incl_element.resize(m_nel_mix,1);
73 for (m = 0; m < m_nel_mix; m++) {
78 if (enm ==
"E" || enm ==
"e") {
88 m_incl_element[m] = 0;
89 for (k = 0; k < m_nsp_mix; k++) {
90 if (m_mix->
nAtoms(k,m) != 0.0) {
91 m_incl_species[k] = 0;
100 if (m_eloc < m_nel_mix) {
101 m_element.push_back(m_eloc);
105 for (m = 0; m < m_nel_mix; m++) {
106 if (m_incl_element[m] == 1 && m != m_eloc) {
108 m_element.push_back(m);
124 for (k = 0; k < m_nsp_mix; k++) {
128 m_incl_species[k] = 0;
132 +
" is excluded since its thermo properties are \n"
133 "not valid at this temperature, but it has "
134 "non-zero moles in the initial state.");
141 for (k = 0; k < m_nsp_mix; k++) {
142 if (m_incl_species[k] ==1) {
144 m_species.push_back(k);
149 m_work.resize(m_nsp);
150 m_work2.resize(m_nsp);
151 m_work3.resize(m_nsp_mix);
152 m_mu.resize(m_nsp_mix);
155 m_moles.resize(m_nsp);
156 m_lastmoles.resize(m_nsp);
157 m_dxi.resize(
nFree());
161 for (ik = 0; ik < m_nsp; ik++) {
166 m_deltaG_RT.resize(
nFree(), 0.0);
168 m_majorsp.resize(m_nsp);
169 m_sortindex.resize(m_nsp,0);
170 m_lastsort.resize(m_nel);
171 m_solnrxn.resize(
nFree());
172 m_A.
resize(m_nel, m_nsp, 0.0);
174 m_order.resize(
std::max(m_nsp, m_nel), 0);
175 for (k = 0; k < m_nsp; k++) {
196 for (k = 0; k < m_nsp; k++) {
197 m_moles[k] += m_work[k];
198 m_lastmoles[k] = m_moles[k];
200 m_dsoln.push_back(1);
202 m_dsoln.push_back(0);
215 doublereal MultiPhaseEquil::equilibrate(
int XY, doublereal err,
216 int maxsteps,
int loglevel)
225 for (i = 0; i < maxsteps; i++) {
227 iterstr =
"iteration "+
int2str(i);
245 throw CanteraError(
"MultiPhaseEquil::equilibrate",
246 "no convergence in " +
int2str(maxsteps) +
247 " iterations. Error = " +
fp2str(error()));
259 void MultiPhaseEquil::updateMixMoles()
261 fill(m_work3.begin(), m_work3.end(), 0.0);
263 for (k = 0; k < m_nsp; k++) {
264 m_work3[m_species[k]] = m_moles[k];
275 fill(m_work3.begin(), m_work3.end(), 0.0);
277 for (k = 0; k < m_nsp; k++) {
278 m_work3[m_species[k]] = (m_moles[k] > 0.0 ? m_moles[k] : 0.0);
297 double not_mu = 1.0e12;
307 double delta_xi, dxi_min = 1.0e10;
326 for (j = 0; j <
nFree(); j++) {
329 for (ik = 0; ik < m_nsp; ik++) {
330 dg_rt += mu(ik) * m_N(ik,j);
334 idir = (dg_rt < 0.0 ? 1 : -1);
336 for (ik = 0; ik < m_nsp; ik++) {
345 delta_xi = fabs(0.99*moles(ik)/nu);
348 if (!redo && delta_xi < 1.0e-10 && ik < m_nel) {
350 addLogEntry(
"component too small",speciesName(ik));
354 if (delta_xi < dxi_min) {
360 for (ik = 0; ik < m_nsp; ik++) {
361 moles(ik) += m_N(ik, j) * idir*dxi_min;
367 for (ik = 0; ik < m_nsp; ik++)
368 if (moles(ik) != 0.0) {
405 if (order.size() != m_nsp) {
406 for (k = 0; k < m_nsp; k++) {
410 for (k = 0; k < m_nsp; k++) {
411 m_order[k] = order[k];
415 index_t nRows = m_nel;
416 index_t nColumns = m_nsp;
420 for (m = 0; m < nRows; m++) {
421 for (k = 0; k < nColumns; k++) {
422 m_A(m, k) = m_mix->
nAtoms(m_species[m_order[k]], m_element[m]);
427 for (m = 0; m < nRows; m++) {
429 bool isZeroRow =
true;
430 for (k = m; k < nColumns; k++) {
431 if (fabs(m_A(m,k)) > sqrt(TINY)) {
438 index_t n = nRows - 1;
439 bool foundSwapCandidate =
false;
441 for (k = m; k < nColumns; k++) {
442 if (fabs(m_A(n,k)) > sqrt(TINY)) {
443 foundSwapCandidate =
true;
447 if (foundSwapCandidate) {
453 for (k = 0; k < nColumns; k++) {
454 std::swap(m_A(n,k), m_A(m,k));
466 if (m < nColumns && m_A(m,m) == 0.0) {
474 doublereal maxmoles = -999.0;
476 for (k = m+1; k < nColumns; k++) {
477 if (m_A(m,k) != 0.0) {
478 if (fabs(m_moles[m_order[k]]) > maxmoles) {
480 maxmoles = fabs(m_moles[m_order[k]]);
487 for (
size_t n = 0; n < nRows; n++) {
488 std::swap(m_A(n, m), m_A(n, kmax));
492 std::swap(m_order[m], m_order[kmax]);
497 for (k = 0; k < nColumns; k++) {
503 for (
size_t n = m+1; n < m_nel; n++) {
504 fctr = m_A(n,m)/m_A(m,m);
505 for (k = 0; k < m_nsp; k++) {
506 m_A(n,k) -= m_A(m,k)*fctr;
513 for (m =
std::min(nRows,nColumns)-1; m > 0; m--) {
514 for (
size_t n = m-1; n !=
npos; n--) {
515 if (m_A(n,m) != 0.0) {
517 for (k = m; k < m_nsp; k++) {
518 m_A(n,k) -= fctr*m_A(m,k);
525 for (
size_t n = 0; n < m_nsp; n++) {
527 for (k = 0; k <
nFree(); k++) {
528 m_N(n, k) = -m_A(n, k + m_nel);
531 for (k = 0; k <
nFree(); k++) {
534 m_N(n, n - m_nel) = 1.0;
539 for (j = 0; j <
nFree(); j++) {
540 m_solnrxn[j] =
false;
541 for (k = 0; k < m_nsp; k++) {
557 copy(x.begin(), x.end(), m_work2.begin());
559 for (k = 0; k < m_nsp; k++) {
560 x[m_order[k]] = m_work2[k];
564 #if defined(WITH_HTML_LOGS)
565 void MultiPhaseEquil::printInfo(
int loglevel)
572 for (m = 0; m < m_nel; m++) {
583 for (m = m_nel; m < m_nsp; m++) {
595 for (k = 0; k <
nFree(); k++) {
609 string sr =
"", sp =
"";
614 for (i = 0; i < m_nsp; i++) {
616 k = m_species[m_order[i]];
626 return sr +
" <=> " + sp;
630 void MultiPhaseEquil::step(doublereal omega,
vector_fp& deltaN,
638 throw CanteraError(
"step",
"negative omega");
641 for (ik = 0; ik < m_nel; ik++) {
643 m_lastmoles[k] = m_moles[k];
650 m_moles[k] += omega * deltaN[k];
653 for (ik = m_nel; ik < m_nsp; ik++) {
655 m_lastmoles[k] = m_moles[k];
657 m_moles[k] += omega * deltaN[k];
659 m_moles[k] = fabs(m_moles[k])*
std::min(10.0,
660 exp(-m_deltaG_RT[ik - m_nel]));
692 doublereal FCTR = 0.99;
693 const doublereal MAJOR_THRESHOLD = 1.0e-12;
695 doublereal omega = 1.0, omax, omegamax = 1.0;
696 for (ik = 0; ik < m_nsp; ik++) {
700 if (m_moles[k] < MAJOR_THRESHOLD) {
710 if (m_dsoln[k] == 1) {
712 if ((m_moles[k] > MAJOR_THRESHOLD) || (ik < m_nel)) {
713 if (m_moles[k] < MAJOR_THRESHOLD) {
716 omax = m_moles[k]*FCTR/(fabs(m_work[k]) + TINY);
717 if (m_work[k] < 0.0 && omax < omegamax) {
719 if (omegamax < 1.0e-5) {
725 m_majorsp[k] =
false;
728 if (m_work[k] < 0.0 && m_moles[k] > 0.0) {
729 omax = -m_moles[k]/m_work[k];
730 if (omax < omegamax) {
732 if (omegamax < 1.0e-5) {
737 if (m_moles[k] < -
Tiny) {
750 step(omegamax, m_work);
754 doublereal not_mu = 1.0e12;
756 doublereal grad1 = 0.0;
757 for (k = 0; k < m_nsp; k++) {
758 grad1 += m_work[k] * m_mu[m_species[k]];
763 omega *= fabs(grad0) / (grad1 + fabs(grad0));
764 for (k = 0; k < m_nsp; k++) {
765 m_moles[k] = m_lastmoles[k];
785 index_t j, k, ik, kc, ip;
786 doublereal stoich, nmoles, csum, term1, fctr, rfctr;
788 const doublereal TINY = 1.0e-20;
789 doublereal grad = 0.0;
793 doublereal not_mu = 1.0e12;
796 for (j = 0; j <
nFree(); j++) {
799 getStoichVector(j, nu);
802 doublereal dg_rt = 0.0;
803 for (k = 0; k < m_nsp; k++) {
804 dg_rt += m_mu[m_species[k]] * nu[k];
808 m_deltaG_RT[j] = dg_rt;
816 if (m_moles[k] <= 0.0 && dg_rt > 0.0) {
821 }
else if (!m_solnrxn[j]) {
827 for (k = 0; k < m_nel; k++) {
830 nmoles = fabs(m_mix->
speciesMoles(m_species[kc])) + TINY;
831 csum += stoich*stoich*m_dsoln[kc]/nmoles;
835 kc = m_order[j + m_nel];
836 nmoles = fabs(m_mix->
speciesMoles(m_species[kc])) + TINY;
837 term1 = m_dsoln[kc]/nmoles;
840 doublereal sum = 0.0, psum;
841 for (ip = 0; ip < m_np; ip++) {
845 for (k = 0; k < m_nsp; k++) {
850 psum += stoich * stoich;
853 sum -= psum / (fabs(m_mix->
phaseMoles(ip)) + TINY);
856 rfctr = term1 + csum + sum;
857 if (fabs(rfctr) < TINY) {
860 fctr = 1.0/(term1 + csum + sum);
863 dxi[j] = -fctr*dg_rt;
866 for (m = 0; m < m_nel; m++) {
867 if (m_moles[m_order[m]] <= 0.0 && (m_N(m, j)*dxi[j] < 0.0)) {
871 grad += dxi[j]*dg_rt;
877 void MultiPhaseEquil::computeN()
880 std::vector<std::pair<double, size_t> > moleFractions(m_nsp);
881 for (
size_t k = 0; k < m_nsp; k++) {
883 moleFractions[k].first = - m_mix->
speciesMoles(m_species[k]);
884 moleFractions[k].second = k;
886 std::sort(moleFractions.begin(), moleFractions.end());
887 for (
size_t k = 0; k < m_nsp; k++) {
888 m_sortindex[k] = moleFractions[k].second;
892 for (
size_t m = 0; m < m_nel; m++) {
894 for (
size_t ik = 0; ik < m_nsp; ik++) {
896 if (m_mix->
nAtoms(m_species[k],m_element[m]) != 0) {
901 for (
size_t ij = 0; ij < m_nel; ij++) {
902 if (k == m_order[ij]) {
906 if (!ok || m_force) {
914 doublereal MultiPhaseEquil::error()
916 doublereal err, maxerr = 0.0;
919 for (
size_t j = 0; j <
nFree(); j++) {
920 size_t ik = j + m_nel;
924 if (!isStoichPhase(ik) && fabs(moles(ik)) <=
SmallNumber) {
930 if (isStoichPhase(ik) && moles(ik) <= 0.0 &&
931 m_deltaG_RT[j] >= 0.0) {
934 err = fabs(m_deltaG_RT[j]);
943 double MultiPhaseEquil::phaseMoles(index_t iph)
const
948 void MultiPhaseEquil::reportCSV(
const std::string& reportFile)
956 size_t nphase = m_np;
958 FILE* FP = fopen(reportFile.c_str(),
"w");
960 printf(
"Failure to open file\n");
965 vector<double> mf(m_nsp_mix, 1.0);
966 vector<double> fe(m_nsp_mix, 0.0);
968 std::vector<double> VolPM;
969 std::vector<double> activity;
970 std::vector<double> ac;
971 std::vector<double> mu;
972 std::vector<double> mu0;
973 std::vector<double> molalities;
977 for (
size_t iphase = 0; iphase < nphase; iphase++) {
979 ThermoPhase& tref = m_mix->
phase(iphase);
981 VolPM.resize(nSpecies, 0.0);
982 tref.getMoleFractions(&mf[istart]);
983 tref.getPartialMolarVolumes(
DATA_PTR(VolPM));
986 double TMolesPhase = phaseMoles(iphase);
987 double VolPhaseVolumes = 0.0;
988 for (k = 0; k < nSpecies; k++) {
989 VolPhaseVolumes += VolPM[k] * mf[istart + k];
991 VolPhaseVolumes *= TMolesPhase;
992 vol += VolPhaseVolumes;
994 fprintf(FP,
"--------------------- VCS_MULTIPHASE_EQUIL FINAL REPORT"
995 " -----------------------------\n");
996 fprintf(FP,
"Temperature = %11.5g kelvin\n", Temp);
997 fprintf(FP,
"Pressure = %11.5g Pascal\n", pres);
998 fprintf(FP,
"Total Volume = %11.5g m**3\n", vol);
1002 for (
size_t iphase = 0; iphase < nphase; iphase++) {
1005 ThermoPhase& tref = m_mix->
phase(iphase);
1006 ThermoPhase* tp = &tref;
1008 string phaseName = tref.name();
1010 double TMolesPhase = phaseMoles(iphase);
1012 nSpecies = tref.nSpecies();
1013 activity.resize(nSpecies, 0.0);
1014 ac.resize(nSpecies, 0.0);
1016 mu0.resize(nSpecies, 0.0);
1017 mu.resize(nSpecies, 0.0);
1018 VolPM.resize(nSpecies, 0.0);
1019 molalities.resize(nSpecies, 0.0);
1021 int actConvention = tp->activityConvention();
1022 tp->getActivities(
DATA_PTR(activity));
1023 tp->getActivityCoefficients(
DATA_PTR(ac));
1024 tp->getStandardChemPotentials(
DATA_PTR(mu0));
1026 tp->getPartialMolarVolumes(
DATA_PTR(VolPM));
1027 tp->getChemPotentials(
DATA_PTR(mu));
1028 double VolPhaseVolumes = 0.0;
1029 for (k = 0; k < nSpecies; k++) {
1030 VolPhaseVolumes += VolPM[k] * mf[istart + k];
1032 VolPhaseVolumes *= TMolesPhase;
1033 vol += VolPhaseVolumes;
1034 if (actConvention == 1) {
1035 MolalityVPSSTP* mTP =
static_cast<MolalityVPSSTP*
>(tp);
1036 mTP->getMolalities(
DATA_PTR(molalities));
1037 tp->getChemPotentials(
DATA_PTR(mu));
1040 fprintf(FP,
" Name, Phase, PhaseMoles, Mole_Fract, "
1041 "Molalities, ActCoeff, Activity,"
1042 "ChemPot_SS0, ChemPot, mole_num, PMVol, Phase_Volume\n");
1044 fprintf(FP,
" , , (kmol), , "
1046 " (kJ/gmol), (kJ/gmol), (kmol), (m**3/kmol), (m**3)\n");
1048 for (k = 0; k < nSpecies; k++) {
1049 sName = tp->speciesName(k);
1050 fprintf(FP,
"%12s, %11s, %11.3e, %11.3e, %11.3e, %11.3e, %11.3e,"
1051 "%11.3e, %11.3e, %11.3e, %11.3e, %11.3e\n",
1053 phaseName.c_str(), TMolesPhase,
1054 mf[istart + k], molalities[k], ac[k], activity[k],
1055 mu0[k]*1.0E-6, mu[k]*1.0E-6,
1056 mf[istart + k] * TMolesPhase,
1057 VolPM[k], VolPhaseVolumes);
1062 fprintf(FP,
" Name, Phase, PhaseMoles, Mole_Fract, "
1063 "Molalities, ActCoeff, Activity,"
1064 " ChemPotSS0, ChemPot, mole_num, PMVol, Phase_Volume\n");
1066 fprintf(FP,
" , , (kmol), , "
1068 " (kJ/gmol), (kJ/gmol), (kmol), (m**3/kmol), (m**3)\n");
1070 for (k = 0; k < nSpecies; k++) {
1071 molalities[k] = 0.0;
1073 for (k = 0; k < nSpecies; k++) {
1074 sName = tp->speciesName(k);
1075 fprintf(FP,
"%12s, %11s, %11.3e, %11.3e, %11.3e, %11.3e, %11.3e, "
1076 "%11.3e, %11.3e,% 11.3e, %11.3e, %11.3e\n",
1078 phaseName.c_str(), TMolesPhase,
1079 mf[istart + k], molalities[k], ac[k],
1080 activity[k], mu0[k]*1.0E-6, mu[k]*1.0E-6,
1081 mf[istart + k] * TMolesPhase,
1082 VolPM[k], VolPhaseVolumes);
1089 tp->getChemPotentials(&(fe[istart]));
1090 for (k = 0; k < nSpecies; k++) {