21 MultiPhase::MultiPhase() :
35 for (
size_t n = 0; n < mix.
nPhases(); n++) {
43 for (
size_t n = 0; n < phases.size(); n++) {
53 "phases cannot be added after init() has been called.");
57 throw CanteraError(
"MultiPhase::addPhase",
"Phase '{}'' is not " 58 "compatible with MultiPhase equilibrium solver", p->
name());
75 for (
size_t m = 0; m < p->
nElements(); m++) {
87 if (ename ==
"E" || ename ==
"e") {
129 for (
size_t m = 0; m <
m_nel; m++) {
132 for (
size_t ip = 0; ip <
nPhases(); ip++) {
136 for (
size_t kp = 0; kp < nsp; kp++) {
137 if (mlocal !=
npos) {
192 doublereal sum = 0.0;
193 for (
size_t i = 0; i <
nPhases(); i++) {
194 double phasesum = 0.0;
195 size_t nsp =
m_phase[i]->nSpecies();
196 for (
size_t ik = 0; ik < nsp; ik++) {
207 doublereal sum = 0.0;
208 for (
size_t i = 0; i <
nPhases(); i++) {
232 doublereal phasesum = 0.0;
233 size_t nsp =
m_phase[p]->nSpecies();
234 for (
size_t ik = 0; ik < nsp; ik++) {
238 return Faraday*phasesum*
m_moles[p];
245 for (
size_t i = 0; i <
nPhases(); i++) {
246 m_phase[i]->getChemPotentials(mu + loc);
252 doublereal* mu,
bool standard)
const 257 for (
size_t i = 0; i <
nPhases(); i++) {
260 m_phase[i]->getChemPotentials(mu + loc);
262 m_phase[i]->getStandardChemPotentials(mu + loc);
282 doublereal sum = 0.0;
284 for (
size_t i = 0; i <
nPhases(); i++) {
294 doublereal sum = 0.0;
296 for (
size_t i = 0; i <
nPhases(); i++) {
306 doublereal sum = 0.0;
308 for (
size_t i = 0; i <
nPhases(); i++) {
318 doublereal sum = 0.0;
320 for (
size_t i = 0; i <
nPhases(); i++) {
330 doublereal sum = 0.0;
332 for (
size_t i = 0; i <
nPhases(); i++) {
348 for (
size_t k = 0; k < p->
nSpecies(); k++) {
357 for (
size_t k = 0; k < kk; k++) {
374 doublereal* dtmp = molNum;
375 for (
size_t ip = 0; ip <
nPhases(); ip++) {
376 doublereal phasemoles =
m_moles[ip];
379 for (
size_t ik = 0; ik < nsp; ik++) {
380 *(dtmp++) *= phasemoles;
392 for (
size_t ip = 0; ip <
nPhases(); ip++) {
395 double phasemoles = 0.0;
396 for (
size_t ik = 0; ik < nsp; ik++) {
402 if (phasemoles > 0.0) {
419 tmpMoles[indexS] += addedMoles;
420 tmpMoles[indexS] = std::max(tmpMoles[indexS], 0.0);
445 for (
size_t eGlobal = 0; eGlobal <
m_nel; eGlobal++) {
454 for (
size_t eGlobal = 0; eGlobal <
m_nel; eGlobal++) {
457 for (
size_t ip = 0; ip <
nPhases(); ip++) {
460 doublereal phasemoles =
m_moles[ip];
461 for (
size_t ik = 0; ik < nspPhase; ik++) {
462 size_t kGlobal = loc + ik;
464 for (
size_t eGlobal = 0; eGlobal <
m_nel; eGlobal++) {
475 for (
size_t i = 0; i <
nPhases(); i++) {
476 double vol = 1.0/
m_phase[i]->molarDensity();
483 int maxsteps,
int maxiter,
487 doublereal dta = 0.0;
495 return e.equilibrate(XY, err, maxsteps, loglevel);
496 }
else if (XY == HP) {
499 double Thigh = 2.0*
m_Tmax;
501 for (
int n = 0; n < maxiter; n++) {
509 e.equilibrate(TP, err, maxsteps, loglevel);
529 double cpb = (Hhigh - Hlow)/(Thigh - Tlow);
530 dt = (h0 - hnow)/cpb;
532 double dtmax = 0.5*fabs(Thigh - Tlow);
537 double tnew = sqrt(Tlow*Thigh);
541 double herr = fabs((h0 - hnow)/h0);
546 double tnew =
m_temp + dt;
562 double tnew = 0.5*(
m_temp + Thigh);
563 if (fabs(tnew -
m_temp) < 1.0) {
570 throw CanteraError(
"MultiPhase::equilibrate_MultiPhaseEquil",
571 "No convergence for T");
572 }
else if (XY == SP) {
575 double Thigh = 1.0e6;
576 for (
int n = 0; n < maxiter; n++) {
580 e.equilibrate(TP, err, maxsteps, loglevel);
583 Tlow = std::max(Tlow,
m_temp);
585 Thigh = std::min(Thigh,
m_temp);
587 double dt = (s0 - snow)*
m_temp/
cp();
588 double dtmax = 0.5*fabs(Thigh - Tlow);
589 dtmax = (dtmax > 500.0 ? 500.0 : dtmax);
597 double tnew =
m_temp + dt;
609 double tnew = 0.5*(
m_temp + Thigh);
614 throw CanteraError(
"MultiPhase::equilibrate_MultiPhaseEquil",
615 "No convergence for T");
616 }
else if (XY == TV) {
619 for (
int n = 0; n < maxiter; n++) {
624 e.equilibrate(TP, err, maxsteps, loglevel);
626 double verr = fabs((v0 - vnow)/v0);
633 double dVdP = (
volume() - vnow)/(0.01*pnow);
637 throw CanteraError(
"MultiPhase::equilibrate_MultiPhaseEquil",
644 double rtol,
int max_steps,
int max_iter,
645 int estimate_equil,
int log_level)
651 double initial_T =
m_temp;
654 if (solver ==
"auto" || solver ==
"vcs") {
656 debuglog(
"Trying VCS equilibrium solver\n", log_level);
658 int ret = eqsolve.
equilibrate(ixy, estimate_equil, log_level-1,
662 "VCS solver failed. Return code: {}", ret);
664 debuglog(
"VCS solver succeeded\n", log_level);
666 }
catch (std::exception& err) {
667 debuglog(
"VCS solver failed.\n", log_level);
674 if (solver ==
"auto") {
681 if (solver ==
"auto" || solver ==
"gibbs") {
683 debuglog(
"Trying MultiPhaseEquil (Gibbs) equilibrium solver\n",
687 debuglog(
"MultiPhaseEquil solver succeeded\n", log_level);
689 }
catch (std::exception& err) {
690 debuglog(
"MultiPhaseEquil solver failed.\n", log_level);
701 if (solver !=
"auto") {
703 "Invalid solver specified: '" + solver +
"'");
737 for (
size_t e = 0; e <
m_nel; e++) {
782 for (
int iph = 0; iph < (int)
nPhases(); iph++) {
783 if (
m_phase[iph]->
id() == pName) {
818 for (
size_t ip = 0; ip <
nPhases(); ip++) {
829 for (
size_t p = 0; p <
nPhases(); p++) {
std::map< std::string, doublereal > compositionMap
Map connecting a string name with a double.
bool tempOK(size_t p) const
Return true if the phase p has valid thermo data for the current temperature.
size_t nElements() const
Number of elements.
std::vector< size_t > m_spphase
Mapping between the global species number and the phase ID.
doublereal m_temp
Current value of the temperature (kelvin)
void checkPhaseArraySize(size_t mm) const
Check that an array size is at least nPhases() Throws an exception if mm is less than nPhases()...
doublereal cp() const
Heat capacity at constant pressure [J/K].
size_t nPhases() const
Number of phases.
doublereal temperature() const
Temperature (K).
void addPhases(std::vector< ThermoPhase *> &phases, const vector_fp &phaseMoles)
Add a vector of phases to the mixture.
bool m_init
True if the init() routine has been called, and the MultiPhase frozen.
doublereal speciesMoles(size_t kGlob) const
Returns the moles of global species k. units = kmol.
doublereal charge() const
Total charge summed over all phases (Coulombs).
virtual void setState_TPX(doublereal t, doublereal p, const doublereal *x)
Set the temperature (K), pressure (Pa), and mole fractions.
void setPhaseMoles(const size_t n, const doublereal moles)
Set the number of moles of phase with index n.
void getMoleFractions(doublereal *const x) const
Returns the global Species mole fractions.
std::map< std::string, size_t > m_enamemap
Returns the global element index, given the element string name.
const size_t npos
index returned by functions to indicate "no position"
void checkSpeciesIndex(size_t k) const
Check that the specified species index is in range.
size_t elementIndex(const std::string &name) const
Return the index of element named 'name'.
vector_int m_atomicNumber
Atomic number of each global element.
void getElemAbundances(doublereal *elemAbundances) const
Retrieves a vector of element abundances.
size_t nSpecies() const
Number of species, summed over all phases.
void addPhase(ThermoPhase *p, doublereal moles)
Add a phase to the mixture.
size_t nSpecies() const
Returns the number of species in the phase.
doublereal volume() const
The total mixture volume [m^3].
virtual doublereal minTemp(size_t k=npos) const
Minimum temperature for which the thermodynamic data for the species or phase are valid...
doublereal pressure() const
Pressure [Pa].
void equilibrate(const std::string &XY, const std::string &solver="auto", double rtol=1e-9, int max_steps=50000, int max_iter=100, int estimate_equil=0, int log_level=0)
Equilibrate a MultiPhase object.
Headers for the MultiPhase object that is used to set up multiphase equilibrium problems (see Equilfu...
const doublereal Undef
Fairly random number to be used to initialize variables against to see if they are subsequently defin...
vector_fp m_moles
Vector of the number of moles in each phase.
Base class for a phase with thermodynamic properties.
int atomicNumber(size_t m) const
Atomic number of element m.
virtual bool compatibleWithMultiPhase() const
Indicates whether this phase type can be used with class MultiPhase for equilibrium calculations...
std::vector< std::string > m_enames
String names of the global elements.
void addSpeciesMoles(const int indexS, const doublereal addedMoles)
Adds moles of a certain species to the mixture.
vector_fp m_moleFractions
Locally stored vector of mole fractions of all species comprising the MultiPhase object.
Cantera's Interface to the Multiphase chemical equilibrium solver.
size_t elementIndex(const std::string &name) const
Returns the index of the element with name name.
void uploadMoleFractionsFromPhases()
Update the locally-stored composition within this object to match the current compositions of the pha...
A class for multiphase mixtures.
void resize(size_t n, size_t m, doublereal v=0.0)
Resize the matrix.
std::vector< bool > m_temp_OK
Vector of bools indicating whether temperatures are ok for phases.
void setState_TPMoles(const doublereal T, const doublereal Pres, const doublereal *Moles)
Set the state of the underlying ThermoPhase objects in one call.
size_t m_eloc
Global ID of the element corresponding to the electronic charge.
doublereal moleFraction(const size_t kGlob) const
Returns the mole fraction of global species k.
std::string speciesName(size_t k) const
Name of the species with index k.
DenseMatrix m_atoms
Global Stoichiometric Coefficient array.
void calcElemAbundances() const
Calculate the element abundance vector.
std::vector< size_t > m_spstart
Vector of ints containing of first species index in the global list of species for each phase...
void setPhaseMoleFractions(const size_t n, const doublereal *const x)
Set the Mole fractions of the nth phase.
size_t m_nsp
Number of distinct species in all of the 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...
doublereal maxTemp() const
Maximum temperature for which all solution phases have valid thermo data.
size_t m_nel
Number of distinct elements in all of the phases.
doublereal nAtoms(const size_t kGlob, const size_t mGlob) const
Returns the Number of atoms of global element mGlob in global species kGlob.
const U & getValue(const std::map< T, U > &m, const T &key, const U &default_val)
Const accessor for a value in a std::map.
void setTemperature(const doublereal T)
Set the temperature [K].
Base class for exceptions thrown by Cantera classes.
void checkPhaseIndex(size_t m) const
Check that the specified phase index is in range Throws an exception if m is greater than nPhases() ...
void checkSpeciesArraySize(size_t kk) const
Check that an array size is at least nSpecies().
double equilibrate_MultiPhaseEquil(int XY, doublereal err, int maxsteps, int maxiter, int loglevel)
Set the mixture to a state of chemical equilibrium using the MultiPhaseEquil solver.
int equilibrate(int XY, int estimateEquil=0, int printLvl=0, doublereal err=1.0e-6, int maxsteps=VCS_MAXSTEPS, int loglevel=-99)
Equilibrate the solution using the current element abundances stored in the MultiPhase object...
void getChemPotentials(doublereal *mu) const
Returns a vector of Chemical potentials.
vector_fp m_elemAbundances
Vector of element abundances.
size_t speciesPhaseIndex(const size_t kGlob) const
Returns the phase index of the Kth "global" species.
void getMoles(doublereal *molNum) const
Get the mole numbers of all species in the multiphase object.
std::vector< ThermoPhase * > m_phase
Vector of the ThermoPhase pointers.
doublereal phaseCharge(size_t p) const
Charge (Coulombs) of phase with index p.
void getMoleFractions(doublereal *const x) const
Get the species mole fraction vector.
void setMolesByName(const compositionMap &xMap)
Set the number of moles of species in the mixture.
virtual doublereal pressure() const
Return the thermodynamic pressure (Pa).
void setPressure(doublereal P)
Set the pressure [Pa].
doublereal elementMoles(size_t m) const
Total moles of global element m, summed over all phases.
void debuglog(const std::string &msg, int loglevel)
Write a message to the log only if loglevel > 0.
Interface class for the vcsnonlinear solver.
int phaseIndex(const std::string &pName) const
Returns the index, given the phase name.
void setState_TP(const doublereal T, const doublereal Pres)
Set the state of the underlying ThermoPhase objects in one call.
compositionMap parseCompString(const std::string &ss, const std::vector< std::string > &names)
Parse a composition string into a map consisting of individual key:composition pairs.
doublereal entropy() const
The entropy of the mixture [J/K].
void updatePhases() const
Set the states of the phase objects to the locally-stored state within this MultiPhase object...
std::string id() const
Return the string id for the phase.
std::vector< double > vector_fp
Turn on the use of stl vectors for the basic array type within cantera Vector of doubles.
doublereal minTemp() const
Minimum temperature for which all solution phases have valid thermo data.
bool solutionSpecies(size_t kGlob) const
Return true is species kGlob is a species in a multicomponent solution phase.
doublereal IntEnergy() const
The internal energy of the mixture [J].
void init()
Process phases and build atomic composition array.
std::string name() const
Return the name of the phase.
Contains declarations for string manipulation functions within Cantera.
void getValidChemPotentials(doublereal not_mu, doublereal *mu, bool standard=false) const
Returns a vector of Valid chemical potentials.
std::vector< std::string > m_snames
Vector of species names in the problem.
std::string elementName(size_t m) const
Returns the name of the global element m.
std::string phaseName(const size_t iph) const
Returns the name of the n'th phase.
An array index is out of range.
void checkElementArraySize(size_t mm) const
Check that an array size is at least nElements().
int _equilflag(const char *xy)
map property strings to integers
thermo_t & phase(size_t n)
Return a reference to phase n.
doublereal m_press
Current value of the pressure (Pa)
doublereal phaseMoles(const size_t n) const
Return the number of moles in phase n.
doublereal gibbs() const
The Gibbs function of the mixture [J].
Namespace for the Cantera kernel.
virtual doublereal maxTemp(size_t k=npos) const
Maximum temperature for which the thermodynamic data for the species are valid.
doublereal m_Tmax
Minimum temperature for which thermo parameterizations are valid.
doublereal nAtoms(size_t k, size_t m) const
Number of atoms of element m in species k.
doublereal m_Tmin
Minimum temperature for which thermo parameterizations are valid.
void checkElementIndex(size_t m) const
Check that the specified element index is in range.
void setMoles(const doublereal *n)
Sets all of the global species mole numbers.
std::string elementName(size_t m) const
Name of the element with index m.
doublereal enthalpy() const
The enthalpy of the mixture [J].
std::string speciesName(const size_t kGlob) const
Name of species with global index kGlob.