12#include <boost/algorithm/string.hpp>
33 m_kappa[k] = 0.37464 + 1.54226*w - 0.26992*w*w;
35 m_kappa[k] = 0.374642 + 1.487503*w - 0.164423*w*w + 0.016666*w*w*w;
42 double sqt_alpha = 1 + m_kappa[k] * (1 - sqt_T_r);
43 m_alpha[k] = sqt_alpha*sqt_alpha;
46 double aAlpha_k = a*m_alpha[k];
47 m_aAlpha_binary(k,k) = aAlpha_k;
50 for (
size_t j = 0; j <
m_kk; j++) {
54 double a0kj = sqrt(m_a_coeffs(j,j) * a);
55 double aAlpha_j = a*m_alpha[j];
56 double a_Alpha = sqrt(aAlpha_j*aAlpha_k);
57 if (m_a_coeffs(j, k) == 0) {
58 m_a_coeffs(j, k) = a0kj;
59 m_aAlpha_binary(j, k) = a_Alpha;
60 m_a_coeffs(k, j) = a0kj;
61 m_aAlpha_binary(k, j) = a_Alpha;
68 const string& species_j,
double a0)
73 m_a_coeffs(ki, kj) = m_a_coeffs(kj, ki) = a0;
77 double alpha_ij = m_alpha[ki] * m_alpha[kj];
78 m_aAlpha_binary(ki, kj) = m_aAlpha_binary(kj, ki) = a0*alpha_ij;
111 double denom = mv * mv + 2 * mv *
m_b -
m_b *
m_b;
125 double vpb2 = mv + (1 +
Sqrt2) *
m_b;
126 double vmb2 = mv + (1 -
Sqrt2) *
m_b;
127 double vmb = mv -
m_b;
130 for (
size_t k = 0; k <
m_kk; k++) {
132 for (
size_t i = 0; i <
m_kk; i++) {
138 double denom2 =
m_b * (mv * mv + 2 * mv *
m_b -
m_b *
m_b);
140 for (
size_t k = 0; k <
m_kk; k++) {
142 ac[k] = (-RT_ * log(pres * mv/ RT_) + RT_ * log(mv / vmb)
143 + RT_ * m_b_coeffs[k] / vmb
144 - (num /denom) * log(vpb2/vmb2)
148 for (
size_t k = 0; k <
m_kk; k++) {
149 ac[k] = exp(ac[k]/ RT_);
159 for (
size_t k = 0; k <
m_kk; k++) {
161 mu[k] += RT_ * (log(xx));
165 double vmb = mv -
m_b;
166 double vpb2 = mv + (1 +
Sqrt2) *
m_b;
167 double vmb2 = mv + (1 -
Sqrt2) *
m_b;
169 for (
size_t k = 0; k <
m_kk; k++) {
171 for (
size_t i = 0; i <
m_kk; i++) {
178 double denom2 =
m_b * (mv * mv + 2 * mv *
m_b -
m_b *
m_b);
180 for (
size_t k = 0; k <
m_kk; k++) {
183 mu[k] += RT_ * log(pres/refP) - RT_ * log(pres * mv / RT_)
184 + RT_ * log(mv / vmb) + RT_ * m_b_coeffs[k] / vmb
185 - (num /denom) * log(vpb2/vmb2)
194 scale(hbar.begin(), hbar.end(), hbar.begin(),
RT());
196 tmp.resize(
m_kk,0.0);
201 double vmb = mv -
m_b;
202 double vpb2 = mv + (1 +
Sqrt2) *
m_b;
203 double vmb2 = mv + (1 -
Sqrt2) *
m_b;
206 for (
size_t k = 0; k <
m_kk; k++) {
209 for (
size_t i = 0; i <
m_kk; i++) {
210 double grad_aAlpha = m_dalphadT[i]/m_alpha[i] + m_dalphadT[k]/m_alpha[k];
216 double denom = mv * mv + 2 * mv *
m_b -
m_b *
m_b;
217 double denom2 = denom * denom;
219 for (
size_t k = 0; k <
m_kk; k++) {
220 m_dpdni[k] = RT_ / vmb + RT_ * m_b_coeffs[k] / (vmb * vmb)
229 for (
size_t k = 0; k <
m_kk; k++) {
230 fac4 = T*tmp[k] -2 *
m_pp[k];
231 double hE_v = mv *
m_dpdni[k] - RT_
232 - m_b_coeffs[k] / fac3 * log(vpb2 / vmb2) * fac
233 + (mv * m_b_coeffs[k]) / (
m_b * denom) * fac
234 + 1 / (2 *
Sqrt2 *
m_b) * log(vpb2 / vmb2) * fac4;
235 hbar[k] = hbar[k] + hE_v;
246 for (
size_t k = 0; k <
m_kk; k++) {
247 sbar[k] = (sbar[k] -
m_workS[k])/T;
257 for (
size_t k = 0; k <
m_kk; k++) {
258 ubar[k] = ubar[k] - p*
m_workS[k];
270 for (
size_t k = 0; k <
m_kk; k++) {
272 for (
size_t i = 0; i <
m_kk; i++) {
278 double vmb = mv -
m_b;
279 double vpb = mv +
m_b;
280 double fac = mv * mv + 2 * mv *
m_b -
m_b *
m_b;
281 double fac2 = fac * fac;
284 for (
size_t k = 0; k <
m_kk; k++) {
285 double num = RT_ + RT_ *
m_b/ vmb + RT_ * m_b_coeffs[k] / vmb
286 + RT_ *
m_b * m_b_coeffs[k] /(vmb * vmb) - 2 * mv *
m_pp[k] / fac
290 vbar[k] = num / denom;
298 }
else if (a <= 0.0) {
310 m_b_coeffs.push_back(0.0);
312 m_kappa.push_back(0.0);
314 m_alpha.push_back(0.0);
315 m_dalphadT.push_back(0.0);
316 m_d2alphadT2.push_back(0.0);
318 m_partialMolarVolumes.push_back(0.0);
329 std::unordered_map<string, AnyMap*> dbSpecies;
334 if (m_a_coeffs(k, k) != 0.0) {
337 bool foundCoeffs =
false;
338 if (data.hasKey(
"equation-of-state") &&
339 data[
"equation-of-state"].hasMapWhere(
"model",
"Peng-Robinson"))
343 auto eos = data[
"equation-of-state"].getMapWhere(
344 "model",
"Peng-Robinson");
345 if (eos.hasKey(
"a") && eos.hasKey(
"b") && eos.hasKey(
"acentric-factor")) {
346 double a0 = eos.convert(
"a",
"Pa*m^6/kmol^2");
347 double b = eos.convert(
"b",
"m^3/kmol");
349 double w = eos[
"acentric-factor"].asDouble();
354 if (eos.hasKey(
"binary-a")) {
357 for (
auto& [name2, coeff] : binary_a) {
358 double a0 = units.
convert(coeff,
"Pa*m^6/kmol^2");
369 double Tc = NAN, Pc = NAN, omega_ac = NAN;
370 if (data.hasKey(
"critical-parameters")) {
373 auto& critProps = data[
"critical-parameters"].as<
AnyMap>();
374 Tc = critProps.
convert(
"critical-temperature",
"K");
375 Pc = critProps.convert(
"critical-pressure",
"Pa");
376 omega_ac = critProps[
"acentric-factor"].asDouble();
380 if (critPropsDb.
empty()) {
382 dbSpecies = critPropsDb[
"species"].asMap(
"name");
386 auto ucName = boost::algorithm::to_upper_copy(
name);
387 if (dbSpecies.count(ucName)) {
388 auto& spec = *dbSpecies.at(ucName);
389 auto& critProps = spec[
"critical-parameters"].as<
AnyMap>();
390 Tc = critProps.
convert(
"critical-temperature",
"K");
391 Pc = critProps.convert(
"critical-pressure",
"Pa");
392 omega_ac = critProps[
"acentric-factor"].asDouble();
404 "No Peng-Robinson model parameters or critical properties found for "
405 "species '{}'",
name);
417 auto& eosNode = speciesNode[
"equation-of-state"].getMapWhere(
418 "model",
"Peng-Robinson",
true);
419 eosNode[
"a"].setQuantity(m_a_coeffs(k, k),
"Pa*m^6/kmol^2");
420 eosNode[
"b"].setQuantity(m_b_coeffs[k],
"m^3/kmol");
423 auto& critProps = speciesNode[
"critical-parameters"];
426 critProps[
"critical-temperature"].setQuantity(Tc,
"K");
427 critProps[
"critical-pressure"].setQuantity(Pc,
"Pa");
435 auto& eosNode = speciesNode[
"equation-of-state"].getMapWhere(
436 "model",
"Peng-Robinson",
true);
439 bin_a[
species].setQuantity(coeff,
"Pa*m^6/kmol^2");
441 eosNode[
"binary-a"] = std::move(bin_a);
448 double hh =
m_b / molarV;
451 double vpb = molarV + (1.0 +
Sqrt2) *
m_b;
452 double vmb = molarV + (1.0 -
Sqrt2) *
m_b;
453 double fac = alpha_1 / (2.0 *
Sqrt2 *
m_b);
454 double sresid_mol_R = log(zz*(1.0 - hh)) + fac * log(vpb / vmb) /
GasConstant;
464 double vpb = molarV + (1 +
Sqrt2) *
m_b;
465 double vmb = molarV + (1 -
Sqrt2) *
m_b;
466 double fac = 1 / (2.0 *
Sqrt2 *
m_b);
473 double v =
m_b * 1.1;
474 double atmp, btmp, aAlphatmp;
476 double pres = std::max(
psatEst(T), presGuess);
478 bool foundLiq =
false;
480 while (m < 100 && !foundLiq) {
481 int nsol =
solveCubic(T, pres, atmp, btmp, aAlphatmp, Vroot);
482 if (nsol == 1 || nsol == 2) {
509 if (rhoGuess == -1.0) {
510 if (phaseRequested >= FLUID_LIQUID_0) {
512 rhoGuess = mmw / lqvol;
519 double volGuess = mmw / rhoGuess;
522 double molarVolLast = m_Vroot[0];
524 if (phaseRequested >= FLUID_LIQUID_0) {
525 molarVolLast = m_Vroot[0];
526 }
else if (phaseRequested == FLUID_GAS || phaseRequested == FLUID_SUPERCRIT) {
527 molarVolLast = m_Vroot[2];
529 if (volGuess > m_Vroot[1]) {
530 molarVolLast = m_Vroot[2];
532 molarVolLast = m_Vroot[0];
535 }
else if (m_NSolns == 1) {
536 if (phaseRequested == FLUID_GAS || phaseRequested == FLUID_SUPERCRIT
537 || phaseRequested == FLUID_UNDEFINED)
539 molarVolLast = m_Vroot[0];
543 }
else if (m_NSolns == -1) {
544 if (phaseRequested >= FLUID_LIQUID_0 || phaseRequested == FLUID_UNDEFINED
545 || phaseRequested == FLUID_SUPERCRIT)
547 molarVolLast = m_Vroot[0];
548 }
else if (T > tcrit) {
549 molarVolLast = m_Vroot[0];
554 molarVolLast = m_Vroot[0];
557 return mmw / molarVolLast;
562 double denom = molarVol * molarVol + 2 * molarVol *
m_b -
m_b *
m_b;
563 double vpb = molarVol +
m_b;
564 double vmb = molarVol -
m_b;
593 double vmb = mv -
m_b;
594 double denom = mv * mv + 2 * mv *
m_b -
m_b *
m_b;
603 for (
size_t j = 0; j <
m_kk; j++) {
605 double sqt_alpha = 1 + m_kappa[j] * (1 - sqrt(temp / critTemp_j));
606 m_alpha[j] = sqt_alpha*sqt_alpha;
610 for (
size_t i = 0; i <
m_kk; i++) {
611 for (
size_t j = 0; j <
m_kk; j++) {
612 m_aAlpha_binary(i, j) = sqrt(m_alpha[i] * m_alpha[j]) * m_a_coeffs(i,j);
623 for (
size_t i = 0; i <
m_kk; i++) {
625 for (
size_t j = 0; j <
m_kk; j++) {
634 double daAlphadT = 0.0, k, Tc, sqtTr, coeff1, coeff2;
635 for (
size_t i = 0; i <
m_kk; i++) {
639 coeff1 = 1 / (Tc*sqtTr);
642 m_dalphadT[i] = coeff1 * (k*k*coeff2 - k);
645 for (
size_t i = 0; i <
m_kk; i++) {
646 for (
size_t j = 0; j <
m_kk; j++) {
648 * m_aAlpha_binary(i, j)
649 * (m_dalphadT[i] / m_alpha[i] + m_dalphadT[j] / m_alpha[j]);
657 for (
size_t i = 0; i <
m_kk; i++) {
660 double coeff1 = 1 / (Tcrit_i*sqt_Tr);
661 double coeff2 = sqt_Tr - 1;
663 double k = m_kappa[i];
664 m_dalphadT[i] = coeff1 * (k*k*coeff2 - k);
665 m_d2alphadT2[i] = (k*k + k) * coeff1 / (2*sqt_Tr*sqt_Tr*Tcrit_i);
669 double d2aAlphadT2 = 0.0;
670 for (
size_t i = 0; i <
m_kk; i++) {
671 double alphai = m_alpha[i];
672 for (
size_t j = 0; j <
m_kk; j++) {
673 double alphaj = m_alpha[j];
674 double alphaij = alphai * alphaj;
675 double term1 = m_d2alphadT2[i] / alphai + m_d2alphadT2[j] / alphaj;
676 double term2 = 2 * m_dalphadT[i] * m_dalphadT[j] / alphaij;
677 double term3 = m_dalphadT[i] / alphai + m_dalphadT[j] / alphaj;
679 * m_aAlpha_binary(i, j)
680 * (term1 + term2 - 0.5 * term3 * term3);
686void PengRobinson::calcCriticalConditions(
double& pc,
double& tc,
double& vc)
const
706 span<double> Vroot)
const
711 double aAlpha_p = aAlpha / pres;
713 double bn = (b - RT_p);
714 double cn = -(2 * RT_p * b - aAlpha_p + 3 * bsqr);
715 double dn = (bsqr * RT_p + bsqr * b - aAlpha_p * b);
722 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'.
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...
double critPressure() const override
Critical pressure (Pa).
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.
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].
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...