12#include <boost/algorithm/string.hpp>
27 double a0,
double a1,
double b)
35 size_t counter = k +
m_kk * k;
36 a_coeff_vec(0, counter) = a0;
37 a_coeff_vec(1, counter) = a1;
40 for (
size_t j = 0; j <
m_kk; j++) {
46 if (isnan(a_coeff_vec(0, j +
m_kk * j))) {
49 }
else if (isnan(a_coeff_vec(0, j +
m_kk * k))) {
52 double a0kj = sqrt(a_coeff_vec(0, j +
m_kk * j) * a0);
53 double a1kj = sqrt(a_coeff_vec(1, j +
m_kk * j) * a1);
54 a_coeff_vec(0, j +
m_kk * k) = a0kj;
55 a_coeff_vec(1, j +
m_kk * k) = a1kj;
56 a_coeff_vec(0, k +
m_kk * j) = a0kj;
57 a_coeff_vec(1, k +
m_kk * j) = a1kj;
60 a_coeff_vec.
getRow(0, a_vec_Curr_);
76 size_t counter1 = ki +
m_kk * kj;
77 size_t counter2 = kj +
m_kk * ki;
78 a_coeff_vec(0, counter1) = a_coeff_vec(0, counter2) = a0;
79 a_coeff_vec(1, counter1) = a_coeff_vec(1, counter2) = a1;
80 a_vec_Curr_[counter1] = a_vec_Curr_[counter2] = a0;
89 double sqt = sqrt(TKelvin);
94 double dadt = da_dt();
95 double fac = TKelvin * dadt - 3.0 *
m_a_current / 2.0;
97 +1.0/(
m_b_current * sqt) * log(vpb/mv) * (-0.5 * dadt));
105 double sqt = sqrt(TKelvin);
109 double dadt = da_dt();
110 double fac = TKelvin * dadt - 3.0 *
m_a_current / 2.0;
111 return (cvref - 1.0/(2.0 *
m_b_current * TKelvin * sqt) * log(vpb/mv)*fac
112 +1.0/(
m_b_current * sqt) * log(vpb/mv)*(-0.5*dadt));
140 for (
size_t k = 0; k <
m_kk; k++) {
142 for (
size_t i = 0; i <
m_kk; i++) {
143 size_t counter = k +
m_kk*i;
149 for (
size_t k = 0; k <
m_kk; k++) {
150 ac[k] = (-
RT() * log(pres * mv /
RT())
151 +
RT() * log(mv / vmb)
152 +
RT() * b_vec_Curr_[k] / vmb
158 for (
size_t k = 0; k <
m_kk; k++) {
159 ac[k] = exp(ac[k]/
RT());
168 for (
size_t k = 0; k <
m_kk; k++) {
170 mu[k] +=
RT()*(log(xx));
178 for (
size_t k = 0; k <
m_kk; k++) {
180 for (
size_t i = 0; i <
m_kk; i++) {
181 size_t counter = k +
m_kk*i;
188 for (
size_t k = 0; k <
m_kk; k++) {
189 mu[k] += (
RT() * log(pres/refP) -
RT() * log(pres * mv /
RT())
190 +
RT() * log(mv / vmb)
191 +
RT() * b_vec_Curr_[k] / vmb
203 scale(hbar.begin(), hbar.end(), hbar.begin(),
RT());
208 double sqt = sqrt(TKelvin);
211 for (
size_t k = 0; k <
m_kk; k++) {
213 for (
size_t i = 0; i <
m_kk; i++) {
214 size_t counter = k +
m_kk*i;
218 for (
size_t k = 0; k <
m_kk; k++) {
219 dpdni_[k] =
RT()/vmb +
RT() * b_vec_Curr_[k] / (vmb * vmb) - 2.0 *
m_pp[k] / (sqt * mv * vpb)
220 +
m_a_current * b_vec_Curr_[k]/(sqt * mv * vpb * vpb);
222 double dadt = da_dt();
223 double fac = TKelvin * dadt - 3.0 *
m_a_current / 2.0;
225 for (
size_t k = 0; k <
m_kk; k++) {
227 for (
size_t i = 0; i <
m_kk; i++) {
228 size_t counter = k +
m_kk*i;
235 for (
size_t k = 0; k <
m_kk; k++) {
238 + b_vec_Curr_[k] / vpb / (
m_b_current * sqt) * fac);
239 hbar[k] = hbar[k] + hE_v;
240 hbar[k] -= fac2 *
dpdni_[k];
249 double sqt = sqrt(TKelvin);
253 for (
size_t k = 0; k <
m_kk; k++) {
257 for (
size_t k = 0; k <
m_kk; k++) {
259 for (
size_t i = 0; i <
m_kk; i++) {
260 size_t counter = k +
m_kk*i;
264 for (
size_t k = 0; k <
m_kk; k++) {
266 for (
size_t i = 0; i <
m_kk; i++) {
267 size_t counter = k +
m_kk*i;
272 double dadt = da_dt();
276 for (
size_t k = 0; k <
m_kk; k++) {
284 - 1.0 / (
m_b_current * sqt) * b_vec_Curr_[k] / vpb * fac
290 for (
size_t k = 0; k <
m_kk; k++) {
291 sbar[k] -= -m_partialMolarVolumes[k] *
dpdT_;
301 for (
size_t k = 0; k <
nSpecies(); k++) {
302 ubar[k] -= p * m_partialMolarVolumes[k];
308 checkArraySize(
"RedlichKwongMFTP::getPartialMolarIntEnergies_TV", utilde.size(),
m_kk);
312 double sqt = sqrt(T);
316 double dadt = da_dt();
318 for (
size_t k = 0; k <
m_kk; k++) {
327 for (
size_t k = 0; k <
m_kk; k++) {
330 for (
size_t i = 0; i <
m_kk; i++) {
331 size_t counter = k +
m_kk * i;
338 double logv = log(vpb / mv);
339 double F = T * dadt - 1.5 * a;
340 double C = -logv / b + 1.0 / vpb;
341 double pref = 1.0 / (b * sqt);
343 for (
size_t k = 0; k <
m_kk; k++) {
345 double ures = pref * (logv * Sk + b_vec_Curr_[k] * F * C);
356 double sqt = sqrt(T);
360 double dadt = da_dt();
363 for (
size_t k = 0; k <
m_kk; k++) {
375 double dvdT_P = -pT / pV;
382 double pTT = (-d2adt2 + dadt / T - 3.0 * a / (4.0 * T * T)) / (sqt * v * vpb);
384 + (dadt - a / (2.0 * T)) * (2.0 * v + b) / (sqt * v * v * vpb * vpb);
385 double pVV = 2.0 *
GasConstant * T / (vmb * vmb * vmb)
386 - 2.0 * a * (3.0 * v * v + 3.0 * b * v + b * b)
387 / (sqt * v * v * v * vpb * vpb * vpb);
388 double d2vdT2_P = -(pTT + 2.0 * pTV * dvdT_P + pVV * dvdT_P * dvdT_P) / pV;
391 for (
size_t k = 0; k <
m_kk; k++) {
394 for (
size_t i = 0; i <
m_kk; i++) {
395 size_t counter = k +
m_kk * i;
401 double F = T * dadt - 1.5 * a;
402 double dF_dT = -0.5 * dadt + T * d2adt2;
403 double logvpv = log(vpb / v);
404 double termLPrime = -b * dvdT_P / (v * vpb);
406 for (
size_t k = 0; k <
m_kk; k++) {
407 double bk = b_vec_Curr_[k];
410 double Sk = 2.0 * T * dAk - 3.0 * Ak;
411 double dSk_dT = -dAk;
414 double dpdni =
RT() / vmb +
RT() * bk / (vmb * vmb)
415 - 2.0 * Ak / (sqt * v * vpb)
416 + a * bk / (sqt * v * vpb * vpb);
422 - 2.0 *
GasConstant * T * bk * dvdT_P / (vmb * vmb * vmb)
423 - 2.0 * dAk / (sqt * v * vpb)
424 + Ak / (T * sqt * v * vpb)
425 + 2.0 * Ak * (2.0 * v + b) * dvdT_P / (sqt * v * v * vpb * vpb)
426 + bk * dadt / (sqt * v * vpb * vpb)
427 - bk * a / (2.0 * T * sqt * v * vpb * vpb)
428 - bk * a * (3.0 * v + b) * dvdT_P / (sqt * v * v * vpb * vpb * vpb);
431 double dterm1 = -bk / (b * b * sqt)
432 * (termLPrime * F + logvpv * dF_dT - logvpv * F / (2.0 * T));
433 double dterm2 = 1.0 / (b * sqt)
434 * (termLPrime * Sk + logvpv * dSk_dT - logvpv * Sk / (2.0 * T));
435 double dterm3 = bk / (b * sqt)
436 * (dF_dT / vpb - F * dvdT_P / (vpb * vpb) - F / (2.0 * T * vpb));
440 + (dvdT_P + T * d2vdT2_P) * dpdni
441 + T * dvdT_P * d_dpdni_dT
442 + dterm1 + dterm2 + dterm3;
452 double sqt = sqrt(T);
456 double dadt = da_dt();
458 for (
size_t k = 0; k <
m_kk; k++) {
467 for (
size_t k = 0; k <
m_kk; k++) {
470 for (
size_t i = 0; i <
m_kk; i++) {
471 size_t counter = k +
m_kk * i;
481 double logv = log(vpb / mv);
482 double F = T * dadt - 1.5 * a;
483 double C = -logv / b + 1.0 / vpb;
484 double pref = 1.0 / (b * sqt);
486 for (
size_t k = 0; k <
m_kk; k++) {
488 double ures = pref * (logv * Sk + b_vec_Curr_[k] * F * C);
489 double dresdT = pref * (-logv *
m_dAkdT[k] - 0.5 * b_vec_Curr_[k] * dadt * C)
491 cvtilde[k] += dresdT;
498 for (
size_t k = 0; k <
m_kk; k++) {
500 for (
size_t i = 0; i <
m_kk; i++) {
501 size_t counter = k +
m_kk*i;
505 for (
size_t k = 0; k <
m_kk; k++) {
507 for (
size_t i = 0; i <
m_kk; i++) {
508 size_t counter = k +
m_kk*i;
517 for (
size_t k = 0; k <
m_kk; k++) {
520 - 2.0 *
m_pp[k] / (sqt * vpb)
525 vbar[k] = num / denom;
533 a_vec_Curr_.resize(
m_kk *
m_kk, 0.0);
537 b_vec_Curr_.push_back(NAN);
543 m_partialMolarVolumes.push_back(0.0);
553 std::unordered_map<string, AnyMap*> dbSpecies;
558 if (!isnan(a_coeff_vec(0, k +
m_kk * k))) {
561 bool foundCoeffs =
false;
563 if (data.hasKey(
"equation-of-state") &&
564 data[
"equation-of-state"].hasMapWhere(
"model",
"Redlich-Kwong"))
568 auto eos = data[
"equation-of-state"].getMapWhere(
569 "model",
"Redlich-Kwong");
571 if (eos.hasKey(
"a") && eos.hasKey(
"b")) {
572 double a0 = 0, a1 = 0;
573 if (eos[
"a"].isScalar()) {
574 a0 = eos.convert(
"a",
"Pa*m^6/kmol^2*K^0.5");
576 auto avec = eos[
"a"].asVector<
AnyValue>(2);
577 a0 = eos.units().convert(avec[0],
"Pa*m^6/kmol^2*K^0.5");
578 a1 = eos.units().convert(avec[1],
"Pa*m^6/kmol^2/K^0.5");
580 double b = eos.convert(
"b",
"m^3/kmol");
586 if (eos.hasKey(
"binary-a")) {
589 for (
auto& [name2, coeff] : binary_a) {
590 double a0 = 0, a1 = 0;
591 if (coeff.isScalar()) {
592 a0 = units.
convert(coeff,
"Pa*m^6/kmol^2*K^0.5");
594 auto avec = coeff.asVector<
AnyValue>(2);
595 a0 = units.convert(avec[0],
"Pa*m^6/kmol^2*K^0.5");
596 a1 = units.convert(avec[1],
"Pa*m^6/kmol^2/K^0.5");
608 double Tc = NAN, Pc = NAN;
609 if (data.hasKey(
"critical-parameters")) {
612 auto& critProps = data[
"critical-parameters"].as<
AnyMap>();
613 Tc = critProps.
convert(
"critical-temperature",
"K");
614 Pc = critProps.convert(
"critical-pressure",
"Pa");
618 if (critPropsDb.
empty()) {
620 dbSpecies = critPropsDb[
"species"].asMap(
"name");
624 auto ucName = boost::algorithm::to_upper_copy(
name);
625 if (dbSpecies.count(ucName)) {
626 auto& spec = *dbSpecies.at(ucName);
627 auto& critProps = spec[
"critical-parameters"].as<
AnyMap>();
628 Tc = critProps.
convert(
"critical-temperature",
"K");
629 Pc = critProps.convert(
"critical-pressure",
"Pa");
642 "No critical property or Redlich-Kwong parameters found "
643 "for species {}.",
name);
649 AnyMap& speciesNode)
const
654 auto& eosNode = speciesNode[
"equation-of-state"].getMapWhere(
655 "model",
"Redlich-Kwong",
true);
657 size_t counter = k +
m_kk * k;
658 if (a_coeff_vec(1, counter) != 0.0) {
659 vector<AnyValue> coeffs(2);
660 coeffs[0].setQuantity(a_coeff_vec(0, counter),
"Pa*m^6/kmol^2*K^0.5");
661 coeffs[1].setQuantity(a_coeff_vec(1, counter),
"Pa*m^6/kmol^2/K^0.5");
662 eosNode[
"a"] = std::move(coeffs);
664 eosNode[
"a"].setQuantity(a_coeff_vec(0, counter),
665 "Pa*m^6/kmol^2*K^0.5");
667 eosNode[
"b"].setQuantity(b_vec_Curr_[k],
"m^3/kmol");
669 auto& critProps = speciesNode[
"critical-parameters"];
670 double a = a_coeff_vec(0, k +
m_kk * k);
671 double b = b_vec_Curr_[k];
674 critProps[
"critical-temperature"].setQuantity(Tc,
"K");
675 critProps[
"critical-pressure"].setQuantity(Pc,
"Pa");
679 auto& eosNode = speciesNode[
"equation-of-state"].getMapWhere(
680 "model",
"Redlich-Kwong",
true);
683 if (coeffs.second == 0) {
684 bin_a[name2].setQuantity(coeffs.first,
"Pa*m^6/kmol^2*K^0.5");
686 vector<AnyValue> C(2);
687 C[0].setQuantity(coeffs.first,
"Pa*m^6/kmol^2*K^0.5");
688 C[1].setQuantity(coeffs.second,
"Pa*m^6/kmol^2/K^0.5");
689 bin_a[name2] = std::move(C);
692 eosNode[
"binary-a"] = std::move(bin_a);
701 double molarV = mmw / rho;
704 double dadt = da_dt();
706 double sqT = sqrt(T);
717 double molarV = mmw / rho;
720 double dadt = da_dt();
722 double sqT = sqrt(T);
733 double pres = std::max(
psatEst(TKelvin), presGuess);
735 bool foundLiq =
false;
737 while (m < 100 && !foundLiq) {
738 int nsol =
solveCubic(TKelvin, pres, atmp, btmp, Vroot);
739 if (nsol == 1 || nsol == 2) {
766 if (rhoguess == -1.0) {
767 if (phaseRequested >= FLUID_LIQUID_0) {
769 rhoguess = mmw / lqvol;
777 double volguess = mmw / rhoguess;
780 double molarVolLast = Vroot_[0];
782 if (phaseRequested >= FLUID_LIQUID_0) {
783 molarVolLast = Vroot_[0];
784 }
else if (phaseRequested == FLUID_GAS || phaseRequested == FLUID_SUPERCRIT) {
785 molarVolLast = Vroot_[2];
787 if (volguess > Vroot_[1]) {
788 molarVolLast = Vroot_[2];
790 molarVolLast = Vroot_[0];
793 }
else if (NSolns_ == 1) {
794 if (phaseRequested == FLUID_GAS || phaseRequested == FLUID_SUPERCRIT || phaseRequested == FLUID_UNDEFINED) {
795 molarVolLast = Vroot_[0];
799 }
else if (NSolns_ == -1) {
800 if (phaseRequested >= FLUID_LIQUID_0 || phaseRequested == FLUID_UNDEFINED || phaseRequested == FLUID_SUPERCRIT) {
801 molarVolLast = Vroot_[0];
802 }
else if (TKelvin > tcrit) {
803 molarVolLast = Vroot_[0];
808 molarVolLast = Vroot_[0];
811 return mmw / molarVolLast;
816 double sqt = sqrt(TKelvin);
822 double dpdv = (-
GasConstant * TKelvin / (vmb * vmb)
858 double sqt = sqrt(TKelvin);
861 double dadt = da_dt();
870 for (
size_t i = 0; i <
m_kk; i++) {
871 for (
size_t j = 0; j <
m_kk; j++) {
872 size_t counter = i *
m_kk + j;
873 a_vec_Curr_[counter] = a_coeff_vec(0,counter) + a_coeff_vec(1,counter) * temp;
880 for (
size_t i = 0; i <
m_kk; i++) {
882 for (
size_t j = 0; j <
m_kk; j++) {
888 fmt::memory_buffer b;
889 for (
size_t k = 0; k <
m_kk; k++) {
890 if (isnan(b_vec_Curr_[k])) {
898 throw CanteraError(
"RedlichKwongMFTP::updateMixingExpressions",
899 "Missing Redlich-Kwong coefficients for species: {}", to_string(b));
908 for (
size_t i = 0; i <
m_kk; i++) {
910 for (
size_t j = 0; j <
m_kk; j++) {
911 size_t counter = i *
m_kk + j;
912 double a_vec_Curr = a_coeff_vec(0,counter) + a_coeff_vec(1,counter) * temp;
917 for (
size_t i = 0; i <
m_kk; i++) {
919 for (
size_t j = 0; j <
m_kk; j++) {
920 size_t counter = i *
m_kk + j;
921 double a_vec_Curr = a_coeff_vec(0,counter);
928double RedlichKwongMFTP::da_dt()
const
932 for (
size_t i = 0; i <
m_kk; i++) {
933 for (
size_t j = 0; j <
m_kk; j++) {
934 size_t counter = i *
m_kk + j;
942void RedlichKwongMFTP::calcCriticalConditions(
double& pc,
double& tc,
double& vc)
const
946 for (
size_t i = 0; i <
m_kk; i++) {
947 for (
size_t j = 0; j <
m_kk; j++) {
948 size_t counter = i +
m_kk * j;
972 double sqrttc, f, dfdt, deltatc;
978 for (
int j = 0; j < 10; j++) {
982 deltatc = - f / dfdt;
986 throw CanteraError(
"RedlichKwongMFTP::calcCriticalConditions",
996 span<double> Vroot)
const
1002 double sqt = sqrt(T);
1003 double cn = - (
GasConstant * T * b / pres - a/(pres * sqt) + b * b);
1004 double dn = - (a * b / (pres * sqt));
1008 double tc = pow(tmp, pp);
1012 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 getGibbs_ref(span< double > g) const override
Returns the vector of the Gibbs function of the reference state at the current temperature of the sol...
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.
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.
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...