15 MultiPhaseEquil::MultiPhaseEquil(
MultiPhase* mix,
bool start,
int loglevel) : m_mix(mix)
29 m_incl_species.resize(m_nsp_mix,1);
30 m_incl_element.resize(m_nel_mix,1);
31 for (m = 0; m < m_nel_mix; m++) {
36 if (enm ==
"E" || enm ==
"e") {
46 m_incl_element[m] = 0;
47 for (k = 0; k < m_nsp_mix; k++) {
48 if (m_mix->
nAtoms(k,m) != 0.0) {
49 m_incl_species[k] = 0;
58 if (m_eloc < m_nel_mix) {
59 m_element.push_back(m_eloc);
63 for (m = 0; m < m_nel_mix; m++) {
64 if (m_incl_element[m] == 1 && m != m_eloc) {
66 m_element.push_back(m);
82 for (k = 0; k < m_nsp_mix; k++) {
86 m_incl_species[k] = 0;
90 +
" is excluded since its thermo properties are \n"
91 "not valid at this temperature, but it has "
92 "non-zero moles in the initial state.");
99 for (k = 0; k < m_nsp_mix; k++) {
100 if (m_incl_species[k] ==1) {
102 m_species.push_back(k);
107 m_work.resize(m_nsp);
108 m_work2.resize(m_nsp);
109 m_work3.resize(m_nsp_mix);
110 m_mu.resize(m_nsp_mix);
113 m_moles.resize(m_nsp);
114 m_lastmoles.resize(m_nsp);
115 m_dxi.resize(
nFree());
119 for (ik = 0; ik < m_nsp; ik++) {
124 m_deltaG_RT.resize(
nFree(), 0.0);
126 m_majorsp.resize(m_nsp);
127 m_sortindex.resize(m_nsp,0);
128 m_lastsort.resize(m_nel);
129 m_solnrxn.resize(
nFree());
130 m_A.
resize(m_nel, m_nsp, 0.0);
132 m_order.resize(std::max(m_nsp, m_nel), 0);
133 for (k = 0; k < m_nsp; k++) {
154 for (k = 0; k < m_nsp; k++) {
155 m_moles[k] += m_work[k];
156 m_lastmoles[k] = m_moles[k];
158 m_dsoln.push_back(1);
160 m_dsoln.push_back(0);
172 doublereal MultiPhaseEquil::equilibrate(
int XY, doublereal err,
173 int maxsteps,
int loglevel)
178 for (i = 0; i < maxsteps; i++) {
185 throw CanteraError(
"MultiPhaseEquil::equilibrate",
186 "no convergence in " +
int2str(maxsteps) +
187 " iterations. Error = " +
fp2str(error()));
193 void MultiPhaseEquil::updateMixMoles()
195 fill(m_work3.begin(), m_work3.end(), 0.0);
197 for (k = 0; k < m_nsp; k++) {
198 m_work3[m_species[k]] = m_moles[k];
205 fill(m_work3.begin(), m_work3.end(), 0.0);
207 for (k = 0; k < m_nsp; k++) {
208 m_work3[m_species[k]] = (m_moles[k] > 0.0 ? m_moles[k] : 0.0);
217 double not_mu = 1.0e12;
224 double delta_xi, dxi_min = 1.0e10;
240 for (j = 0; j <
nFree(); j++) {
243 for (ik = 0; ik < m_nsp; ik++) {
244 dg_rt += mu(ik) * m_N(ik,j);
248 idir = (dg_rt < 0.0 ? 1 : -1);
250 for (ik = 0; ik < m_nsp; ik++) {
259 delta_xi = fabs(0.99*moles(ik)/nu);
262 if (!redo && delta_xi < 1.0e-10 && ik < m_nel) {
265 dxi_min = std::min(dxi_min, delta_xi);
269 for (ik = 0; ik < m_nsp; ik++) {
270 moles(ik) += m_N(ik, j) * idir*dxi_min;
285 if (order.size() != m_nsp) {
286 for (k = 0; k < m_nsp; k++) {
290 for (k = 0; k < m_nsp; k++) {
291 m_order[k] = order[k];
295 size_t nRows = m_nel;
296 size_t nColumns = m_nsp;
300 for (m = 0; m < nRows; m++) {
301 for (k = 0; k < nColumns; k++) {
302 m_A(m, k) = m_mix->
nAtoms(m_species[m_order[k]], m_element[m]);
307 for (m = 0; m < nRows; m++) {
309 bool isZeroRow =
true;
310 for (k = m; k < nColumns; k++) {
311 if (fabs(m_A(m,k)) > sqrt(
Tiny)) {
318 size_t n = nRows - 1;
319 bool foundSwapCandidate =
false;
321 for (k = m; k < nColumns; k++) {
322 if (fabs(m_A(n,k)) > sqrt(
Tiny)) {
323 foundSwapCandidate =
true;
327 if (foundSwapCandidate) {
333 for (k = 0; k < nColumns; k++) {
334 std::swap(m_A(n,k), m_A(m,k));
346 if (m < nColumns && m_A(m,m) == 0.0) {
354 doublereal maxmoles = -999.0;
356 for (k = m+1; k < nColumns; k++) {
357 if (m_A(m,k) != 0.0) {
358 if (fabs(m_moles[m_order[k]]) > maxmoles) {
360 maxmoles = fabs(m_moles[m_order[k]]);
367 for (
size_t n = 0; n < nRows; n++) {
368 std::swap(m_A(n, m), m_A(n, kmax));
372 std::swap(m_order[m], m_order[kmax]);
377 for (k = 0; k < nColumns; k++) {
383 for (
size_t n = m+1; n < m_nel; n++) {
384 fctr = m_A(n,m)/m_A(m,m);
385 for (k = 0; k < m_nsp; k++) {
386 m_A(n,k) -= m_A(m,k)*fctr;
393 for (m = std::min(nRows,nColumns)-1; m > 0; m--) {
394 for (
size_t n = m-1; n !=
npos; n--) {
395 if (m_A(n,m) != 0.0) {
397 for (k = m; k < m_nsp; k++) {
398 m_A(n,k) -= fctr*m_A(m,k);
405 for (
size_t n = 0; n < m_nsp; n++) {
407 for (k = 0; k <
nFree(); k++) {
408 m_N(n, k) = -m_A(n, k + m_nel);
411 for (k = 0; k <
nFree(); k++) {
414 m_N(n, n - m_nel) = 1.0;
419 for (j = 0; j <
nFree(); j++) {
420 m_solnrxn[j] =
false;
421 for (k = 0; k < m_nsp; k++) {
432 copy(x.begin(), x.end(), m_work2.begin());
434 for (k = 0; k < m_nsp; k++) {
435 x[m_order[k]] = m_work2[k];
439 void MultiPhaseEquil::step(doublereal omega,
vector_fp& deltaN,
444 throw CanteraError(
"MultiPhaseEquil::step",
"negative omega");
447 for (ik = 0; ik < m_nel; ik++) {
449 m_lastmoles[k] = m_moles[k];
450 m_moles[k] += omega * deltaN[k];
453 for (ik = m_nel; ik < m_nsp; ik++) {
455 m_lastmoles[k] = m_moles[k];
457 m_moles[k] += omega * deltaN[k];
459 m_moles[k] = fabs(m_moles[k])*std::min(10.0,
460 exp(-m_deltaG_RT[ik - m_nel]));
481 doublereal FCTR = 0.99;
482 const doublereal MAJOR_THRESHOLD = 1.0e-12;
484 doublereal omega = 1.0, omax, omegamax = 1.0;
485 for (ik = 0; ik < m_nsp; ik++) {
489 if (m_moles[k] < MAJOR_THRESHOLD) {
499 if (m_dsoln[k] == 1) {
501 if ((m_moles[k] > MAJOR_THRESHOLD) || (ik < m_nel)) {
502 if (m_moles[k] < MAJOR_THRESHOLD) {
505 omax = m_moles[k]*FCTR/(fabs(m_work[k]) +
Tiny);
506 if (m_work[k] < 0.0 && omax < omegamax) {
508 if (omegamax < 1.0e-5) {
514 m_majorsp[k] =
false;
517 if (m_work[k] < 0.0 && m_moles[k] > 0.0) {
518 omax = -m_moles[k]/m_work[k];
519 if (omax < omegamax) {
521 if (omegamax < 1.0e-5) {
531 step(omegamax, m_work);
535 doublereal not_mu = 1.0e12;
537 doublereal grad1 = 0.0;
538 for (k = 0; k < m_nsp; k++) {
539 grad1 += m_work[k] * m_mu[m_species[k]];
544 omega *= fabs(grad0) / (grad1 + fabs(grad0));
545 for (k = 0; k < m_nsp; k++) {
546 m_moles[k] = m_lastmoles[k];
555 size_t j, k, ik, kc, ip;
556 doublereal stoich, nmoles, csum, term1, fctr, rfctr;
558 doublereal grad = 0.0;
562 doublereal not_mu = 1.0e12;
565 for (j = 0; j <
nFree(); j++) {
568 getStoichVector(j, nu);
571 doublereal dg_rt = 0.0;
572 for (k = 0; k < m_nsp; k++) {
573 dg_rt += m_mu[m_species[k]] * nu[k];
577 m_deltaG_RT[j] = dg_rt;
585 if (m_moles[k] <= 0.0 && dg_rt > 0.0) {
590 }
else if (!m_solnrxn[j]) {
596 for (k = 0; k < m_nel; k++) {
600 csum += stoich*stoich*m_dsoln[kc]/nmoles;
604 kc = m_order[j + m_nel];
606 term1 = m_dsoln[kc]/nmoles;
609 doublereal sum = 0.0, psum;
610 for (ip = 0; ip < m_np; ip++) {
614 for (k = 0; k < m_nsp; k++) {
618 psum += stoich * stoich;
624 rfctr = term1 + csum + sum;
625 if (fabs(rfctr) <
Tiny) {
628 fctr = 1.0/(term1 + csum + sum);
631 dxi[j] = -fctr*dg_rt;
634 for (m = 0; m < m_nel; m++) {
635 if (m_moles[m_order[m]] <= 0.0 && (m_N(m, j)*dxi[j] < 0.0)) {
639 grad += dxi[j]*dg_rt;
645 void MultiPhaseEquil::computeN()
648 std::vector<std::pair<double, size_t> > moleFractions(m_nsp);
649 for (
size_t k = 0; k < m_nsp; k++) {
651 moleFractions[k].first = - m_mix->
speciesMoles(m_species[k]);
652 moleFractions[k].second = k;
654 std::sort(moleFractions.begin(), moleFractions.end());
655 for (
size_t k = 0; k < m_nsp; k++) {
656 m_sortindex[k] = moleFractions[k].second;
660 for (
size_t m = 0; m < m_nel; m++) {
662 for (
size_t ik = 0; ik < m_nsp; ik++) {
664 if (m_mix->
nAtoms(m_species[k],m_element[m]) != 0) {
669 for (
size_t ij = 0; ij < m_nel; ij++) {
670 if (k == m_order[ij]) {
674 if (!ok || m_force) {
682 doublereal MultiPhaseEquil::error()
684 doublereal err, maxerr = 0.0;
687 for (
size_t j = 0; j <
nFree(); j++) {
688 size_t ik = j + m_nel;
692 if (!isStoichPhase(ik) && fabs(moles(ik)) <=
SmallNumber) {
698 if (isStoichPhase(ik) && moles(ik) <= 0.0 &&
699 m_deltaG_RT[j] >= 0.0) {
702 err = fabs(m_deltaG_RT[j]);
704 maxerr = std::max(maxerr, err);
709 double MultiPhaseEquil::phaseMoles(
size_t iph)
const
714 void MultiPhaseEquil::reportCSV(
const std::string& reportFile)
722 size_t nphase = m_np;
724 FILE* FP = fopen(reportFile.c_str(),
"w");
726 throw CanteraError(
"MultiPhaseEquil::reportCSV",
"Failure to open file");
730 vector<double> mf(m_nsp_mix, 1.0);
731 vector<double> fe(m_nsp_mix, 0.0);
733 std::vector<double> VolPM;
734 std::vector<double> activity;
735 std::vector<double> ac;
736 std::vector<double> mu;
737 std::vector<double> mu0;
738 std::vector<double> molalities;
742 for (
size_t iphase = 0; iphase < nphase; iphase++) {
744 ThermoPhase& tref = m_mix->
phase(iphase);
746 VolPM.resize(nSpecies, 0.0);
747 tref.getMoleFractions(&mf[istart]);
748 tref.getPartialMolarVolumes(
DATA_PTR(VolPM));
750 double TMolesPhase = phaseMoles(iphase);
751 double VolPhaseVolumes = 0.0;
752 for (k = 0; k < nSpecies; k++) {
753 VolPhaseVolumes += VolPM[k] * mf[istart + k];
755 VolPhaseVolumes *= TMolesPhase;
756 vol += VolPhaseVolumes;
758 fprintf(FP,
"--------------------- VCS_MULTIPHASE_EQUIL FINAL REPORT"
759 " -----------------------------\n");
760 fprintf(FP,
"Temperature = %11.5g kelvin\n", Temp);
761 fprintf(FP,
"Pressure = %11.5g Pascal\n", pres);
762 fprintf(FP,
"Total Volume = %11.5g m**3\n", vol);
764 for (
size_t iphase = 0; iphase < nphase; iphase++) {
767 ThermoPhase& tref = m_mix->
phase(iphase);
768 ThermoPhase* tp = &tref;
770 string phaseName = tref.name();
771 double TMolesPhase = phaseMoles(iphase);
772 nSpecies = tref.nSpecies();
773 activity.resize(nSpecies, 0.0);
774 ac.resize(nSpecies, 0.0);
776 mu0.resize(nSpecies, 0.0);
777 mu.resize(nSpecies, 0.0);
778 VolPM.resize(nSpecies, 0.0);
779 molalities.resize(nSpecies, 0.0);
781 int actConvention = tp->activityConvention();
782 tp->getActivities(
DATA_PTR(activity));
783 tp->getActivityCoefficients(
DATA_PTR(ac));
784 tp->getStandardChemPotentials(
DATA_PTR(mu0));
786 tp->getPartialMolarVolumes(
DATA_PTR(VolPM));
787 tp->getChemPotentials(
DATA_PTR(mu));
788 double VolPhaseVolumes = 0.0;
789 for (k = 0; k < nSpecies; k++) {
790 VolPhaseVolumes += VolPM[k] * mf[istart + k];
792 VolPhaseVolumes *= TMolesPhase;
793 vol += VolPhaseVolumes;
794 if (actConvention == 1) {
795 MolalityVPSSTP* mTP =
static_cast<MolalityVPSSTP*
>(tp);
796 mTP->getMolalities(
DATA_PTR(molalities));
797 tp->getChemPotentials(
DATA_PTR(mu));
800 fprintf(FP,
" Name, Phase, PhaseMoles, Mole_Fract, "
801 "Molalities, ActCoeff, Activity,"
802 "ChemPot_SS0, ChemPot, mole_num, PMVol, Phase_Volume\n");
804 fprintf(FP,
" , , (kmol), , "
806 " (kJ/gmol), (kJ/gmol), (kmol), (m**3/kmol), (m**3)\n");
808 for (k = 0; k < nSpecies; k++) {
809 sName = tp->speciesName(k);
810 fprintf(FP,
"%12s, %11s, %11.3e, %11.3e, %11.3e, %11.3e, %11.3e,"
811 "%11.3e, %11.3e, %11.3e, %11.3e, %11.3e\n",
813 phaseName.c_str(), TMolesPhase,
814 mf[istart + k], molalities[k], ac[k], activity[k],
815 mu0[k]*1.0E-6, mu[k]*1.0E-6,
816 mf[istart + k] * TMolesPhase,
817 VolPM[k], VolPhaseVolumes);
822 fprintf(FP,
" Name, Phase, PhaseMoles, Mole_Fract, "
823 "Molalities, ActCoeff, Activity,"
824 " ChemPotSS0, ChemPot, mole_num, PMVol, Phase_Volume\n");
826 fprintf(FP,
" , , (kmol), , "
828 " (kJ/gmol), (kJ/gmol), (kmol), (m**3/kmol), (m**3)\n");
830 for (k = 0; k < nSpecies; k++) {
833 for (k = 0; k < nSpecies; k++) {
834 sName = tp->speciesName(k);
835 fprintf(FP,
"%12s, %11s, %11.3e, %11.3e, %11.3e, %11.3e, %11.3e, "
836 "%11.3e, %11.3e,% 11.3e, %11.3e, %11.3e\n",
838 phaseName.c_str(), TMolesPhase,
839 mf[istart + k], molalities[k], ac[k],
840 activity[k], mu0[k]*1.0E-6, mu[k]*1.0E-6,
841 mf[istart + k] * TMolesPhase,
842 VolPM[k], VolPhaseVolumes);
std::string int2str(const int n, const std::string &fmt)
Convert an int to a string using a format converter.
bool tempOK(size_t p) const
Return true if the phase p has valid thermo data for the current temperature.
doublereal speciesMoles(size_t kGlob) const
Returns the moles of global species k. units = kmol.
doublereal computeReactionSteps(vector_fp &dxi)
Compute the change in extent of reaction for each reaction.
const size_t npos
index returned by functions to indicate "no position"
std::string speciesName(const size_t kGlob) const
Name of species with global index kGlob.
size_t nFree() const
Number of degrees of freedom.
doublereal phaseMoles(const size_t n) const
Return the number of moles in phase n.
doublereal temperature() const
Temperature [K].
void getComponents(const std::vector< size_t > &order)
This method finds a set of component species and a complete set of formation reactions for the non-co...
void getMoleFractions(doublereal *const x) const
Get the species mole fraction vector.
Header for intermediate ThermoPhase object for phases which employ molality based activity coefficien...
int setInitialMoles(int loglevel=0)
Estimate the initial mole numbers.
Base class for a phase with thermodynamic properties.
size_t speciesPhaseIndex(const size_t kGlob) const
Returns the phase index of the Kth "global" species.
A class for multiphase mixtures.
doublereal stepComposition(int loglevel=0)
Take one step in composition, given the gradient of G at the starting point, and a vector of reaction...
void resize(size_t n, size_t m, doublereal v=0.0)
Resize the matrix.
void unsort(vector_fp &x)
Re-arrange a vector of species properties in sorted form (components first) into unsorted, sequential form.
std::string elementName(size_t m) const
Returns the name of the global element m.
void getValidChemPotentials(doublereal not_mu, doublereal *mu, bool standard=false) const
Returns a vector of Valid chemical potentials.
std::string fp2str(const double x, const std::string &fmt)
Convert a double into a c++ string.
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.
Base class for exceptions thrown by Cantera classes.
doublereal nAtoms(const size_t kGlob, const size_t mGlob) const
Returns the Number of atoms of global element mGlob in global species kGlob.
size_t nSpecies() const
Returns the number of species in the phase.
const doublereal SmallNumber
smallest number to compare to zero.
std::vector< double > vector_fp
Turn on the use of stl vectors for the basic array type within cantera Vector of doubles.
const doublereal Tiny
Small number to compare differences of mole fractions against.
const doublereal GasConstant
Universal Gas Constant. [J/kmol/K].
doublereal elementMoles(size_t m) const
Total moles of global element m, summed over all phases.
Contains declarations for string manipulation functions within Cantera.
#define DATA_PTR(vec)
Creates a pointer to the start of the raw data for a vector.
bool solutionSpecies(size_t kGlob) const
Return true is species kGlob is a species in a multicomponent solution phase.
thermo_t & phase(size_t n)
Return a reference to phase n.
doublereal pressure() const
Pressure [Pa].
size_t nSpecies() const
Number of species, summed over all phases.
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...
void finish()
Clean up the composition.
void setMoles(const doublereal *n)
Sets all of the global species mole numbers.
size_t nElements() const
Number of elements.
size_t nPhases() const
Number of phases.