32 for (
size_t k = 0; k <
nSpecies(); k++) {
52 return Units(1.0, 0, -
static_cast<double>(
nDim()), 0, 0, 0, 1);
63 for (
size_t k = 0; k <
nSpecies(); k++) {
71 for (
size_t k = 0; k <
m_kk; k++) {
72 lnac[k] = std::log(lnac[k]);
80 for (
size_t k = 0; k <
m_kk; k++) {
147 AnyMap state = input_state;
150 if (state.
hasKey(
"mass-fractions")) {
151 state[
"Y"] = state[
"mass-fractions"];
152 state.
erase(
"mass-fractions");
154 if (state.
hasKey(
"mole-fractions")) {
155 state[
"X"] = state[
"mole-fractions"];
156 state.
erase(
"mole-fractions");
158 if (state.
hasKey(
"temperature")) {
159 state[
"T"] = state[
"temperature"];
161 if (state.
hasKey(
"pressure")) {
162 state[
"P"] = state[
"pressure"];
164 if (state.
hasKey(
"enthalpy")) {
165 state[
"H"] = state[
"enthalpy"];
167 if (state.
hasKey(
"int-energy")) {
168 state[
"U"] = state[
"int-energy"];
170 if (state.
hasKey(
"internal-energy")) {
171 state[
"U"] = state[
"internal-energy"];
173 if (state.
hasKey(
"specific-volume")) {
174 state[
"V"] = state[
"specific-volume"];
176 if (state.
hasKey(
"entropy")) {
177 state[
"S"] = state[
"entropy"];
179 if (state.
hasKey(
"density")) {
180 state[
"D"] = state[
"density"];
182 if (state.
hasKey(
"vapor-fraction")) {
183 state[
"Q"] = state[
"vapor-fraction"];
188 if (state[
"X"].is<string>()) {
194 }
else if (state.
hasKey(
"Y")) {
195 if (state[
"Y"].is<string>()) {
203 if (state.
size() == 0) {
206 double T = state.
convert(
"T",
"K");
207 double P = state.
convert(
"P",
"Pa");
243 }
else if (state.
hasKey(
"T")) {
245 }
else if (state.
hasKey(
"P")) {
249 "'state' did not specify a recognized set of properties.\n"
250 "Keys provided were: {}", input_state.
keys_str());
263 double rtol,
bool doUV)
273 "Input specific volume is too small or negative. v = {}", v);
279 "Input pressure is too small or negative. p = {}", p);
292 }
else if (Tnew < Tmin) {
306 bool ignoreBounds =
false;
309 bool unstablePhase =
false;
312 double Tunstable = -1.0;
313 bool unstablePhaseNew =
false;
316 for (
int n = 0; n < 500; n++) {
321 unstablePhase =
true;
325 dt =
clip((Htarget - Hold)/cpd, -100.0, 100.0);
332 if ((dt > 0.0 && unstablePhase) || (dt <= 0.0 && !unstablePhase)) {
333 if (Hbot < Htarget && Tnew < (0.75 * Tbot + 0.25 * Told)) {
334 dt = 0.75 * (Tbot - Told);
337 }
else if (Htop > Htarget && Tnew > (0.75 * Ttop + 0.25 * Told)) {
338 dt = 0.75 * (Ttop - Told);
343 if (Tnew > Tmax && !ignoreBounds) {
346 if (Hmax >= Htarget) {
347 if (Htop < Htarget) {
356 if (Tnew < Tmin && !ignoreBounds) {
359 if (Hmin <= Htarget) {
360 if (Hbot > Htarget) {
373 for (
int its = 0; its < 10; its++) {
375 if (Tnew < Told / 3.0) {
377 dt = -2.0 * Told / 3.0;
388 unstablePhaseNew =
true;
391 unstablePhaseNew =
false;
394 if (unstablePhase ==
false && unstablePhaseNew ==
true) {
399 if (Hnew == Htarget) {
401 }
else if (Hnew > Htarget && (Htop < Htarget || Hnew < Htop)) {
404 }
else if (Hnew < Htarget && (Hbot > Htarget || Hnew > Hbot)) {
409 double Herr = Htarget - Hnew;
410 double acpd = std::max(fabs(cpd), 1.0E-5);
411 double denom = std::max(fabs(Htarget), acpd * Tnew);
412 double HConvErr = fabs((Herr)/denom);
413 if (HConvErr < rtol || fabs(dt/Tnew) < rtol) {
421 string ErrString =
"No convergence in 500 iterations\n";
423 ErrString += fmt::format(
424 "\tTarget Internal Energy = {}\n"
425 "\tCurrent Specific Volume = {}\n"
426 "\tStarting Temperature = {}\n"
427 "\tCurrent Temperature = {}\n"
428 "\tCurrent Internal Energy = {}\n"
429 "\tCurrent Delta T = {}\n",
430 Htarget, v, Tinit, Tnew, Hnew, dt);
432 ErrString += fmt::format(
433 "\tTarget Enthalpy = {}\n"
434 "\tCurrent Pressure = {}\n"
435 "\tStarting Temperature = {}\n"
436 "\tCurrent Temperature = {}\n"
437 "\tCurrent Enthalpy = {}\n"
438 "\tCurrent Delta T = {}\n",
439 Htarget, p, Tinit, Tnew, Hnew, dt);
442 ErrString += fmt::format(
443 "\t - The phase became unstable (Cp < 0) T_unstable_last = {}\n",
447 throw CanteraError(
"ThermoPhase::setState_HPorUV (UV)", ErrString);
449 throw CanteraError(
"ThermoPhase::setState_HPorUV (HP)", ErrString);
465 double rtol,
bool doSV)
473 "Input specific volume is too small or negative. v = {}", v);
479 "Input pressure is too small or negative. p = {}", p);
492 }
else if (Tnew < Tmin) {
506 bool ignoreBounds =
false;
509 bool unstablePhase =
false;
510 double Tunstable = -1.0;
511 bool unstablePhaseNew =
false;
514 for (
int n = 0; n < 500; n++) {
519 unstablePhase =
true;
523 dt =
clip((Starget - Sold)*Told/cpd, -100.0, 100.0);
527 if ((dt > 0.0 && unstablePhase) || (dt <= 0.0 && !unstablePhase)) {
528 if (Sbot < Starget && Tnew < Tbot) {
529 dt = 0.75 * (Tbot - Told);
532 }
else if (Stop > Starget && Tnew > Ttop) {
533 dt = 0.75 * (Ttop - Told);
538 if (Tnew > Tmax && !ignoreBounds) {
541 if (Smax >= Starget) {
542 if (Stop < Starget) {
550 }
else if (Tnew < Tmin && !ignoreBounds) {
553 if (Smin <= Starget) {
554 if (Sbot > Starget) {
567 for (
int its = 0; its < 10; its++) {
573 unstablePhaseNew =
true;
576 unstablePhaseNew =
false;
579 if (unstablePhase ==
false && unstablePhaseNew ==
true) {
584 if (Snew == Starget) {
586 }
else if (Snew > Starget && (Stop < Starget || Snew < Stop)) {
589 }
else if (Snew < Starget && (Sbot > Starget || Snew > Sbot)) {
594 double Serr = Starget - Snew;
595 double acpd = std::max(fabs(cpd), 1.0E-5);
596 double denom = std::max(fabs(Starget), acpd * Tnew);
597 double SConvErr = fabs((Serr * Tnew)/denom);
598 if (SConvErr < rtol || fabs(dt/Tnew) < rtol) {
606 string ErrString =
"No convergence in 500 iterations\n";
608 ErrString += fmt::format(
609 "\tTarget Entropy = {}\n"
610 "\tCurrent Specific Volume = {}\n"
611 "\tStarting Temperature = {}\n"
612 "\tCurrent Temperature = {}\n"
613 "\tCurrent Entropy = {}\n"
614 "\tCurrent Delta T = {}\n",
615 Starget, v, Tinit, Tnew, Snew, dt);
617 ErrString += fmt::format(
618 "\tTarget Entropy = {}\n"
619 "\tCurrent Pressure = {}\n"
620 "\tStarting Temperature = {}\n"
621 "\tCurrent Temperature = {}\n"
622 "\tCurrent Entropy = {}\n"
623 "\tCurrent Delta T = {}\n",
624 Starget, p, Tinit, Tnew, Snew, dt);
627 ErrString += fmt::format(
"\t - The phase became unstable (Cp < 0) T_unstable_last = {}\n",
631 throw CanteraError(
"ThermoPhase::setState_SPorSV (SV)", ErrString);
633 throw CanteraError(
"ThermoPhase::setState_SPorSV (SP)", ErrString);
646 for (
size_t k = 0; k !=
m_kk; ++k) {
650 o2req += x *
nAtoms(k, iC);
653 o2req += x *
nAtoms(k, iS);
656 o2req += x * 0.25 *
nAtoms(k, iH);
661 "No composition specified");
671 for (
size_t k = 0; k !=
m_kk; ++k) {
677 "No composition specified");
679 return 0.5 * o2pres / sum;
695 parseCompString(fuelComp.find(
":") != string::npos ? fuelComp : fuelComp+
":1.0"),
696 parseCompString(oxComp.find(
":") != string::npos ? oxComp : oxComp+
":1.0"),
701 const double* oxComp,
704 vector<double> fuel, ox;
705 if (basis == ThermoBasis::molar) {
710 fuelComp = fuel.data();
717 if (o2_required_fuel < 0.0 || o2_required_ox > 0.0) {
719 "Fuel composition contains too much oxygen or "
720 "oxidizer contains not enough oxygen. "
721 "Fuel and oxidizer composition mixed up?");
724 if (o2_required_ox == 0.0) {
725 return std::numeric_limits<double>::infinity();
728 return o2_required_fuel / (-o2_required_ox);
736 "Equivalence ratio phi must be >= 0");
741 vector<double> fuel, ox;
742 if (basis == ThermoBasis::molar) {
747 fuelComp = fuel.data();
753 double sum_f = std::accumulate(fuelComp, fuelComp+
m_kk, 0.0);
754 double sum_o = std::accumulate(oxComp, oxComp+
m_kk, 0.0);
756 vector<double> y(
m_kk);
757 for (
size_t k = 0; k !=
m_kk; ++k) {
758 y[k] = phi * fuelComp[k]/sum_f + AFR_st * oxComp[k]/sum_o;
769 parseCompString(fuelComp.find(
":") != string::npos ? fuelComp : fuelComp+
":1.0"),
770 parseCompString(oxComp.find(
":") != string::npos ? oxComp : oxComp+
":1.0"),
787 if (o2_present == 0.0) {
788 return std::numeric_limits<double>::infinity();
791 return o2_required / o2_present;
807 parseCompString(fuelComp.find(
":") != string::npos ? fuelComp : fuelComp+
":1.0"),
808 parseCompString(oxComp.find(
":") != string::npos ? oxComp : oxComp+
":1.0"),
813 const double* oxComp,
823 return std::numeric_limits<double>::infinity();
826 vector<double> fuel, ox;
827 if (basis == ThermoBasis::molar) {
832 fuelComp = fuel.data();
838 return std::max(Z / (1.0 - Z) * AFR_st, 0.0);
853 parseCompString(fuelComp.find(
":") != string::npos ? fuelComp : fuelComp+
":1.0"),
854 parseCompString(oxComp.find(
":") != string::npos ? oxComp : oxComp+
":1.0"),
861 if (mixFrac < 0.0 || mixFrac > 1.0) {
863 "Mixture fraction must be between 0 and 1");
866 vector<double> fuel, ox;
867 if (basis == ThermoBasis::molar) {
872 fuelComp = fuel.data();
876 double sum_yf = std::accumulate(fuelComp, fuelComp+
m_kk, 0.0);
877 double sum_yo = std::accumulate(oxComp, oxComp+
m_kk, 0.0);
879 if (sum_yf == 0.0 || sum_yo == 0.0) {
881 "No fuel and/or oxidizer composition specified");
886 vector<double> y(
m_kk);
888 for (
size_t k = 0; k !=
m_kk; ++k) {
889 y[k] = mixFrac * fuelComp[k]/sum_yf + (1.0-mixFrac) * oxComp[k]/sum_yo;
899 const string& element)
const
910 parseCompString(fuelComp.find(
":") != string::npos ? fuelComp : fuelComp+
":1.0"),
911 parseCompString(oxComp.find(
":") != string::npos ? oxComp : oxComp+
":1.0"),
918 vector<double> fuel, ox;
919 if (basis == ThermoBasis::molar) {
924 fuelComp = fuel.data();
928 if (element ==
"Bilger")
934 if (o2_required_fuel < 0.0 || o2_required_ox > 0.0) {
936 "Fuel composition contains too much oxygen or "
937 "oxidizer contains not enough oxygen. "
938 "Fuel and oxidizer composition mixed up?");
941 double denominator = o2_required_fuel - o2_required_ox;
943 if (denominator == 0.0) {
945 "Fuel and oxidizer have the same composition");
948 double Z = (o2_required_mix - o2_required_ox) / denominator;
950 return std::min(std::max(Z, 0.0), 1.0);
953 double sum_yf = std::accumulate(fuelComp, fuelComp+
m_kk, 0.0);
954 double sum_yo = std::accumulate(oxComp, oxComp+
m_kk, 0.0);
956 if (sum_yf == 0.0 || sum_yo == 0.0) {
958 "No fuel and/or oxidizer composition specified");
961 auto elementalFraction = [
this](
size_t m,
const double* y) {
963 for (
size_t k = 0; k !=
m_kk; ++k) {
970 double Z_m_fuel = elementalFraction(m, fuelComp)/sum_yf;
971 double Z_m_ox = elementalFraction(m, oxComp)/sum_yo;
974 if (Z_m_fuel == Z_m_ox) {
976 "Fuel and oxidizer have the same composition for element {}",
979 double Z = (Z_m_mix - Z_m_ox) / (Z_m_fuel - Z_m_ox);
980 return std::min(std::max(Z, 0.0), 1.0);
997 if (inputFile.empty()) {
1001 size_t dot = inputFile.find_last_of(
".");
1004 extension = inputFile.substr(
dot+1);
1006 if (extension ==
"xml" || extension ==
"cti") {
1008 "The CTI and XML formats are no longer supported.");
1012 auto& phase = root[
"phases"].getMapWhere(
"name",
id);
1021 "Missing species thermo data");
1033 "Temperature ({}), pressure ({}) and vapor fraction ({}) "
1034 "are inconsistent, above the critical temperature.",
1040 if (std::abs(Psat / P - 1) < 1e-6) {
1042 }
else if ((Q == 0 && P >= Psat) || (Q == 1 && P <= Psat)) {
1046 "Temperature ({}), pressure ({}) and vapor fraction ({}) "
1047 "are inconsistent.\nPsat at this T: {}\n"
1048 "Consider specifying the state using two fully independent "
1049 "properties (for example, temperature and density)",
1056 if (!spec->thermo) {
1058 "Species {} has no thermo data", spec->name);
1062 spec->thermo->validate(spec->name);
1070 if (!spec->thermo) {
1072 "Species {} has no thermo data", spec->name);
1077 "New species '{}' does not match existing species '{}' at index {}",
1080 spec->thermo->validate(spec->name);
1101 phaseNode[
"name"] =
name();
1104 for (
size_t i = 0; i <
nElements(); i++) {
1112 if (stateVars.count(
"T")) {
1116 if (stateVars.count(
"D")) {
1117 state[
"density"].setQuantity(
density(),
"kg/m^3");
1118 }
else if (stateVars.count(
"P")) {
1119 state[
"P"].setQuantity(
pressure(),
"Pa");
1122 if (stateVars.count(
"Y")) {
1123 map<string, double> Y;
1124 for (
size_t k = 0; k <
m_kk; k++) {
1132 }
else if (stateVars.count(
"X")) {
1133 map<string, double> X;
1134 for (
size_t k = 0; k <
m_kk; k++) {
1144 phaseNode[
"state"] = std::move(state);
1148 phaseNode[
"__type__"] =
"Phase";
1168 double rtol,
int max_steps,
int max_iter,
1169 int estimate_equil,
int log_level)
1171 if (solver ==
"auto" || solver ==
"element_potential") {
1172 vector<double> initial_state;
1174 debuglog(
"Trying ChemEquil solver\n", log_level);
1179 int ret = E.
equilibrate(*
this, XY.c_str(), log_level-1);
1182 "ChemEquil solver failed. Return code: {}", ret);
1184 debuglog(
"ChemEquil solver succeeded\n", log_level);
1186 }
catch (std::exception& err) {
1187 debuglog(
"ChemEquil solver failed.\n", log_level);
1190 if (solver ==
"auto") {
1197 if (solver ==
"auto" || solver ==
"vcs" || solver ==
"gibbs") {
1201 M.
equilibrate(XY, solver, rtol, max_steps, max_iter,
1202 estimate_equil, log_level);
1206 if (solver !=
"auto") {
1208 "Invalid solver specified: '{}'", solver);
1214 for (
size_t m = 0; m <
m_kk; m++) {
1215 for (
size_t k = 0; k <
m_kk; k++) {
1216 dlnActCoeffdlnN[ld * k + m] = 0.0;
1222void ThermoPhase::getdlnActCoeffdlnN_numderiv(
const size_t ld,
1223 double*
const dlnActCoeffdlnN)
1225 double deltaMoles_j = 0.0;
1229 vector<double> ActCoeff_Base(
m_kk);
1231 vector<double> Xmol_Base(
m_kk);
1235 vector<double> ActCoeff(
m_kk);
1236 vector<double> Xmol(
m_kk);
1237 double v_totalMoles = 1.0;
1238 double TMoles_base = v_totalMoles;
1241 for (
size_t j = 0; j <
m_kk; j++) {
1247 double moles_j_base = v_totalMoles * Xmol_Base[j];
1248 deltaMoles_j = 1.0E-7 * moles_j_base + v_totalMoles * 1.0E-13 + 1.0E-150;
1252 v_totalMoles = TMoles_base + deltaMoles_j;
1253 for (
size_t k = 0; k <
m_kk; k++) {
1254 Xmol[k] = Xmol_Base[k] * TMoles_base / v_totalMoles;
1256 Xmol[j] = (moles_j_base + deltaMoles_j) / v_totalMoles;
1264 double*
const lnActCoeffCol = dlnActCoeffdlnN + ld * j;
1265 for (
size_t k = 0; k <
m_kk; k++) {
1266 lnActCoeffCol[k] = (2*moles_j_base + deltaMoles_j) *(ActCoeff[k] - ActCoeff_Base[k]) /
1267 ((ActCoeff[k] + ActCoeff_Base[k]) * deltaMoles_j);
1270 v_totalMoles = TMoles_base;
1279 if (
type() ==
"none") {
1281 "Not implemented for thermo model 'none'");
1284 fmt::memory_buffer b;
1286 int name_width = 18;
1288 string blank_leader = fmt::format(
"{:{}}",
"", name_width);
1290 string string_property = fmt::format(
"{{:>{}}} {{}}\n", name_width);
1292 string one_property = fmt::format(
"{{:>{}}} {{:<.5g}} {{}}\n", name_width);
1294 string two_prop_header =
"{} {:^15} {:^15}\n";
1295 string kg_kmol_header = fmt::format(
1296 two_prop_header, blank_leader,
"1 kg",
"1 kmol"
1298 string Y_X_header = fmt::format(
1299 two_prop_header, blank_leader,
"mass frac. Y",
"mole frac. X"
1301 string two_prop_sep = fmt::format(
1302 "{} {:-^15} {:-^15}\n", blank_leader,
"",
""
1304 string two_property = fmt::format(
1305 "{{:>{}}} {{:15.5g}} {{:15.5g}} {{}}\n", name_width
1308 string three_prop_header = fmt::format(
1309 "{} {:^15} {:^15} {:^15}\n", blank_leader,
"mass frac. Y",
1310 "mole frac. X",
"chem. pot. / RT"
1312 string three_prop_sep = fmt::format(
1313 "{} {:-^15} {:-^15} {:-^15}\n", blank_leader,
"",
"",
""
1315 string three_property = fmt::format(
1316 "{{:>{}}} {{:15.5g}} {{:15.5g}} {{:15.5g}}\n", name_width
1332 fmt_append(b, one_property,
"potential", phi,
"V");
1356 "heat capacity c_v",
"<not implemented>");
1360 vector<double> x(
m_kk);
1361 vector<double> y(
m_kk);
1362 vector<double> mu(
m_kk);
1367 double xMinor = 0.0;
1368 double yMinor = 0.0;
1373 for (
size_t k = 0; k <
m_kk; k++) {
1374 if (abs(x[k]) >= threshold) {
1390 for (
size_t k = 0; k <
m_kk; k++) {
1391 if (abs(x[k]) >= threshold) {
1401 string minor = fmt::format(
"[{:+5d} minor]", nMinor);
1402 fmt_append(b, two_property, minor, yMinor, xMinor,
"");
1405 return to_string(b) + err.
what();
1407 return to_string(b);
Headers for the MultiPhase object that is used to set up multiphase equilibrium problems (see Chemica...
Pure Virtual Base class for individual species reference state thermodynamic managers and text for th...
Declaration for class Cantera::Species.
Headers for the factory class that can create known ThermoPhase objects (see Thermodynamic Properties...
Header file for class ThermoPhase, the base class for phases with thermodynamic properties,...
A map of string keys to values whose type can vary at runtime.
size_t size() const
Returns the number of elements in this map.
bool hasKey(const string &key) const
Returns true if the map contains an item named key.
double convert(const string &key, const string &units) const
Convert the item stored by the given key to the units specified in units.
void setFlowStyle(bool flow=true)
Use "flow" style when outputting this AnyMap to YAML.
void erase(const string &key)
Erase the value held by key.
static AnyMap fromYamlFile(const string &name, const string &parent_name="")
Create an AnyMap from a YAML file.
void update(const AnyMap &other, bool keepExisting=true)
Add items from other to this AnyMap.
static bool addOrderingRules(const string &objectType, const vector< vector< string > > &specs)
Add global rules for setting the order of elements when outputting AnyMap objects to YAML.
string keys_str() const
Return a string listing the keys in this AnyMap, for use in error messages, for example.
Base class for exceptions thrown by Cantera classes.
const char * what() const override
Get a description of the error.
Class ChemEquil implements a chemical equilibrium solver for single-phase solutions.
int equilibrate(ThermoPhase &s, const char *XY, int loglevel=0)
Equilibrate a phase, holding the elemental composition fixed at the initial value found within the Th...
EquilOpt options
Options controlling how the calculation is carried out.
double relTolerance
Relative tolerance.
int maxIterations
Maximum number of iterations.
string canonicalize(const string &name)
Get the canonical name registered for a type.
A class for multiphase mixtures.
void init()
Process phases and build atomic composition array.
void addPhase(ThermoPhase *p, double moles)
Add a phase to the mixture.
A species thermodynamic property manager for a phase.
bool ready(size_t nSpecies)
Check if data for all species (0 through nSpecies-1) has been installed.
virtual void install_STIT(size_t index, shared_ptr< SpeciesThermoInterpType > stit)
Install a new species thermodynamic property parameterization for one species.
virtual void modifySpecies(size_t index, shared_ptr< SpeciesThermoInterpType > spec)
Modify the species thermodynamic property parameterization for a species.
virtual void resetHf298(const size_t k)
Restore the original heat of formation of one or more species.
An error indicating that an unimplemented function has been called.
double massFraction(size_t k) const
Return the mass fraction of a single species.
virtual bool addSpecies(shared_ptr< Species > spec)
Add a Species to this Phase.
virtual void setMoleFractions(const double *const x)
Set the mole fractions to the specified values.
void assertCompressible(const string &setter) const
Ensure that phase is compressible.
void restoreState(const vector< double > &state)
Restore a state saved on a previous call to saveState.
size_t nSpecies() const
Returns the number of species in the phase.
virtual void setMassFractions_NoNorm(const double *const y)
Set the mass fractions to the specified values without normalizing.
virtual map< string, size_t > nativeState() const
Return a map of properties defining the native state of a substance.
size_t m_kk
Number of species in the phase.
virtual void modifySpecies(size_t k, shared_ptr< Species > spec)
Modify the thermodynamic data associated with a species.
size_t nDim() const
Returns the number of spatial dimensions (1, 2, or 3)
void setState_TD(double t, double rho)
Set the internally stored temperature (K) and density (kg/m^3)
double temperature() const
Temperature (K).
virtual void setPressure(double p)
Set the internally stored pressure (Pa) at constant temperature and composition.
double meanMolecularWeight() const
The mean molecular weight. Units: (kg/kmol)
void moleFractionsToMassFractions(const double *X, double *Y) const
Converts a mixture composition from mass fractions to mole fractions.
void saveState(vector< double > &state) const
Save the current internal state of the phase.
size_t elementIndex(const string &name) const
Return the index of element named 'name'.
void setMassFractionsByName(const Composition &yMap)
Set the species mass fractions by name.
string speciesName(size_t k) const
Name of the species with index k.
virtual void setDensity(const double density_)
Set the internally stored density (kg/m^3) of the phase.
vector< double > getCompositionFromMap(const Composition &comp) const
Converts a Composition to a vector with entries for each species Species that are not specified are s...
void getMoleFractions(double *const x) const
Get the species mole fraction vector.
void setMoleFractionsByName(const Composition &xMap)
Set the species mole fractions by name.
const double * massFractions() const
Return a const pointer to the mass fraction array.
const vector< double > & molecularWeights() const
Return a const reference to the internal vector of molecular weights.
double moleFraction(size_t k) const
Return the mole fraction of a single species.
const vector< string > & elementNames() const
Return a read-only reference to the vector of element names.
virtual double density() const
Density (kg/m^3).
double nAtoms(size_t k, size_t m) const
Number of atoms of element m in species k.
virtual void setTemperature(double temp)
Set the internally stored temperature of the phase (K).
size_t nElements() const
Number of elements.
virtual void setMassFractions(const double *const y)
Set the mass fractions to the specified values and normalize them.
const vector< string > & speciesNames() const
Return a const reference to the vector of species names.
double molecularWeight(size_t k) const
Molecular weight of species k.
virtual void invalidateCache()
Invalidate any cached values which are normally updated only when a change in state is detected.
void getMassFractions(double *const y) const
Get the species mass fractions.
virtual double pressure() const
Return the thermodynamic pressure (Pa).
string elementName(size_t m) const
Name of the element with index m.
double charge(size_t k) const
Dimensionless electrical charge of a single molecule of species k The charge is normalized by the the...
string name() const
Return the name of the phase.
static ThermoFactory * factory()
Static function that creates a static instance of the factory.
int m_ssConvention
Contains the standard state convention.
virtual double critTemperature() const
Critical temperature (K).
virtual void setState_HP(double h, double p, double tol=1e-9)
Set the internally stored specific enthalpy (J/kg) and pressure (Pa) of the phase.
double electricPotential() const
Returns the electric potential of this phase (V).
virtual void setState_UV(double u, double v, double tol=1e-9)
Set the specific internal energy (J/kg) and specific volume (m^3/kg).
virtual double cp_mole() const
Molar heat capacity at constant pressure. Units: J/kmol/K.
double equivalenceRatio() const
Compute the equivalence ratio for the current mixture from available oxygen and required oxygen.
virtual void setParameters(const AnyMap &phaseNode, const AnyMap &rootNode=AnyMap())
Set equation of state parameters from an AnyMap phase description.
virtual void getParameters(AnyMap &phaseNode) const
Store the parameters of a ThermoPhase object such that an identical one could be reconstructed using ...
virtual double enthalpy_mole() const
Molar enthalpy. Units: J/kmol.
virtual void setState_TP(double t, double p)
Set the temperature (K) and pressure (Pa)
virtual double standardConcentration(size_t k=0) const
Return the standard concentration for the kth species.
virtual void setState_TV(double t, double v, double tol=1e-9)
Set the temperature (K) and specific volume (m^3/kg).
virtual double logStandardConc(size_t k=0) const
Natural logarithm of the standard concentration of the kth species.
double o2Present(const double *y) const
Helper function for computing the amount of oxygen available in the current mixture.
virtual void setState_PV(double p, double v, double tol=1e-9)
Set the pressure (Pa) and specific volume (m^3/kg).
virtual void setState(const AnyMap &state)
Set the state using an AnyMap containing any combination of properties supported by the thermodynamic...
virtual double minTemp(size_t k=npos) const
Minimum temperature for which the thermodynamic data for the species or phase are valid.
virtual void setState_TPX(double t, double p, const double *x)
Set the temperature (K), pressure (Pa), and mole fractions.
void setState_SPorSV(double s, double p, double tol=1e-9, bool doSV=false)
Carry out work in SP and SV calculations.
double RT() const
Return the Gas Constant multiplied by the current temperature.
virtual double critPressure() const
Critical pressure (Pa).
virtual void setState_TPY(double t, double p, const double *y)
Set the internally stored temperature (K), pressure (Pa), and mass fractions of the phase.
double m_tlast
last value of the temperature processed by reference state
virtual void setState_ST(double s, double t, double tol=1e-9)
Set the specific entropy (J/kg/K) and temperature (K).
void setState_HPorUV(double h, double p, double tol=1e-9, bool doUV=false)
Carry out work in HP and UV calculations.
double gibbs_mass() const
Specific Gibbs function. Units: J/kg.
virtual void getActivityConcentrations(double *c) const
This method returns an array of generalized concentrations.
double stoichAirFuelRatio(const double *fuelComp, const double *oxComp, ThermoBasis basis=ThermoBasis::molar) const
Compute the stoichiometric air to fuel ratio (kg oxidizer / kg fuel) given fuel and oxidizer composit...
string type() const override
String indicating the thermodynamic model implemented.
AnyMap parameters(bool withInput=true) const
Returns the parameters of a ThermoPhase object such that an identical one could be reconstructed usin...
virtual string report(bool show_thermo=true, double threshold=-1e-14) const
returns a summary of the state of the phase as a string
virtual double maxTemp(size_t k=npos) const
Maximum temperature for which the thermodynamic data for the species are valid.
double mixtureFraction(const double *fuelComp, const double *oxComp, ThermoBasis basis=ThermoBasis::molar, const string &element="Bilger") const
Compute the mixture fraction = kg fuel / (kg oxidizer + kg fuel) for the current mixture given fuel a...
double o2Required(const double *y) const
Helper function for computing the amount of oxygen required for complete oxidation.
void getElectrochemPotentials(double *mu) const
Get the species electrochemical potentials.
virtual void getdlnActCoeffdlnN(const size_t ld, double *const dlnActCoeffdlnN)
Get the array of derivatives of the log activity coefficients with respect to the log of the species ...
virtual void getActivityCoefficients(double *ac) const
Get the array of non-dimensional molar-based activity coefficients at the current solution temperatur...
virtual string phaseOfMatter() const
String indicating the mechanical phase of the matter in this Phase.
virtual void setState_Tsat(double t, double x)
Set the state to a saturated system at a particular temperature.
virtual double entropy_mole() const
Molar entropy. Units: J/kmol/K.
double cv_mass() const
Specific heat at constant volume. Units: J/kg/K.
virtual int activityConvention() const
This method returns the convention used in specification of the activities, of which there are curren...
virtual void initThermo()
Initialize the ThermoPhase object after all species have been set up.
double entropy_mass() const
Specific entropy. Units: J/kg/K.
virtual MultiSpeciesThermo & speciesThermo(int k=-1)
Return a changeable reference to the calculation manager for species reference-state thermodynamic pr...
virtual void setState_UP(double u, double p, double tol=1e-9)
Set the specific internal energy (J/kg) and pressure (Pa).
void initThermoFile(const string &inputFile, const string &id)
Initialize a ThermoPhase object using an input file.
virtual void setState_SP(double s, double p, double tol=1e-9)
Set the specific entropy (J/kg/K) and pressure (Pa).
virtual int standardStateConvention() const
This method returns the convention used in specification of the standard state, of which there are cu...
void modifySpecies(size_t k, shared_ptr< Species > spec) override
Modify the thermodynamic data associated with a species.
virtual void setState_SH(double s, double h, double tol=1e-9)
Set the specific entropy (J/kg/K) and the specific enthalpy (J/kg)
void invalidateCache() override
Invalidate any cached values which are normally updated only when a change in state is detected.
virtual void getActivities(double *a) const
Get the array of non-dimensional activities at the current solution temperature, pressure,...
void setMixtureFraction(double mixFrac, const double *fuelComp, const double *oxComp, ThermoBasis basis=ThermoBasis::molar)
Set the mixture composition according to the mixture fraction = kg fuel / (kg oxidizer + kg fuel)
virtual void resetHf298(const size_t k=npos)
Restore the original heat of formation of one or more species.
virtual void getChemPotentials(double *mu) const
Get the species chemical potentials. Units: J/kmol.
double cp_mass() const
Specific heat at constant pressure. Units: J/kg/K.
virtual void setState_TH(double t, double h, double tol=1e-9)
Set the temperature (K) and the specific enthalpy (J/kg)
virtual void getLnActivityCoefficients(double *lnac) const
Get the array of non-dimensional molar-based ln activity coefficients at the current solution tempera...
double intEnergy_mass() const
Specific internal energy. Units: J/kg.
virtual Units standardConcentrationUnits() const
Returns the units of the "standard concentration" for this phase.
virtual double cv_mole() const
Molar heat capacity at constant volume. Units: J/kmol/K.
MultiSpeciesThermo m_spthermo
Pointer to the calculation manager for species reference-state thermodynamic properties.
virtual double satPressure(double t)
Return the saturation pressure given the temperature.
bool addSpecies(shared_ptr< Species > spec) override
Add a Species to this Phase.
AnyMap m_input
Data supplied via setParameters.
virtual double intEnergy_mole() const
Molar internal energy. Units: J/kmol.
virtual void setState_DP(double rho, double p)
Set the density (kg/m**3) and pressure (Pa) at constant composition.
void setEquivalenceRatio(double phi, const double *fuelComp, const double *oxComp, ThermoBasis basis=ThermoBasis::molar)
Set the mixture composition according to the equivalence ratio.
void setState_TPQ(double T, double P, double Q)
Set the temperature, pressure, and vapor fraction (quality).
virtual void setState_VH(double v, double h, double tol=1e-9)
Set the specific volume (m^3/kg) and the specific enthalpy (J/kg)
virtual double gibbs_mole() const
Molar Gibbs function. Units: J/kmol.
virtual void setState_SV(double s, double v, double tol=1e-9)
Set the specific entropy (J/kg/K) and specific volume (m^3/kg).
const AnyMap & input() const
Access input data associated with the phase description.
virtual void setState_Psat(double p, double x)
Set the state to a saturated system at a particular pressure.
void setState_conditional_TP(double t, double p, bool set_p)
Helper function used by setState_HPorUV and setState_SPorSV.
double enthalpy_mass() const
Specific enthalpy. Units: J/kg.
A representation of the units associated with a dimensional quantity.
void fmt_append(fmt::memory_buffer &b, Args... args)
Versions 6.2.0 and 6.2.1 of fmtlib do not include this define before they include windows....
Composition parseCompString(const string &ss, const vector< string > &names)
Parse a composition string into a map consisting of individual key:composition pairs.
void equilibrate(const string &XY, const 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 ThermoPhase object.
void debuglog(const string &msg, int loglevel)
Write a message to the log only if loglevel > 0.
double dot(InputIter x_begin, InputIter x_end, InputIter2 y_begin)
Function that calculates a templated inner product.
T clip(const T &value, const T &lower, const T &upper)
Clip value such that lower <= value <= upper.
const double Faraday
Faraday constant [C/kmol].
const double OneAtm
One atmosphere [Pa].
void setupPhase(ThermoPhase &thermo, const AnyMap &phaseNode, const AnyMap &rootNode)
Initialize a ThermoPhase object.
Namespace for the Cantera kernel.
const size_t npos
index returned by functions to indicate "no position"
const int cAC_CONVENTION_MOLAR
Standard state uses the molar convention.
const double SmallNumber
smallest number to compare to zero.
ThermoBasis
Differentiate between mole fractions and mass fractions for input mixture composition.
map< string, double > Composition
Map from string names to doubles.
Contains declarations for string manipulation functions within Cantera.