6#ifndef CT_REDLICHKWONGMFTP_H
7#define CT_REDLICHKWONGMFTP_H
30 const std::string&
id=
"");
44 virtual std::string
type()
const {
45 return "RedlichKwong";
50 virtual doublereal
cp_mole()
const;
51 virtual doublereal
cv_mole()
const;
136 virtual bool addSpecies(shared_ptr<Species> spec);
141 AnyMap& speciesNode)
const;
158 virtual std::vector<double>
getCoeff(
const std::string& iName);
192 const std::string& species_j,
double a0,
double a1);
211 virtual doublereal
sresid()
const;
212 virtual doublereal
hresid()
const;
215 virtual doublereal
liquidVolEst(doublereal TKelvin, doublereal& pres)
const;
216 virtual doublereal
densityCalc(doublereal T, doublereal
pressure,
int phase, doublereal rhoguess);
220 virtual doublereal
dpdVCalc(doublereal TKelvin, doublereal molarVol, doublereal& presCalc)
const;
245 void calculateAB(doublereal temp, doublereal& aCalc, doublereal& bCalc)
const;
249 doublereal da_dt()
const;
251 void calcCriticalConditions(doublereal& pc, doublereal& tc, doublereal& vc)
const;
254 int solveCubic(
double T,
double pres,
double a,
double b,
double Vroot[3])
const;
284 enum class CoeffSource { EoS, CritProps, Database };
290 doublereal Vroot_[3];
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.
std::string name() const
Return the name of the phase.
shared_ptr< Species > species(const std::string &name) const
Return the Species object for the named species.
Implementation of a multi-species Redlich-Kwong equation of state.
virtual bool addSpecies(shared_ptr< Species > spec)
virtual void getActivityCoefficients(doublereal *ac) const
Get the array of non-dimensional activity coefficients at the current solution temperature,...
static const doublereal omega_b
Omega constant for b.
int m_formTempParam
Form of the temperature parameterization.
virtual doublereal hresid() const
Calculate the deviation terms for the total enthalpy of the mixture from the ideal gas mixture.
doublereal m_b_current
Value of b in the equation of state.
doublereal dpdT_
The derivative of the pressure wrt the temperature.
virtual doublereal pressure() const
Return the thermodynamic pressure (Pa).
std::map< std::string, std::map< std::string, std::pair< double, double > > > m_binaryParameters
Explicitly-specified binary interaction parameters.
virtual void getPartialMolarIntEnergies(doublereal *ubar) const
Return an array of partial molar internal energies for the species in the mixture.
virtual doublereal cp_mole() const
Molar heat capacity at constant pressure. Units: J/kmol/K.
virtual std::vector< double > getCoeff(const std::string &iName)
Retrieve a and b coefficients by looking up tabulated critical parameters.
virtual void getPartialMolarEnthalpies(doublereal *hbar) const
Returns an array of partial molar enthalpies for the species in the mixture.
virtual void getPartialMolarEntropies(doublereal *sbar) const
Returns an array of partial molar entropies of the species in the solution.
virtual void getPartialMolarCp(double *cpbar) const
Return an array of partial molar heat capacities for the species in the mixture.
virtual void initThermoXML(XML_Node &phaseNode, const std::string &id)
Import and initialize a ThermoPhase object using an XML tree.
void readXMLCrossFluid(XML_Node &pureFluidParam)
Read the cross species RedlichKwong input parameters.
virtual void getPartialMolarVolumes(doublereal *vbar) const
Return an array of partial molar volumes for the species in the mixture.
vector_fp m_pp
Temporary storage - length = m_kk.
RedlichKwongMFTP(const std::string &infile="", const std::string &id="")
Construct a RedlichKwongMFTP object from an input file.
virtual doublereal cv_mole() const
Molar heat capacity at constant volume. Units: J/kmol/K.
void readXMLPureFluid(XML_Node &pureFluidParam)
Read the pure species RedlichKwong input parameters.
virtual doublereal dpdVCalc(doublereal TKelvin, doublereal molarVol, doublereal &presCalc) const
Calculate the pressure and the pressure derivative given the temperature and the molar volume.
std::vector< CoeffSource > m_coeffSource
For each species, specifies the source of the a and b coefficients.
virtual doublereal sresid() const
Calculate the deviation terms for the total entropy of the mixture from the ideal gas mixture.
virtual void updateMixingExpressions()
Update the a and b parameters.
void setBinaryCoeffs(const std::string &species_i, const std::string &species_j, double a0, double a1)
Set values for the interaction parameter between two species.
void pressureDerivatives() const
Calculate dpdV and dpdT at the current conditions.
static const doublereal omega_vc
Omega constant for the critical molar volume.
virtual doublereal liquidVolEst(doublereal TKelvin, doublereal &pres) const
Estimate for the molar volume of the liquid.
virtual void initThermo()
Initialize the ThermoPhase object after all species have been set up.
virtual doublereal densityCalc(doublereal T, doublereal pressure, int phase, doublereal rhoguess)
Calculates the density given the temperature and the pressure and a guess at the density.
doublereal m_a_current
Value of a in the equation of state.
virtual doublereal densSpinodalGas() const
Return the value of the density at the gas spinodal point (on the gas side) for the current temperatu...
virtual std::string type() const
String indicating the thermodynamic model implemented.
void setSpeciesCoeffs(const std::string &species, double a0, double a1, double b)
Set the pure fluid interaction parameters for a species.
virtual void getSpeciesParameters(const std::string &name, AnyMap &speciesNode) const
Get phase-specific parameters of a Species object such that an identical one could be reconstructed a...
void calculateAB(doublereal temp, doublereal &aCalc, doublereal &bCalc) const
Calculate the a and the b parameters given the temperature.
virtual void getChemPotentials_RT(doublereal *mu) const
Get the array of non-dimensional species chemical potentials.
doublereal dpdV_
The derivative of the pressure wrt the volume.
virtual doublereal densSpinodalLiquid() const
Return the value of the density at the liquid spinodal point (on the liquid side) for the current tem...
virtual void getChemPotentials(doublereal *mu) const
Get the species chemical potentials. Units: J/kmol.
static const doublereal omega_a
Omega constant for a -> value of a in terms of critical properties.
vector_fp dpdni_
Vector of derivatives of pressure wrt mole number.
virtual void setParametersFromXML(const XML_Node &thermoNode)
Set equation of state parameter values from XML entries.
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.
virtual doublereal standardConcentration(size_t k=0) const
Returns the standard concentration , which is used to normalize the generalized concentration.
Class XML_Node is a tree-based representation of the contents of an XML file.
Namespace for the Cantera kernel.
std::vector< double > vector_fp
Turn on the use of stl vectors for the basic array type within cantera Vector of doubles.