12#include <boost/algorithm/string.hpp>
34 m_kappa[k] = 0.37464 + 1.54226*w - 0.26992*w*w;
36 m_kappa[k] = 0.374642 + 1.487503*w - 0.164423*w*w + 0.016666*w*w*w;
43 double sqt_alpha = 1 + m_kappa[k] * (1 - sqt_T_r);
44 m_alpha[k] = sqt_alpha*sqt_alpha;
47 double aAlpha_k = a*m_alpha[k];
48 m_aAlpha_binary(k,k) = aAlpha_k;
51 for (
size_t j = 0; j <
m_kk; j++) {
55 double a0kj = sqrt(m_a_coeffs(j,j) * a);
56 double aAlpha_j = a*m_alpha[j];
57 double a_Alpha = sqrt(aAlpha_j*aAlpha_k);
58 if (m_a_coeffs(j, k) == 0) {
59 m_a_coeffs(j, k) = a0kj;
60 m_aAlpha_binary(j, k) = a_Alpha;
61 m_a_coeffs(k, j) = a0kj;
62 m_aAlpha_binary(k, j) = a_Alpha;
69 const string& species_j,
double a0)
74 m_a_coeffs(ki, kj) = m_a_coeffs(kj, ki) = a0;
78 double alpha_ij = m_alpha[ki] * m_alpha[kj];
79 m_aAlpha_binary(ki, kj) = m_aAlpha_binary(kj, ki) = a0*alpha_ij;
119 double denom = mv * mv + 2 * mv *
m_b -
m_b *
m_b;
132 for (
size_t k = 0; k <
m_kk; k++) {
139 double vpb2 = mv + (1 +
Sqrt2) *
m_b;
140 double vmb2 = mv + (1 -
Sqrt2) *
m_b;
141 double vmb = mv -
m_b;
144 for (
size_t k = 0; k <
m_kk; k++) {
146 for (
size_t i = 0; i <
m_kk; i++) {
152 double denom2 =
m_b * (mv * mv + 2 * mv *
m_b -
m_b *
m_b);
154 for (
size_t k = 0; k <
m_kk; k++) {
156 ac[k] = (-RT_ * log(pres * mv/ RT_) + RT_ * log(mv / vmb)
157 + RT_ * m_b_coeffs[k] / vmb
158 - (num /denom) * log(vpb2/vmb2)
162 for (
size_t k = 0; k <
m_kk; k++) {
163 ac[k] = exp(ac[k]/ RT_);
173 for (
size_t k = 0; k <
m_kk; k++) {
175 mu[k] += RT_ * log(xx);
182 double vmb = mv -
m_b;
183 double vpb2 = mv + (1 +
Sqrt2) *
m_b;
184 double vmb2 = mv + (1 -
Sqrt2) *
m_b;
186 for (
size_t k = 0; k <
m_kk; k++) {
188 for (
size_t i = 0; i <
m_kk; i++) {
194 double denom2 =
m_b * (mv * mv + 2 * mv *
m_b -
m_b *
m_b);
196 for (
size_t k = 0; k <
m_kk; k++) {
199 mu[k] += - RT_ * log(pres * mv / RT_)
200 + RT_ * log(mv / vmb) + RT_ * m_b_coeffs[k] / vmb
201 - (num /denom) * log(vpb2/vmb2)
210 scale(hbar.begin(), hbar.end(), hbar.begin(),
RT());
215 tmp.resize(
m_kk,0.0);
220 double vmb = mv -
m_b;
221 double vpb2 = mv + (1 +
Sqrt2) *
m_b;
222 double vmb2 = mv + (1 -
Sqrt2) *
m_b;
225 for (
size_t k = 0; k <
m_kk; k++) {
228 for (
size_t i = 0; i <
m_kk; i++) {
229 double grad_aAlpha = m_dalphadT[i]/m_alpha[i] + m_dalphadT[k]/m_alpha[k];
235 double denom = mv * mv + 2 * mv *
m_b -
m_b *
m_b;
236 double denom2 = denom * denom;
238 for (
size_t k = 0; k <
m_kk; k++) {
239 m_dpdni[k] = RT_ / vmb + RT_ * m_b_coeffs[k] / (vmb * vmb)
248 for (
size_t k = 0; k <
m_kk; k++) {
249 fac4 = T*tmp[k] -2 *
m_pp[k];
250 double hE_v = mv *
m_dpdni[k] - RT_
251 - m_b_coeffs[k] / fac3 * log(vpb2 / vmb2) * fac
252 + (mv * m_b_coeffs[k]) / (
m_b * denom) * fac
253 + 1 / (2 *
Sqrt2 *
m_b) * log(vpb2 / vmb2) * fac4;
254 hbar[k] = hbar[k] + hE_v;
265 for (
size_t k = 0; k <
m_kk; k++) {
266 sbar[k] = (sbar[k] -
m_workS[k])/T;
276 for (
size_t k = 0; k <
m_kk; k++) {
277 ubar[k] = ubar[k] - p*
m_workS[k];
289 for (
size_t k = 0; k <
m_kk; k++) {
291 for (
size_t i = 0; i <
m_kk; i++) {
297 double vmb = mv -
m_b;
298 double vpb = mv +
m_b;
299 double fac = mv * mv + 2 * mv *
m_b -
m_b *
m_b;
300 double fac2 = fac * fac;
303 for (
size_t k = 0; k <
m_kk; k++) {
304 double num = RT_ + RT_ *
m_b/ vmb + RT_ * m_b_coeffs[k] / vmb
305 + RT_ *
m_b * m_b_coeffs[k] /(vmb * vmb) - 2 * mv *
m_pp[k] / fac
309 vbar[k] = num / denom;
317 }
else if (a <= 0.0) {
329 m_b_coeffs.push_back(0.0);
331 m_kappa.push_back(0.0);
333 m_alpha.push_back(0.0);
334 m_dalphadT.push_back(0.0);
335 m_d2alphadT2.push_back(0.0);
337 m_partialMolarVolumes.push_back(0.0);
348 std::unordered_map<string, AnyMap*> dbSpecies;
353 if (m_a_coeffs(k, k) != 0.0) {
356 bool foundCoeffs =
false;
357 if (data.hasKey(
"equation-of-state") &&
358 data[
"equation-of-state"].hasMapWhere(
"model",
"Peng-Robinson"))
362 auto eos = data[
"equation-of-state"].getMapWhere(
363 "model",
"Peng-Robinson");
364 if (eos.hasKey(
"a") && eos.hasKey(
"b") && eos.hasKey(
"acentric-factor")) {
365 double a0 = eos.convert(
"a",
"Pa*m^6/kmol^2");
366 double b = eos.convert(
"b",
"m^3/kmol");
368 double w = eos[
"acentric-factor"].asDouble();
373 if (eos.hasKey(
"binary-a")) {
376 for (
auto& [name2, coeff] : binary_a) {
377 double a0 = units.
convert(coeff,
"Pa*m^6/kmol^2");
388 double Tc = NAN, Pc = NAN, omega_ac = NAN;
389 if (data.hasKey(
"critical-parameters")) {
392 auto& critProps = data[
"critical-parameters"].as<
AnyMap>();
393 Tc = critProps.
convert(
"critical-temperature",
"K");
394 Pc = critProps.convert(
"critical-pressure",
"Pa");
395 omega_ac = critProps[
"acentric-factor"].asDouble();
399 if (critPropsDb.
empty()) {
401 dbSpecies = critPropsDb[
"species"].asMap(
"name");
405 auto ucName = boost::algorithm::to_upper_copy(
name);
406 if (dbSpecies.count(ucName)) {
407 auto& spec = *dbSpecies.at(ucName);
408 auto& critProps = spec[
"critical-parameters"].as<
AnyMap>();
409 Tc = critProps.
convert(
"critical-temperature",
"K");
410 Pc = critProps.convert(
"critical-pressure",
"Pa");
411 omega_ac = critProps[
"acentric-factor"].asDouble();
423 "No Peng-Robinson model parameters or critical properties found for "
424 "species '{}'",
name);
436 auto& eosNode = speciesNode[
"equation-of-state"].getMapWhere(
437 "model",
"Peng-Robinson",
true);
438 eosNode[
"a"].setQuantity(m_a_coeffs(k, k),
"Pa*m^6/kmol^2");
439 eosNode[
"b"].setQuantity(m_b_coeffs[k],
"m^3/kmol");
442 auto& critProps = speciesNode[
"critical-parameters"];
445 critProps[
"critical-temperature"].setQuantity(Tc,
"K");
446 critProps[
"critical-pressure"].setQuantity(Pc,
"Pa");
454 auto& eosNode = speciesNode[
"equation-of-state"].getMapWhere(
455 "model",
"Peng-Robinson",
true);
458 bin_a[
species].setQuantity(coeff,
"Pa*m^6/kmol^2");
460 eosNode[
"binary-a"] = std::move(bin_a);
470 double hh =
m_b / molarV;
473 double vpb = molarV + (1.0 +
Sqrt2) *
m_b;
474 double vmb = molarV + (1.0 -
Sqrt2) *
m_b;
475 double fac = alpha_1 / (2.0 *
Sqrt2 *
m_b);
476 double sresid_mol_R = log(zz*(1.0 - hh)) + fac * log(vpb / vmb) /
GasConstant;
489 double vpb = molarV + (1 +
Sqrt2) *
m_b;
490 double vmb = molarV + (1 -
Sqrt2) *
m_b;
491 double fac = 1 / (2.0 *
Sqrt2 *
m_b);
498 double v =
m_b * 1.1;
499 double atmp, btmp, aAlphatmp;
501 double pres = std::max(
psatEst(T), presGuess);
503 bool foundLiq =
false;
505 while (m < 100 && !foundLiq) {
506 int nsol =
solveCubic(T, pres, atmp, btmp, aAlphatmp, Vroot);
507 if (nsol == 1 || nsol == 2) {
534 if (rhoGuess == -1.0) {
535 if (phaseRequested >= FLUID_LIQUID_0) {
537 rhoGuess = mmw / lqvol;
544 double volGuess = mmw / rhoGuess;
549 double vmin = std::max(0.0,
m_b * (1.0 + 1e-12));
550 vector<double> physicalRoots;
551 for (
double root : m_Vroot) {
552 if (std::isfinite(
root) &&
root > vmin) {
553 physicalRoots.push_back(
root);
556 std::sort(physicalRoots.begin(), physicalRoots.end());
557 if (physicalRoots.empty()) {
559 }
else if (physicalRoots.size() == 1) {
563 if ((phaseRequested == FLUID_GAS && m_NSolns < 0)
564 || (phaseRequested >= FLUID_LIQUID_0 && m_NSolns > 0))
569 return mmw / physicalRoots[0];
570 }
else if (physicalRoots.size() == 2) {
571 m_Vroot[0] = physicalRoots[0];
572 m_Vroot[1] = 0.5 * (physicalRoots[0] + physicalRoots[1]);
573 m_Vroot[2] = physicalRoots[1];
575 m_Vroot[0] = physicalRoots[0];
576 m_Vroot[1] = physicalRoots[1];
577 m_Vroot[2] = physicalRoots[2];
580 double molarVolLast = m_Vroot[0];
582 if (phaseRequested >= FLUID_LIQUID_0) {
583 molarVolLast = m_Vroot[0];
584 }
else if (phaseRequested == FLUID_GAS || phaseRequested == FLUID_SUPERCRIT) {
585 molarVolLast = m_Vroot[2];
587 if (volGuess > m_Vroot[1]) {
588 molarVolLast = m_Vroot[2];
590 molarVolLast = m_Vroot[0];
593 }
else if (m_NSolns == 1) {
594 if (phaseRequested == FLUID_GAS || phaseRequested == FLUID_SUPERCRIT
595 || phaseRequested == FLUID_UNDEFINED)
597 molarVolLast = m_Vroot[0];
601 }
else if (m_NSolns == -1) {
602 if (phaseRequested >= FLUID_LIQUID_0 || phaseRequested == FLUID_UNDEFINED
603 || phaseRequested == FLUID_SUPERCRIT)
605 molarVolLast = m_Vroot[0];
606 }
else if (T > tcrit) {
607 molarVolLast = m_Vroot[0];
612 molarVolLast = m_Vroot[0];
615 return mmw / molarVolLast;
620 double denom = molarVol * molarVol + 2 * molarVol *
m_b -
m_b *
m_b;
621 double vpb = molarVol +
m_b;
622 double vmb = molarVol -
m_b;
651 double vmb = mv -
m_b;
652 double denom = mv * mv + 2 * mv *
m_b -
m_b *
m_b;
661 for (
size_t j = 0; j <
m_kk; j++) {
663 double sqt_alpha = 1 + m_kappa[j] * (1 - sqrt(temp / critTemp_j));
664 m_alpha[j] = sqt_alpha*sqt_alpha;
668 for (
size_t i = 0; i <
m_kk; i++) {
669 for (
size_t j = 0; j <
m_kk; j++) {
670 m_aAlpha_binary(i, j) = sqrt(m_alpha[i] * m_alpha[j]) * m_a_coeffs(i,j);
681 for (
size_t i = 0; i <
m_kk; i++) {
683 for (
size_t j = 0; j <
m_kk; j++) {
692 double daAlphadT = 0.0, k, Tc, sqtTr, coeff1, coeff2;
693 for (
size_t i = 0; i <
m_kk; i++) {
697 coeff1 = 1 / (Tc*sqtTr);
700 m_dalphadT[i] = coeff1 * (k*k*coeff2 - k);
703 for (
size_t i = 0; i <
m_kk; i++) {
704 for (
size_t j = 0; j <
m_kk; j++) {
706 * m_aAlpha_binary(i, j)
707 * (m_dalphadT[i] / m_alpha[i] + m_dalphadT[j] / m_alpha[j]);
715 for (
size_t i = 0; i <
m_kk; i++) {
718 double coeff1 = 1 / (Tcrit_i*sqt_Tr);
719 double coeff2 = sqt_Tr - 1;
721 double k = m_kappa[i];
722 m_dalphadT[i] = coeff1 * (k*k*coeff2 - k);
723 m_d2alphadT2[i] = (k*k + k) * coeff1 / (2*sqt_Tr*sqt_Tr*Tcrit_i);
727 double d2aAlphadT2 = 0.0;
728 for (
size_t i = 0; i <
m_kk; i++) {
729 double alphai = m_alpha[i];
730 for (
size_t j = 0; j <
m_kk; j++) {
731 double alphaj = m_alpha[j];
732 double alphaij = alphai * alphaj;
733 double term1 = m_d2alphadT2[i] / alphai + m_d2alphadT2[j] / alphaj;
734 double term2 = 2 * m_dalphadT[i] * m_dalphadT[j] / alphaij;
735 double term3 = m_dalphadT[i] / alphai + m_dalphadT[j] / alphaj;
737 * m_aAlpha_binary(i, j)
738 * (term1 + term2 - 0.5 * term3 * term3);
744void PengRobinson::calcCriticalConditions(
double& pc,
double& tc,
double& vc)
const
764 span<double> Vroot)
const
769 double aAlpha_p = aAlpha / pres;
771 double bn = (b - RT_p);
772 double cn = -(2 * RT_p * b - aAlpha_p + 3 * bsqr);
773 double dn = (bsqr * RT_p + bsqr * b - aAlpha_p * b);
780 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.
virtual void resize(size_t n, size_t m, double v=0.0)
Resize the array, and fill the new entries with 'v'.
double critPressure() const override
Critical pressure (Pa).
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.
An error indicating that an unimplemented function has been called.
double thermalExpansionCoeff() const override
Return the volumetric thermal expansion coefficient. Units: 1/K.
void setSpeciesCoeffs(const string &species, double a, double b, double w)
Set the pure fluid interaction parameters for a species.
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
Returns an array of partial molar enthalpies for the species in the mixture.
double soundSpeed() const override
Return the speed of sound. Units: m/s.
double pressure() const override
Return the thermodynamic pressure (Pa).
map< string, map< string, double > > m_binaryParameters
Explicitly-specified binary interaction parameters, to enable serialization.
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...
void getPartialMolarCp(span< double > cpbar) const override
Calculate species-specific molar specific heats.
static const double omega_b
Omega constant: b0 (= omega_b) used in Peng-Robinson equation of state.
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 setBinaryCoeffs(const string &species_i, const string &species_j, double a)
Set values for the interaction parameter between two species.
vector< double > m_dpdni
Vector of derivatives of pressure with respect to mole number.
PengRobinson(const string &infile="", const string &id="")
Construct and initialize a PengRobinson object directly from an input file.
void initThermo() override
Initialize the ThermoPhase object after all species have been set up.
double m_dpdT
The derivative of the pressure with respect to the temperature.
double cv_mole() const override
Molar heat capacity at constant volume and composition [J/kmol/K].
double liquidVolEst(double T, double &pres) const override
Estimate for the molar volume of the liquid.
double daAlpha_dT() const
Calculate temperature derivative .
void calculatePressureDerivatives() const
Calculate and 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.
double m_aAlpha_mix
Value of in the equation of state.
double m_a
Value of in the equation of state.
static const double omega_a
Omega constant: a0 (= omega_a) used in Peng-Robinson equation of state.
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, b, and omega coefficients.
void calculateAB(double &aCalc, double &bCalc, double &aAlpha) const
Calculate the , , and parameters given the temperature.
void updateMixingExpressions() override
Update the , , and parameters.
double cp_mole() const override
Molar heat capacity at constant pressure and composition [J/kmol/K].
double m_b
Value of in the equation of state.
void getPartialMolarVolumes(span< double > vbar) const override
Return an array of partial molar volumes for the species in the mixture.
int solveCubic(double T, double pres, double a, double b, double aAlpha, span< double > Vroot) const
Prepare variables and call the function to solve the cubic 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.
AnyMap getAuxiliaryData() override
Return parameters used by the Peng-Robinson equation of state.
bool addSpecies(shared_ptr< Species > spec) override
Add a Species to this Phase.
double dpdVCalc(double T, double molarVol, double &presCalc) const override
Calculate the pressure and the pressure derivative given the temperature and the molar volume.
void getChemPotentials(span< double > mu) const override
Get the species chemical potentials. Units: J/kmol.
double d2aAlpha_dT2() const
Calculate second derivative .
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.
double m_dpdV
The derivative of the pressure with respect to the volume.
double speciesCritTemperature(double a, double b) const
Calculate species-specific critical temperature.
void getPartialMolarIntEnergies(span< double > ubar) const override
Return an array of partial molar internal energies for the species in the mixture.
vector< double > m_acentric
acentric factor for each species, length m_kk
vector< double > m_workS
Vector of size m_kk, used as a temporary holding area.
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)
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.
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 RT() const
Return the Gas Constant multiplied by the current temperature.
AnyMap parameters(bool withInput=true) const
Returns the parameters of a ThermoPhase object such that an identical one could be reconstructed usin...
virtual AnyMap getAuxiliaryData()
Return intermediate or model-specific parameters used by particular derived classes.
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...
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].
const double Sqrt2
Sqrt(2)
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...