6#ifndef CT_REDLICHKWONGMFTP_H
7#define CT_REDLICHKWONGMFTP_H
31 string type()
const override {
32 return "Redlich-Kwong";
37 double cp_mole()
const override;
38 double cv_mole()
const override;
108 bool addSpecies(shared_ptr<Species> spec)
override;
143 const string& species_j,
double a0,
double a1);
148 double sresid()
const override;
149 double hresid()
const override;
152 double liquidVolEst(
double TKelvin,
double& pres)
const override;
157 double dpdVCalc(
double TKelvin,
double molarVol,
double& presCalc)
const override;
186 void calculateAB(
double temp,
double& aCalc,
double& bCalc)
const;
190 double da_dt()
const;
192 void calcCriticalConditions(
double& pc,
double& tc,
double& vc)
const override;
195 int solveCubic(
double T,
double pres,
double a,
double b,
double Vroot[3])
const;
217 vector<double> a_vec_Curr_;
218 vector<double> b_vec_Curr_;
225 enum class CoeffSource { EoS, CritProps, Database };
231 double Vroot_[3] = {0.0, 0.0, 0.0};
237 mutable vector<double> m_partialMolarVolumes;
Header file for class Cantera::Array2D.
Header file for a derived class of ThermoPhase that handles non-ideal mixtures based on the fugacity ...
A map of string keys to values whose type can vary at runtime.
A class for 2D arrays stored in column-major (Fortran-compatible) form.
This is a filter class for ThermoPhase that implements some preparatory steps for efficiently handlin...
An error indicating that an unimplemented function has been called.
shared_ptr< Species > species(const string &name) const
Return the Species object for the named species.
string name() const
Return the name of the phase.
Implementation of a multi-species Redlich-Kwong equation of state.
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.
double densSpinodalLiquid() const override
Return the value of the density at the liquid spinodal point (on the liquid side) for the current tem...
void getPartialMolarEnthalpies(double *hbar) const override
Returns an array of partial molar enthalpies for the species in the mixture.
void getChemPotentials(double *mu) const override
Get the species chemical potentials. Units: J/kmol.
double sresid() const override
Calculate the deviation terms for the total entropy of the mixture from the ideal gas mixture.
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...
static const double omega_b
Omega constant for b.
vector< double > m_pp
Temporary storage - length = m_kk.
string type() const override
String indicating the thermodynamic model implemented.
void initThermo() override
Initialize the ThermoPhase object after all species have been set up.
double densSpinodalGas() const override
Return the value of the density at the gas spinodal point (on the gas side) for the current temperatu...
void getPartialMolarVolumes(double *vbar) const override
Return an array of partial molar volumes for the species in the mixture.
double cv_mole() const override
Molar heat capacity at constant volume. Units: 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.
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 getPartialMolarIntEnergies(double *ubar) const override
Return an array of partial molar internal energies for the species in the mixture.
void updateMixingExpressions() override
Update the a and b parameters.
double cp_mole() const override
Molar heat capacity at constant pressure. Units: J/kmol/K.
void getPartialMolarCp(double *cpbar) const override
Return an array of partial molar heat capacities for the species in the mixture.
double m_a_current
Value of a in the equation of state.
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 getActivityCoefficients(double *ac) const override
Get the array of non-dimensional activity coefficients at the current solution temperature,...
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.
void getPartialMolarEntropies(double *sbar) const override
Returns an array of partial molar entropies of the species in the solution.
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.
int solveCubic(double T, double pres, double a, double b, double Vroot[3]) const
Prepare variables and call the function to solve the cubic equation of state.
Namespace for the Cantera kernel.