Cantera  3.1.0
Loading...
Searching...
No Matches
RedlichKwongMFTP.h
Go to the documentation of this file.
1//! @file RedlichKwongMFTP.h
2
3// This file is part of Cantera. See License.txt in the top-level directory or
4// at https://cantera.org/license.txt for license and copyright information.
5
6#ifndef CT_REDLICHKWONGMFTP_H
7#define CT_REDLICHKWONGMFTP_H
8
9#include "MixtureFugacityTP.h"
10#include "cantera/base/Array.h"
11
12namespace Cantera
13{
14/**
15 * Implementation of a multi-species Redlich-Kwong equation of state
16 *
17 * @ingroup thermoprops
18 */
20{
21public:
22 //! Construct a RedlichKwongMFTP object from an input file
23 /*!
24 * @param infile Name of the input file containing the phase definition.
25 * If blank, an empty phase will be created.
26 * @param id name (ID) of the phase in the input file. If empty, the
27 * first phase definition in the input file will be used.
28 */
29 explicit RedlichKwongMFTP(const string& infile="", const string& id="");
30
31 string type() const override {
32 return "Redlich-Kwong";
33 }
34
35 //! @name Molar Thermodynamic properties
36 //! @{
37 double cp_mole() const override;
38 double cv_mole() const override;
39 //! @}
40 //! @name Mechanical Properties
41 //! @{
42
43 //! Return the thermodynamic pressure (Pa).
44 /*!
45 * Since the mass density, temperature, and mass fractions are stored,
46 * this method uses these values to implement the
47 * mechanical equation of state @f$ P(T, \rho, Y_1, \dots, Y_K) @f$.
48 *
49 * @f[
50 * P = \frac{RT}{v-b_{mix}} - \frac{a_{mix}}{T^{0.5} v \left( v + b_{mix} \right) }
51 * @f]
52 */
53 double pressure() const override;
54
55 //! @}
56
57public:
58
59 //! Returns the standard concentration @f$ C^0_k @f$, which is used to
60 //! normalize the generalized concentration.
61 /*!
62 * This is defined as the concentration by which the generalized
63 * concentration is normalized to produce the activity. In many cases, this
64 * quantity will be the same for all species in a phase. Since the activity
65 * for an ideal gas mixture is simply the mole fraction, for an ideal gas
66 * @f$ C^0_k = P/\hat R T @f$.
67 *
68 * @param k Optional parameter indicating the species. The default is to
69 * assume this refers to species 0.
70 * @return
71 * Returns the standard Concentration in units of m3 kmol-1.
72 */
73 double standardConcentration(size_t k=0) const override;
74
75 //! Get the array of non-dimensional activity coefficients at the current
76 //! solution temperature, pressure, and solution concentration.
77 /*!
78 * For all objects with the Mixture Fugacity approximation, we define the
79 * standard state as an ideal gas at the current temperature and pressure of
80 * the solution. The activities are based on this standard state.
81 *
82 * @param ac Output vector of activity coefficients. Length: m_kk.
83 */
84 void getActivityCoefficients(double* ac) const override;
85
86 //! @name Partial Molar Properties of the Solution
87 //! @{
88
89 void getChemPotentials(double* mu) const override;
90 void getPartialMolarEnthalpies(double* hbar) const override;
91 void getPartialMolarEntropies(double* sbar) const override;
92 void getPartialMolarIntEnergies(double* ubar) const override;
93 void getPartialMolarCp(double* cpbar) const override {
94 throw NotImplementedError("RedlichKwongMFTP::getPartialMolarCp");
95 }
96 void getPartialMolarVolumes(double* vbar) const override;
97 //! @}
98
99public:
100 //! @name Initialization Methods - For Internal use
101 //!
102 //! The following methods are used in the process of constructing
103 //! the phase and setting its parameters from a specification in an
104 //! input file. They are not normally used in application programs.
105 //! To see how they are used, see importPhase().
106 //! @{
107
108 bool addSpecies(shared_ptr<Species> spec) override;
109 void initThermo() override;
110 void getSpeciesParameters(const string& name, AnyMap& speciesNode) const override;
111
112 //! Set the pure fluid interaction parameters for a species
113 /*!
114 * The "a" parameter for species *i* in the Redlich-Kwong model is assumed
115 * to be a linear function of temperature:
116 * @f[ a = a_0 + a_1 T @f]
117 *
118 * @param species Name of the species
119 * @param a0 constant term in the expression for the "a" parameter
120 * of the specified species [Pa-m^6/kmol^2]
121 * @param a1 temperature-proportional term in the expression for the
122 * "a" parameter of the specified species [Pa-m^6/kmol^2/K]
123 * @param b "b" parameter in the Redlich-Kwong model [m^3/kmol]
124 */
125 void setSpeciesCoeffs(const string& species, double a0, double a1, double b);
126
127 //! Set values for the interaction parameter between two species
128 /*!
129 * The "a" parameter for interactions between species *i* and *j* is
130 * assumed by default to be computed as:
131 * @f[ a_{ij} = \sqrt(a_{i,0} a_{j,0}) + \sqrt(a_{i,1} a_{j,1}) T @f]
132 *
133 * This function overrides the defaults with the specified parameters:
134 * @f[ a_{ij} = a_{ij,0} + a_{ij,1} T @f]
135 *
136 * @param species_i Name of one species
137 * @param species_j Name of the other species
138 * @param a0 constant term in the "a" expression [Pa-m^6/kmol^2]
139 * @param a1 temperature-proportional term in the "a" expression
140 * [Pa-m^6/kmol^2/K]
141 */
142 void setBinaryCoeffs(const string& species_i,
143 const string& species_j, double a0, double a1);
144 //! @}
145
146protected:
147 // Special functions inherited from MixtureFugacityTP
148 double sresid() const override;
149 double hresid() const override;
150
151public:
152 double liquidVolEst(double TKelvin, double& pres) const override;
153 double densityCalc(double T, double pressure, int phase, double rhoguess) override;
154
155 double densSpinodalLiquid() const override;
156 double densSpinodalGas() const override;
157 double dpdVCalc(double TKelvin, double molarVol, double& presCalc) const override;
158
159 double isothermalCompressibility() const override;
160 double thermalExpansionCoeff() const override;
161 double soundSpeed() const override;
162
163 //! Calculate dpdV and dpdT at the current conditions
164 /*!
165 * These are stored internally.
166 */
167 void pressureDerivatives() const;
168
169 //! Update the a and b parameters
170 /*!
171 * The a and the b parameters depend on the mole fraction and the
172 * temperature. This function updates the internal numbers based on the
173 * state of the object.
174 */
175 void updateMixingExpressions() override;
176
177 //! Calculate the a and the b parameters given the temperature
178 /*!
179 * This function doesn't change the internal state of the object, so it is a
180 * const function. It does use the stored mole fractions in the object.
181 *
182 * @param temp Temperature (TKelvin)
183 * @param aCalc (output) Returns the a value
184 * @param bCalc (output) Returns the b value.
185 */
186 void calculateAB(double temp, double& aCalc, double& bCalc) const;
187
188 // Special functions not inherited from MixtureFugacityTP
189
190 double da_dt() const;
191
192 void calcCriticalConditions(double& pc, double& tc, double& vc) const override;
193
194 //! Prepare variables and call the function to solve the cubic equation of state
195 int solveCubic(double T, double pres, double a, double b, double Vroot[3]) const;
196
197protected:
198 //! Form of the temperature parameterization
199 /*!
200 * - 0 = There is no temperature parameterization of a or b
201 * - 1 = The a_ij parameter is a linear function of the temperature
202 */
204
205 //! Value of b in the equation of state
206 /*!
207 * m_b is a function of the temperature and the mole fraction.
208 */
209 double m_b_current = 0.0;
210
211 //! Value of a in the equation of state
212 /*!
213 * a_b is a function of the temperature and the mole fraction.
214 */
215 double m_a_current = 0.0;
216
217 vector<double> a_vec_Curr_;
218 vector<double> b_vec_Curr_;
219
220 Array2D a_coeff_vec;
221
222 //! Explicitly-specified binary interaction parameters
223 map<string, map<string, pair<double, double>>> m_binaryParameters;
224
225 enum class CoeffSource { EoS, CritProps, Database };
226 //! For each species, specifies the source of the a and b coefficients
227 vector<CoeffSource> m_coeffSource;
228
229 int NSolns_ = 0;
230
231 double Vroot_[3] = {0.0, 0.0, 0.0};
232
233 //! Temporary storage - length = m_kk.
234 mutable vector<double> m_pp;
235
236 // Partial molar volumes of the species
237 mutable vector<double> m_partialMolarVolumes;
238
239 //! The derivative of the pressure wrt the volume
240 /*!
241 * Calculated at the current conditions. temperature and mole number kept
242 * constant
243 */
244 mutable double dpdV_ = 0.0;
245
246 //! The derivative of the pressure wrt the temperature
247 /*!
248 * Calculated at the current conditions. Total volume and mole number kept
249 * constant
250 */
251 mutable double dpdT_ = 0.0;
252
253 //! Vector of derivatives of pressure wrt mole number
254 /*!
255 * Calculated at the current conditions. Total volume, temperature and
256 * other mole number kept constant
257 */
258 mutable vector<double> dpdni_;
259
260private:
261 //! Omega constant for a -> value of a in terms of critical properties
262 /*!
263 * this was calculated from a small nonlinear solve
264 */
265 static const double omega_a;
266
267 //! Omega constant for b
268 static const double omega_b;
269
270 //! Omega constant for the critical molar volume
271 static const double omega_vc;
272};
273}
274
275#endif
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.
Definition AnyMap.h:431
A class for 2D arrays stored in column-major (Fortran-compatible) form.
Definition Array.h:32
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.
Definition Phase.cpp:873
string name() const
Return the name of the phase.
Definition Phase.cpp:20
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.
Definition AnyMap.cpp:595