12#include <boost/algorithm/string.hpp>
28 double a0,
double a1,
double b)
36 size_t counter = k +
m_kk * k;
37 a_coeff_vec(0, counter) = a0;
38 a_coeff_vec(1, counter) = a1;
41 for (
size_t j = 0; j <
m_kk; j++) {
47 if (isnan(a_coeff_vec(0, j +
m_kk * j))) {
50 }
else if (isnan(a_coeff_vec(0, j +
m_kk * k))) {
53 double a0kj = sqrt(a_coeff_vec(0, j +
m_kk * j) * a0);
54 double a1kj = sqrt(a_coeff_vec(1, j +
m_kk * j) * a1);
55 a_coeff_vec(0, j +
m_kk * k) = a0kj;
56 a_coeff_vec(1, j +
m_kk * k) = a1kj;
57 a_coeff_vec(0, k +
m_kk * j) = a0kj;
58 a_coeff_vec(1, k +
m_kk * j) = a1kj;
61 a_coeff_vec.
getRow(0, a_vec_Curr_);
77 size_t counter1 = ki +
m_kk * kj;
78 size_t counter2 = kj +
m_kk * ki;
79 a_coeff_vec(0, counter1) = a_coeff_vec(0, counter2) = a0;
80 a_coeff_vec(1, counter1) = a_coeff_vec(1, counter2) = a1;
81 a_vec_Curr_[counter1] = a_vec_Curr_[counter2] = a0;
94 double sqt = sqrt(TKelvin);
98 double dadt = da_dt();
99 double fac = TKelvin * dadt - 3.0 *
m_a_current / 2.0;
101 +1.0/(
m_b_current * sqt) * log(vpb/mv) * (-0.5 * dadt));
113 double sqt = sqrt(TKelvin);
116 double dadt = da_dt();
117 double fac = TKelvin * dadt - 3.0 *
m_a_current / 2.0;
118 return (cvref - 1.0/(2.0 *
m_b_current * TKelvin * sqt) * log(vpb/mv)*fac
119 +1.0/(
m_b_current * sqt) * log(vpb/mv)*(-0.5*dadt));
142 for (
size_t k = 0; k <
m_kk; k++) {
153 for (
size_t k = 0; k <
m_kk; k++) {
155 for (
size_t i = 0; i <
m_kk; i++) {
156 size_t counter = k +
m_kk*i;
162 for (
size_t k = 0; k <
m_kk; k++) {
163 ac[k] = (-
RT() * log(pres * mv /
RT())
164 +
RT() * log(mv / vmb)
165 +
RT() * b_vec_Curr_[k] / vmb
171 for (
size_t k = 0; k <
m_kk; k++) {
172 ac[k] = exp(ac[k]/
RT());
181 for (
size_t k = 0; k <
m_kk; k++) {
183 mu[k] +=
RT() * log(xx);
194 for (
size_t k = 0; k <
m_kk; k++) {
196 for (
size_t i = 0; i <
m_kk; i++) {
197 size_t counter = k +
m_kk*i;
203 for (
size_t k = 0; k <
m_kk; k++) {
204 mu[k] += (-
RT() * log(pres * mv /
RT())
205 +
RT() * log(mv / vmb)
206 +
RT() * b_vec_Curr_[k] / vmb
218 scale(hbar.begin(), hbar.end(), hbar.begin(),
RT());
226 double sqt = sqrt(TKelvin);
229 for (
size_t k = 0; k <
m_kk; k++) {
231 for (
size_t i = 0; i <
m_kk; i++) {
232 size_t counter = k +
m_kk*i;
236 for (
size_t k = 0; k <
m_kk; k++) {
237 dpdni_[k] =
RT()/vmb +
RT() * b_vec_Curr_[k] / (vmb * vmb) - 2.0 *
m_pp[k] / (sqt * mv * vpb)
238 +
m_a_current * b_vec_Curr_[k]/(sqt * mv * vpb * vpb);
240 double dadt = da_dt();
241 double fac = TKelvin * dadt - 3.0 *
m_a_current / 2.0;
243 for (
size_t k = 0; k <
m_kk; k++) {
245 for (
size_t i = 0; i <
m_kk; i++) {
246 size_t counter = k +
m_kk*i;
253 for (
size_t k = 0; k <
m_kk; k++) {
256 + b_vec_Curr_[k] / vpb / (
m_b_current * sqt) * fac);
257 hbar[k] = hbar[k] + hE_v;
258 hbar[k] -= fac2 *
dpdni_[k];
267 double sqt = sqrt(TKelvin);
272 for (
size_t k = 0; k <
m_kk; k++) {
279 for (
size_t k = 0; k <
m_kk; k++) {
281 for (
size_t i = 0; i <
m_kk; i++) {
282 size_t counter = k +
m_kk*i;
286 for (
size_t k = 0; k <
m_kk; k++) {
288 for (
size_t i = 0; i <
m_kk; i++) {
289 size_t counter = k +
m_kk*i;
294 double dadt = da_dt();
298 for (
size_t k = 0; k <
m_kk; k++) {
306 + 1.0 / (
m_b_current * sqt) * b_vec_Curr_[k] / vpb * fac);
311 for (
size_t k = 0; k <
m_kk; k++) {
312 sbar[k] += m_partialMolarVolumes[k] *
dpdT_;
322 for (
size_t k = 0; k <
nSpecies(); k++) {
323 ubar[k] -= p * m_partialMolarVolumes[k];
329 checkArraySize(
"RedlichKwongMFTP::getPartialMolarIntEnergies_TV", utilde.size(),
m_kk);
333 double sqt = sqrt(T);
337 double dadt = da_dt();
339 for (
size_t k = 0; k <
m_kk; k++) {
348 for (
size_t k = 0; k <
m_kk; k++) {
351 for (
size_t i = 0; i <
m_kk; i++) {
352 size_t counter = k +
m_kk * i;
359 double logv = log(vpb / mv);
360 double F = T * dadt - 1.5 * a;
361 double C = -logv / b + 1.0 / vpb;
362 double pref = 1.0 / (b * sqt);
364 for (
size_t k = 0; k <
m_kk; k++) {
366 double ures = pref * (logv * Sk + b_vec_Curr_[k] * F * C);
377 double sqt = sqrt(T);
381 double dadt = da_dt();
384 for (
size_t k = 0; k <
m_kk; k++) {
396 double dvdT_P = -pT / pV;
403 double pTT = (-d2adt2 + dadt / T - 3.0 * a / (4.0 * T * T)) / (sqt * v * vpb);
405 + (dadt - a / (2.0 * T)) * (2.0 * v + b) / (sqt * v * v * vpb * vpb);
406 double pVV = 2.0 *
GasConstant * T / (vmb * vmb * vmb)
407 - 2.0 * a * (3.0 * v * v + 3.0 * b * v + b * b)
408 / (sqt * v * v * v * vpb * vpb * vpb);
409 double d2vdT2_P = -(pTT + 2.0 * pTV * dvdT_P + pVV * dvdT_P * dvdT_P) / pV;
412 for (
size_t k = 0; k <
m_kk; k++) {
415 for (
size_t i = 0; i <
m_kk; i++) {
416 size_t counter = k +
m_kk * i;
422 double F = T * dadt - 1.5 * a;
423 double dF_dT = -0.5 * dadt + T * d2adt2;
424 double logvpv = log(vpb / v);
425 double termLPrime = -b * dvdT_P / (v * vpb);
427 for (
size_t k = 0; k <
m_kk; k++) {
428 double bk = b_vec_Curr_[k];
431 double Sk = 2.0 * T * dAk - 3.0 * Ak;
432 double dSk_dT = -dAk;
435 double dpdni =
RT() / vmb +
RT() * bk / (vmb * vmb)
436 - 2.0 * Ak / (sqt * v * vpb)
437 + a * bk / (sqt * v * vpb * vpb);
443 - 2.0 *
GasConstant * T * bk * dvdT_P / (vmb * vmb * vmb)
444 - 2.0 * dAk / (sqt * v * vpb)
445 + Ak / (T * sqt * v * vpb)
446 + 2.0 * Ak * (2.0 * v + b) * dvdT_P / (sqt * v * v * vpb * vpb)
447 + bk * dadt / (sqt * v * vpb * vpb)
448 - bk * a / (2.0 * T * sqt * v * vpb * vpb)
449 - bk * a * (3.0 * v + b) * dvdT_P / (sqt * v * v * vpb * vpb * vpb);
452 double dterm1 = -bk / (b * b * sqt)
453 * (termLPrime * F + logvpv * dF_dT - logvpv * F / (2.0 * T));
454 double dterm2 = 1.0 / (b * sqt)
455 * (termLPrime * Sk + logvpv * dSk_dT - logvpv * Sk / (2.0 * T));
456 double dterm3 = bk / (b * sqt)
457 * (dF_dT / vpb - F * dvdT_P / (vpb * vpb) - F / (2.0 * T * vpb));
461 + (dvdT_P + T * d2vdT2_P) * dpdni
462 + T * dvdT_P * d_dpdni_dT
463 + dterm1 + dterm2 + dterm3;
473 double sqt = sqrt(T);
477 double dadt = da_dt();
479 for (
size_t k = 0; k <
m_kk; k++) {
488 for (
size_t k = 0; k <
m_kk; k++) {
491 for (
size_t i = 0; i <
m_kk; i++) {
492 size_t counter = k +
m_kk * i;
502 double logv = log(vpb / mv);
503 double F = T * dadt - 1.5 * a;
504 double C = -logv / b + 1.0 / vpb;
505 double pref = 1.0 / (b * sqt);
507 for (
size_t k = 0; k <
m_kk; k++) {
509 double ures = pref * (logv * Sk + b_vec_Curr_[k] * F * C);
510 double dresdT = pref * (-logv *
m_dAkdT[k] - 0.5 * b_vec_Curr_[k] * dadt * C)
512 cvtilde[k] += dresdT;
519 for (
size_t k = 0; k <
m_kk; k++) {
521 for (
size_t i = 0; i <
m_kk; i++) {
522 size_t counter = k +
m_kk*i;
526 for (
size_t k = 0; k <
m_kk; k++) {
528 for (
size_t i = 0; i <
m_kk; i++) {
529 size_t counter = k +
m_kk*i;
538 for (
size_t k = 0; k <
m_kk; k++) {
541 - 2.0 *
m_pp[k] / (sqt * vpb)
546 vbar[k] = num / denom;
554 a_vec_Curr_.resize(
m_kk *
m_kk, 0.0);
558 b_vec_Curr_.push_back(NAN);
564 m_partialMolarVolumes.push_back(0.0);
574 std::unordered_map<string, AnyMap*> dbSpecies;
579 if (!isnan(a_coeff_vec(0, k +
m_kk * k))) {
582 bool foundCoeffs =
false;
584 if (data.hasKey(
"equation-of-state") &&
585 data[
"equation-of-state"].hasMapWhere(
"model",
"Redlich-Kwong"))
589 auto eos = data[
"equation-of-state"].getMapWhere(
590 "model",
"Redlich-Kwong");
592 if (eos.hasKey(
"a") && eos.hasKey(
"b")) {
593 double a0 = 0, a1 = 0;
594 if (eos[
"a"].isScalar()) {
595 a0 = eos.convert(
"a",
"Pa*m^6/kmol^2*K^0.5");
597 auto avec = eos[
"a"].asVector<
AnyValue>(2);
598 a0 = eos.units().convert(avec[0],
"Pa*m^6/kmol^2*K^0.5");
599 a1 = eos.units().convert(avec[1],
"Pa*m^6/kmol^2/K^0.5");
601 double b = eos.convert(
"b",
"m^3/kmol");
607 if (eos.hasKey(
"binary-a")) {
610 for (
auto& [name2, coeff] : binary_a) {
611 double a0 = 0, a1 = 0;
612 if (coeff.isScalar()) {
613 a0 = units.
convert(coeff,
"Pa*m^6/kmol^2*K^0.5");
615 auto avec = coeff.asVector<
AnyValue>(2);
616 a0 = units.convert(avec[0],
"Pa*m^6/kmol^2*K^0.5");
617 a1 = units.convert(avec[1],
"Pa*m^6/kmol^2/K^0.5");
629 double Tc = NAN, Pc = NAN;
630 if (data.hasKey(
"critical-parameters")) {
633 auto& critProps = data[
"critical-parameters"].as<
AnyMap>();
634 Tc = critProps.
convert(
"critical-temperature",
"K");
635 Pc = critProps.convert(
"critical-pressure",
"Pa");
639 if (critPropsDb.
empty()) {
641 dbSpecies = critPropsDb[
"species"].asMap(
"name");
645 auto ucName = boost::algorithm::to_upper_copy(
name);
646 if (dbSpecies.count(ucName)) {
647 auto& spec = *dbSpecies.at(ucName);
648 auto& critProps = spec[
"critical-parameters"].as<
AnyMap>();
649 Tc = critProps.
convert(
"critical-temperature",
"K");
650 Pc = critProps.convert(
"critical-pressure",
"Pa");
663 "No critical property or Redlich-Kwong parameters found "
664 "for species {}.",
name);
670 AnyMap& speciesNode)
const
675 auto& eosNode = speciesNode[
"equation-of-state"].getMapWhere(
676 "model",
"Redlich-Kwong",
true);
678 size_t counter = k +
m_kk * k;
679 if (a_coeff_vec(1, counter) != 0.0) {
680 vector<AnyValue> coeffs(2);
681 coeffs[0].setQuantity(a_coeff_vec(0, counter),
"Pa*m^6/kmol^2*K^0.5");
682 coeffs[1].setQuantity(a_coeff_vec(1, counter),
"Pa*m^6/kmol^2/K^0.5");
683 eosNode[
"a"] = std::move(coeffs);
685 eosNode[
"a"].setQuantity(a_coeff_vec(0, counter),
686 "Pa*m^6/kmol^2*K^0.5");
688 eosNode[
"b"].setQuantity(b_vec_Curr_[k],
"m^3/kmol");
690 auto& critProps = speciesNode[
"critical-parameters"];
691 double a = a_coeff_vec(0, k +
m_kk * k);
692 double b = b_vec_Curr_[k];
695 critProps[
"critical-temperature"].setQuantity(Tc,
"K");
696 critProps[
"critical-pressure"].setQuantity(Pc,
"Pa");
700 auto& eosNode = speciesNode[
"equation-of-state"].getMapWhere(
701 "model",
"Redlich-Kwong",
true);
704 if (coeffs.second == 0) {
705 bin_a[name2].setQuantity(coeffs.first,
"Pa*m^6/kmol^2*K^0.5");
707 vector<AnyValue> C(2);
708 C[0].setQuantity(coeffs.first,
"Pa*m^6/kmol^2*K^0.5");
709 C[1].setQuantity(coeffs.second,
"Pa*m^6/kmol^2/K^0.5");
710 bin_a[name2] = std::move(C);
713 eosNode[
"binary-a"] = std::move(bin_a);
725 double molarV = mmw / rho;
728 double dadt = da_dt();
730 double sqT = sqrt(T);
744 double molarV = mmw / rho;
747 double dadt = da_dt();
749 double sqT = sqrt(T);
760 double pres = std::max(
psatEst(TKelvin), presGuess);
762 bool foundLiq =
false;
764 while (m < 100 && !foundLiq) {
765 int nsol =
solveCubic(TKelvin, pres, atmp, btmp, Vroot);
766 if (nsol == 1 || nsol == 2) {
793 if (rhoguess == -1.0) {
794 if (phaseRequested >= FLUID_LIQUID_0) {
796 rhoguess = mmw / lqvol;
804 double volguess = mmw / rhoguess;
808 double vmin = std::max(0.0,
m_b_current * (1.0 + 1e-12));
809 vector<double> physicalRoots;
810 for (
double root : Vroot_) {
811 if (std::isfinite(
root) &&
root > vmin) {
812 physicalRoots.push_back(
root);
815 std::sort(physicalRoots.begin(), physicalRoots.end());
816 if (physicalRoots.empty()) {
818 }
else if (physicalRoots.size() == 1) {
822 if ((phaseRequested == FLUID_GAS && NSolns_ < 0)
823 || (phaseRequested >= FLUID_LIQUID_0 && NSolns_ > 0))
828 return mmw / physicalRoots[0];
829 }
else if (physicalRoots.size() == 2) {
830 Vroot_[0] = physicalRoots[0];
831 Vroot_[1] = 0.5 * (physicalRoots[0] + physicalRoots[1]);
832 Vroot_[2] = physicalRoots[1];
834 Vroot_[0] = physicalRoots[0];
835 Vroot_[1] = physicalRoots[1];
836 Vroot_[2] = physicalRoots[2];
839 double molarVolLast = Vroot_[0];
841 if (phaseRequested >= FLUID_LIQUID_0) {
842 molarVolLast = Vroot_[0];
843 }
else if (phaseRequested == FLUID_GAS || phaseRequested == FLUID_SUPERCRIT) {
844 molarVolLast = Vroot_[2];
846 if (volguess > Vroot_[1]) {
847 molarVolLast = Vroot_[2];
849 molarVolLast = Vroot_[0];
852 }
else if (NSolns_ == 1) {
853 if (phaseRequested == FLUID_GAS || phaseRequested == FLUID_SUPERCRIT || phaseRequested == FLUID_UNDEFINED) {
854 molarVolLast = Vroot_[0];
858 }
else if (NSolns_ == -1) {
859 if (phaseRequested >= FLUID_LIQUID_0 || phaseRequested == FLUID_UNDEFINED || phaseRequested == FLUID_SUPERCRIT) {
860 molarVolLast = Vroot_[0];
861 }
else if (TKelvin > tcrit) {
862 molarVolLast = Vroot_[0];
867 molarVolLast = Vroot_[0];
870 return mmw / molarVolLast;
875 double sqt = sqrt(TKelvin);
881 double dpdv = (-
GasConstant * TKelvin / (vmb * vmb)
917 double sqt = sqrt(TKelvin);
920 double dadt = da_dt();
929 for (
size_t i = 0; i <
m_kk; i++) {
930 for (
size_t j = 0; j <
m_kk; j++) {
931 size_t counter = i *
m_kk + j;
932 a_vec_Curr_[counter] = a_coeff_vec(0,counter) + a_coeff_vec(1,counter) * temp;
939 for (
size_t i = 0; i <
m_kk; i++) {
941 for (
size_t j = 0; j <
m_kk; j++) {
947 fmt::memory_buffer b;
948 for (
size_t k = 0; k <
m_kk; k++) {
949 if (isnan(b_vec_Curr_[k])) {
957 throw CanteraError(
"RedlichKwongMFTP::updateMixingExpressions",
958 "Missing Redlich-Kwong coefficients for species: {}", to_string(b));
967 for (
size_t i = 0; i <
m_kk; i++) {
969 for (
size_t j = 0; j <
m_kk; j++) {
970 size_t counter = i *
m_kk + j;
971 double a_vec_Curr = a_coeff_vec(0,counter) + a_coeff_vec(1,counter) * temp;
976 for (
size_t i = 0; i <
m_kk; i++) {
978 for (
size_t j = 0; j <
m_kk; j++) {
979 size_t counter = i *
m_kk + j;
980 double a_vec_Curr = a_coeff_vec(0,counter);
987double RedlichKwongMFTP::da_dt()
const
991 for (
size_t i = 0; i <
m_kk; i++) {
992 for (
size_t j = 0; j <
m_kk; j++) {
993 size_t counter = i *
m_kk + j;
1001void RedlichKwongMFTP::calcCriticalConditions(
double& pc,
double& tc,
double& vc)
const
1005 for (
size_t i = 0; i <
m_kk; i++) {
1006 for (
size_t j = 0; j <
m_kk; j++) {
1007 size_t counter = i +
m_kk * j;
1031 double sqrttc, f, dfdt, deltatc;
1037 for (
int j = 0; j < 10; j++) {
1041 deltatc = - f / dfdt;
1044 if (deltatc > 0.1) {
1045 throw CanteraError(
"RedlichKwongMFTP::calcCriticalConditions",
1055 span<double> Vroot)
const
1061 double sqt = sqrt(T);
1062 double cn = - (
GasConstant * T * b / pres - a/(pres * sqt) + b * b);
1063 double dn = - (a * b / (pres * sqt));
1067 double tc = pow(tmp, pp);
1071 int nSolnValues =
MixtureFugacityTP::solveCubic(T, pres, a, b, a, Vroot, an, bn, cn, dn, tc, vc);
Declaration for class Cantera::Species.
Headers for the factory class that can create known ThermoPhase objects (see Thermodynamic Properties...
A map of string keys to values whose type can vary at runtime.
const UnitSystem & units() const
Return the default units that should be used to convert stored values.
bool empty() const
Return boolean indicating whether AnyMap is empty.
double convert(const string &key, const string &units) const
Convert the item stored by the given key to the units specified in units.
static AnyMap fromYamlFile(const string &name, const string &parent_name="")
Create an AnyMap from a YAML file.
A wrapper for a variable whose type is determined at runtime.
void getRow(size_t n, span< double > rw) const
Get the nth row and return it in a vector.
virtual void resize(size_t n, size_t m, double v=0.0)
Resize the array, and fill the new entries with 'v'.
Base class for exceptions thrown by Cantera classes.
void getEntropy_R_ref(span< double > er) const override
Returns the vector of nondimensional entropies of the reference state at the current temperature of t...
double critPressure() const override
Critical pressure (Pa).
vector< double > m_h0_RT
Temporary storage for dimensionless reference state enthalpies.
void getStandardChemPotentials(span< double > mu) const override
Get the array of chemical potentials at unit activity.
double critTemperature() const override
Critical temperature (K).
virtual void _updateReferenceStateThermo() const
Updates the reference state thermodynamic functions at the current T of the solution.
void getEnthalpy_RT_ref(span< double > hrt) const override
Returns the vector of nondimensional enthalpies of the reference state at the current temperature of ...
vector< double > moleFractions_
Storage for the current values of the mole fractions of the species.
int solveCubic(double T, double pres, double a, double b, double aAlpha, span< double > Vroot, double an, double bn, double cn, double dn, double tc, double vc) const
Solve the cubic equation of state.
void setTemperature(const double temp) override
Set the temperature of the phase.
virtual double psatEst(double TKelvin) const
Estimate for the saturation pressure.
void getStandardVolumes(span< double > vol) const override
Get the molar volumes of each species in their standard states at the current T and P of the solution...
vector< double > m_cp0_R
Temporary storage for dimensionless reference state heat capacities.
bool addSpecies(shared_ptr< Species > spec) override
Add a Species to this Phase.
double z() const
Calculate the value of z.
vector< double > m_workS
Vector of size m_kk, used as a temporary holding area.
size_t nSpecies() const
Returns the number of species in the phase.
size_t m_kk
Number of species in the phase.
size_t speciesIndex(const string &name, bool raise=true) const
Returns the index of a species named 'name' within the Phase object.
double temperature() const
Temperature (K).
double meanMolecularWeight() const
The mean molecular weight. Units: (kg/kmol)
string speciesName(size_t k) const
Name of the species with index k.
map< string, shared_ptr< Species > > m_species
Map of Species objects.
double mean_X(span< const double > Q) const
Evaluate the mole-fraction-weighted mean of an array Q.
double moleFraction(size_t k) const
Return the mole fraction of a single species.
virtual double density() const
Density (kg/m^3).
shared_ptr< Species > species(const string &name) const
Return the Species object for the named species.
virtual double molarVolume() const
Molar volume (m^3/kmol).
string name() const
Return the name of the phase.
double dpdT_
The derivative of the pressure wrt the temperature.
void calculateAB(double temp, double &aCalc, double &bCalc) const
Calculate the a and the b parameters given the temperature.
double thermalExpansionCoeff() const override
Return the volumetric thermal expansion coefficient. Units: 1/K.
void setBinaryCoeffs(const string &species_i, const string &species_j, double a0, double a1)
Set values for the interaction parameter between two species.
vector< double > m_dAkdT
Temporary storage for dA_k/dT; length m_kk.
double sresid() const override
Calculate the deviation terms for the total entropy of the mixture from the ideal gas mixture.
void getPartialMolarEnthalpies(span< double > hbar) const override
Return partial molar enthalpies (J/kmol).
double soundSpeed() const override
Return the speed of sound. Units: m/s.
double pressure() const override
Return the thermodynamic pressure (Pa).
int m_formTempParam
Form of the temperature parameterization.
void getSpeciesParameters(const string &name, AnyMap &speciesNode) const override
Get phase-specific parameters of a Species object such that an identical one could be reconstructed a...
RedlichKwongMFTP(const string &infile="", const string &id="")
Construct a RedlichKwongMFTP object from an input file.
void getPartialMolarCp(span< double > cpbar) const override
Get the partial molar heat capacities at constant pressure [J/kmol/K].
static const double omega_b
Omega constant for b.
vector< double > m_pp
Temporary storage - length = m_kk.
void getActivityCoefficients(span< double > ac) const override
Get the array of non-dimensional activity coefficients at the current solution temperature,...
void initThermo() override
Initialize the ThermoPhase object after all species have been set up.
void getPartialMolarIntEnergies_TV(span< double > utilde) const override
Get the species molar internal energies associated with the derivatives of total internal energy at c...
double internalPressure() const override
Return the internal pressure [Pa].
double cv_mole() const override
Molar heat capacity at constant volume and composition [J/kmol/K].
void pressureDerivatives() const
Calculate dpdV and dpdT at the current conditions.
double isothermalCompressibility() const override
Returns the isothermal compressibility. Units: 1/Pa.
double hresid() const override
Calculate the deviation terms for the total enthalpy of the mixture from the ideal gas mixture.
static const double omega_a
Omega constant for a -> value of a in terms of critical properties.
double liquidVolEst(double TKelvin, double &pres) const override
Estimate for the molar volume of the liquid.
void getPartialMolarCv_TV(span< double > cvtilde) const override
Get the species molar heat capacities associated with the constant volume partial molar internal ener...
double dpdVCalc(double TKelvin, double molarVol, double &presCalc) const override
Calculate the pressure and the pressure derivative given the temperature and the molar volume.
static const double omega_vc
Omega constant for the critical molar volume.
vector< CoeffSource > m_coeffSource
For each species, specifies the source of the a and b coefficients.
void updateMixingExpressions() override
Update the a and b parameters.
double cp_mole() const override
Molar heat capacity at constant pressure and composition [J/kmol/K].
void getPartialMolarVolumes(span< double > vbar) const override
Return an array of partial molar volumes for the species in the mixture.
double m_a_current
Value of a in the equation of state.
void getPartialMolarEntropies(span< double > sbar) const override
Returns an array of partial molar entropies of the species in the solution.
double standardConcentration(size_t k=0) const override
Returns the standard concentration , which is used to normalize the generalized concentration.
double dpdV_
The derivative of the pressure wrt the volume.
bool addSpecies(shared_ptr< Species > spec) override
Add a Species to this Phase.
double m_b_current
Value of b in the equation of state.
void getChemPotentials(span< double > mu) const override
Get the species chemical potentials. Units: J/kmol.
int solveCubic(double T, double pres, double a, double b, span< double > Vroot) const
Prepare variables and call the function to solve the cubic equation of state.
double densityCalc(double T, double pressure, int phase, double rhoguess) override
Calculates the density given the temperature and the pressure and a guess at the density.
vector< double > dpdni_
Vector of derivatives of pressure wrt mole number.
map< string, map< string, pair< double, double > > > m_binaryParameters
Explicitly-specified binary interaction parameters.
void setSpeciesCoeffs(const string &species, double a0, double a1, double b)
Set the pure fluid interaction parameters for a species.
void getPartialMolarIntEnergies(span< double > ubar) const override
Return an array of partial molar internal energies for the species in the mixture.
double RT() const
Return the Gas Constant multiplied by the current temperature.
void initThermoFile(const string &inputFile, const string &id)
Initialize a ThermoPhase object using an input file.
shared_ptr< Solution > root() const
Get the Solution object containing this ThermoPhase object and linked Kinetics and Transport objects.
virtual void getSpeciesParameters(const string &name, AnyMap &speciesNode) const
Get phase-specific parameters of a Species object such that an identical one could be reconstructed a...
virtual double refPressure() const
Returns the reference pressure in Pa.
double convert(double value, const string &src, const string &dest) const
Convert value from the units of src to the units of dest.
void scale(InputIter begin, InputIter end, OutputIter out, S scale_factor)
Multiply elements of an array by a scale factor.
const double GasConstant
Universal Gas Constant [J/kmol/K].
Namespace for the Cantera kernel.
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.
Various templated functions that carry out common vector and polynomial operations (see Templated Arr...