Cantera  2.5.1
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 
12 namespace Cantera
13 {
14 /**
15  * Implementation of a multi-species Redlich-Kwong equation of state
16  *
17  * @ingroup thermoprops
18  */
20 {
21 public:
22  //! @name Constructors and Duplicators
23  //! @{
24 
25  //! Base constructor.
27 
28  //! Construct a RedlichKwongMFTP object from an input file
29  /*!
30  * @param inputFile Name of the input file containing the phase definition
31  * @param id name (ID) of the phase in the input file. If empty, the
32  * first phase definition in the input file will be used.
33  */
34  RedlichKwongMFTP(const std::string& infile, const std::string& id="");
35 
36  //! Construct and initialize a RedlichKwongMFTP object directly from an
37  //! XML database
38  /*!
39  * @param phaseRef XML phase node containing the description of the phase
40  * @param id id attribute containing the name of the phase. (default
41  * is the empty string)
42  *
43  * @deprecated The XML input format is deprecated and will be removed in
44  * Cantera 3.0.
45  */
46  RedlichKwongMFTP(XML_Node& phaseRef, const std::string& id = "");
47 
48  virtual std::string type() const {
49  return "RedlichKwong";
50  }
51 
52  //! @name Molar Thermodynamic properties
53  //! @{
54 
55  virtual doublereal enthalpy_mole() const;
56  virtual doublereal entropy_mole() const;
57  virtual doublereal cp_mole() const;
58  virtual doublereal cv_mole() const;
59 
60  //! @}
61  //! @name Mechanical Properties
62  //! @{
63 
64  //! Return the thermodynamic pressure (Pa).
65  /*!
66  * Since the mass density, temperature, and mass fractions are stored,
67  * this method uses these values to implement the
68  * mechanical equation of state \f$ P(T, \rho, Y_1, \dots, Y_K) \f$.
69  *
70  * \f[
71  * P = \frac{RT}{v-b_{mix}} - \frac{a_{mix}}{T^{0.5} v \left( v + b_{mix} \right) }
72  * \f]
73  */
74  virtual doublereal pressure() const;
75 
76  // @}
77 
78 protected:
79  /**
80  * Calculate the density of the mixture using the partial molar volumes and
81  * mole fractions as input
82  *
83  * The formula for this is
84  *
85  * \f[
86  * \rho = \frac{\sum_k{X_k W_k}}{\sum_k{X_k V_k}}
87  * \f]
88  *
89  * where \f$X_k\f$ are the mole fractions, \f$W_k\f$ are the molecular
90  * weights, and \f$V_k\f$ are the pure species molar volumes.
91  *
92  * Note, the basis behind this formula is that in an ideal solution the
93  * partial molar volumes are equal to the species standard state molar
94  * volumes. The species molar volumes may be functions of temperature and
95  * pressure.
96  */
97  virtual void calcDensity();
98 
99  virtual void setTemperature(const doublereal temp);
100  virtual void compositionChanged();
101 
102 public:
103  virtual void getActivityConcentrations(doublereal* c) const;
104 
105  //! Returns the standard concentration \f$ C^0_k \f$, which is used to
106  //! normalize the generalized concentration.
107  /*!
108  * This is defined as the concentration by which the generalized
109  * concentration is normalized to produce the activity. In many cases, this
110  * quantity will be the same for all species in a phase. Since the activity
111  * for an ideal gas mixture is simply the mole fraction, for an ideal gas
112  * \f$ C^0_k = P/\hat R T \f$.
113  *
114  * @param k Optional parameter indicating the species. The default is to
115  * assume this refers to species 0.
116  * @return
117  * Returns the standard Concentration in units of m3 kmol-1.
118  */
119  virtual doublereal standardConcentration(size_t k=0) const;
120 
121  //! Get the array of non-dimensional activity coefficients at the current
122  //! solution temperature, pressure, and solution concentration.
123  /*!
124  * For all objects with the Mixture Fugacity approximation, we define the
125  * standard state as an ideal gas at the current temperature and pressure of
126  * the solution. The activities are based on this standard state.
127  *
128  * @param ac Output vector of activity coefficients. Length: m_kk.
129  */
130  virtual void getActivityCoefficients(doublereal* ac) const;
131 
132  /// @name Partial Molar Properties of the Solution
133  //@{
134 
135  //! Get the array of non-dimensional species chemical potentials.
136  //! These are partial molar Gibbs free energies.
137  /*!
138  * \f$ \mu_k / \hat R T \f$.
139  * Units: unitless
140  *
141  * We close the loop on this function, here, calling getChemPotentials() and
142  * then dividing by RT. No need for child classes to handle.
143  *
144  * @param mu Output vector of non-dimensional species chemical potentials
145  * Length: m_kk.
146  */
147  virtual void getChemPotentials_RT(doublereal* mu) const;
148 
149  virtual void getChemPotentials(doublereal* mu) const;
150  virtual void getPartialMolarEnthalpies(doublereal* hbar) const;
151  virtual void getPartialMolarEntropies(doublereal* sbar) const;
152  virtual void getPartialMolarIntEnergies(doublereal* ubar) const;
153  virtual void getPartialMolarCp(doublereal* cpbar) const;
154  virtual void getPartialMolarVolumes(doublereal* vbar) const;
155 
156  //@}
157  /// @name Critical State Properties.
158  //@{
159 
160  virtual doublereal critTemperature() const;
161  virtual doublereal critPressure() const;
162  virtual doublereal critVolume() const;
163  virtual doublereal critCompressibility() const;
164  virtual doublereal critDensity() const;
165 
166 public:
167  //@}
168  //! @name Initialization Methods - For Internal use
169  /*!
170  * The following methods are used in the process of constructing
171  * the phase and setting its parameters from a specification in an
172  * input file. They are not normally used in application programs.
173  * To see how they are used, see importPhase().
174  */
175  //@{
176 
177  virtual bool addSpecies(shared_ptr<Species> spec);
178  virtual void setParametersFromXML(const XML_Node& thermoNode);
179  virtual void initThermoXML(XML_Node& phaseNode, const std::string& id);
180  virtual void initThermo();
181 
182  //! Retrieve a and b coefficients by looking up tabulated critical parameters
183  /*!
184  * If pureFluidParameters are not provided for any species in the phase,
185  * consult the critical properties tabulated in /thermo/critProperties.xml.
186  * If the species is found there, calculate pure fluid parameters a_k and b_k as:
187  * \f[ a_k = 0.4278*R**2*T_c^2.5/P_c \f]
188  *
189  * and:
190  * \f[ b_k = 0.08664*R*T_c/P_c \f]
191  *
192  * @param iName Name of the species
193  */
194  virtual std::vector<double> getCoeff(const std::string& iName);
195 
196  //! Set the pure fluid interaction parameters for a species
197  /*!
198  * The "a" parameter for species *i* in the Redlich-Kwong model is assumed
199  * to be a linear function of temperature:
200  * \f[ a = a_0 + a_1 T \f]
201  *
202  * @param species Name of the species
203  * @param a0 constant term in the expression for the "a" parameter
204  * of the specified species [Pa-m^6/kmol^2]
205  * @param a1 temperature-proportional term in the expression for the
206  * "a" parameter of the specified species [Pa-m^6/kmol^2/K]
207  * @param b "b" parameter in the Redlich-Kwong model [m^3/kmol]
208  */
209  void setSpeciesCoeffs(const std::string& species, double a0, double a1,
210  double b);
211 
212  //! Set values for the interaction parameter between two species
213  /*!
214  * The "a" parameter for interactions between species *i* and *j* is
215  * assumed by default to be computed as:
216  * \f[ a_{ij} = \sqrt(a_{i,0} a_{j,0}) + \sqrt(a_{i,1} a_{j,1}) T \f]
217  *
218  * This function overrides the defaults with the specified parameters:
219  * \f[ a_{ij} = a_{ij,0} + a_{ij,1} T \f]
220  *
221  * @param species_i Name of one species
222  * @param species_j Name of the other species
223  * @param a0 constant term in the "a" expression [Pa-m^6/kmol^2]
224  * @param a1 temperature-proportional term in the "a" expression
225  * [Pa-m^6/kmol^2/K]
226  */
227  void setBinaryCoeffs(const std::string& species_i,
228  const std::string& species_j, double a0, double a1);
229 
230 private:
231  //! Read the pure species RedlichKwong input parameters
232  /*!
233  * @param pureFluidParam XML_Node for the pure fluid parameters
234  */
235  void readXMLPureFluid(XML_Node& pureFluidParam);
236 
237  //! Read the cross species RedlichKwong input parameters
238  /*!
239  * @param pureFluidParam XML_Node for the cross fluid parameters
240  */
241  void readXMLCrossFluid(XML_Node& pureFluidParam);
242 
243  // @}
244 
245 protected:
246  // Special functions inherited from MixtureFugacityTP
247  virtual doublereal sresid() const;
248  virtual doublereal hresid() const;
249 
250 public:
251  virtual doublereal liquidVolEst(doublereal TKelvin, doublereal& pres) const;
252  virtual doublereal densityCalc(doublereal TKelvin, doublereal pressure, int phase, doublereal rhoguess);
253 
254  virtual doublereal densSpinodalLiquid() const;
255  virtual doublereal densSpinodalGas() const;
256  virtual doublereal pressureCalc(doublereal TKelvin, doublereal molarVol) const;
257  virtual doublereal dpdVCalc(doublereal TKelvin, doublereal molarVol, doublereal& presCalc) const;
258 
259  //! Calculate dpdV and dpdT at the current conditions
260  /*!
261  * These are stored internally.
262  */
263  void pressureDerivatives() const;
264 
265  virtual void updateMixingExpressions();
266 
267  //! Update the a and b parameters
268  /*!
269  * The a and the b parameters depend on the mole fraction and the
270  * temperature. This function updates the internal numbers based on the
271  * state of the object.
272  */
273  void updateAB();
274 
275  //! Calculate the a and the b parameters given the temperature
276  /*!
277  * This function doesn't change the internal state of the object, so it is a
278  * const function. It does use the stored mole fractions in the object.
279  *
280  * @param temp Temperature (TKelvin)
281  * @param aCalc (output) Returns the a value
282  * @param bCalc (output) Returns the b value.
283  */
284  void calculateAB(doublereal temp, doublereal& aCalc, doublereal& bCalc) const;
285 
286  // Special functions not inherited from MixtureFugacityTP
287 
288  doublereal da_dt() const;
289 
290  void calcCriticalConditions(doublereal a, doublereal b, doublereal a0_coeff, doublereal aT_coeff,
291  doublereal& pc, doublereal& tc, doublereal& vc) const;
292 
293  //! Solve the cubic equation of state
294  /*!
295  * The R-K equation of state may be solved via the following formula:
296  *
297  * V**3 - V**2(RT/P) - V(RTb/P - a/(P T**.5) + b*b) - (a b / (P T**.5)) = 0
298  *
299  * Returns the number of solutions found. If it only finds the liquid
300  * branch solution, it will return a -1 or a -2 instead of 1 or 2. If it
301  * returns 0, then there is an error.
302  */
303  int NicholsSolve(double TKelvin, double pres, doublereal a, doublereal b,
304  doublereal Vroot[3]) const;
305 
306 protected:
307  //! Form of the temperature parameterization
308  /*!
309  * - 0 = There is no temperature parameterization of a or b
310  * - 1 = The a_ij parameter is a linear function of the temperature
311  */
313 
314  //! Value of b in the equation of state
315  /*!
316  * m_b is a function of the temperature and the mole fraction.
317  */
318  doublereal m_b_current;
319 
320  //! Value of a in the equation of state
321  /*!
322  * a_b is a function of the temperature and the mole fraction.
323  */
324  doublereal m_a_current;
325 
326  vector_fp a_vec_Curr_;
327  vector_fp b_vec_Curr_;
328 
329  Array2D a_coeff_vec;
330 
331  int NSolns_;
332 
333  doublereal Vroot_[3];
334 
335  //! Temporary storage - length = m_kk.
336  mutable vector_fp m_pp;
337 
338  //! Temporary storage - length = m_kk.
339  mutable vector_fp m_tmpV;
340 
341  // Partial molar volumes of the species
342  mutable vector_fp m_partialMolarVolumes;
343 
344  //! The derivative of the pressure wrt the volume
345  /*!
346  * Calculated at the current conditions. temperature and mole number kept
347  * constant
348  */
349  mutable doublereal dpdV_;
350 
351  //! The derivative of the pressure wrt the temperature
352  /*!
353  * Calculated at the current conditions. Total volume and mole number kept
354  * constant
355  */
356  mutable doublereal dpdT_;
357 
358  //! Vector of derivatives of pressure wrt mole number
359  /*!
360  * Calculated at the current conditions. Total volume, temperature and
361  * other mole number kept constant
362  */
363  mutable vector_fp dpdni_;
364 
365 public:
366  //! Omega constant for a -> value of a in terms of critical properties
367  /*!
368  * this was calculated from a small nonlinear solve
369  */
370  static const doublereal omega_a;
371 
372  //! Omega constant for b
373  static const doublereal omega_b;
374 
375  //! Omega constant for the critical molar volume
376  static const doublereal omega_vc;
377 };
378 }
379 
380 #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 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...
shared_ptr< Species > species(const std::string &name) const
Return the Species object for the named species.
Definition: Phase.cpp:980
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,...
void updateAB()
Update the a and b parameters.
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).
virtual void getPartialMolarIntEnergies(doublereal *ubar) const
Return an array of partial molar internal energies for the species in the mixture.
virtual doublereal critVolume() const
Critical volume (m3/kmol).
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 initThermoXML(XML_Node &phaseNode, const std::string &id)
Import and initialize a ThermoPhase object using an XML tree.
virtual doublereal critPressure() const
Critical pressure (Pa).
virtual doublereal enthalpy_mole() const
Molar enthalpy. Units: J/kmol.
virtual doublereal densityCalc(doublereal TKelvin, doublereal pressure, int phase, doublereal rhoguess)
Calculates the density given the temperature and the pressure and a guess at the density.
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.
virtual doublereal cv_mole() const
Molar heat capacity at constant volume. Units: J/kmol/K.
virtual void getActivityConcentrations(doublereal *c) const
This method returns an array of generalized concentrations.
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.
virtual doublereal sresid() const
Calculate the deviation terms for the total entropy of the mixture from the ideal gas mixture.
virtual void getPartialMolarCp(doublereal *cpbar) const
Return an array of partial molar heat capacities for the species in the mixture.
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.
virtual doublereal entropy_mole() const
Molar entropy. Units: J/kmol/K.
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 critTemperature() const
Critical temperature (K).
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.
RedlichKwongMFTP()
Base constructor.
virtual void setTemperature(const doublereal temp)
Set the temperature of the phase.
int NicholsSolve(double TKelvin, double pres, doublereal a, doublereal b, doublereal Vroot[3]) const
Solve the cubic equation of state.
virtual doublereal pressureCalc(doublereal TKelvin, doublereal molarVol) const
Calculate the pressure given the temperature and the molar volume.
virtual void compositionChanged()
Apply changes to the state which are needed after the composition changes.
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 critCompressibility() const
Critical compressibility (unitless).
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.
virtual doublereal critDensity() const
Critical density (kg/m3).
vector_fp m_tmpV
Temporary storage - length = m_kk.
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.
virtual void calcDensity()
Calculate the density of the mixture using the partial molar volumes and mole fractions as input.
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.
Definition: xml.h:104
std::vector< double > vector_fp
Turn on the use of stl vectors for the basic array type within cantera Vector of doubles.
Definition: ct_defs.h:180
Namespace for the Cantera kernel.
Definition: AnyMap.cpp:264