Cantera  2.2.1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
PureFluidPhase.h
Go to the documentation of this file.
1 /**
2  * @file PureFluidPhase.h
3  *
4  * Header for a ThermoPhase class for a pure fluid phase consisting of
5  * gas, liquid, mixed-gas-liquid and supercrit fluid (see \ref thermoprops
6  * and class \link Cantera::PureFluidPhase PureFluidPhase\endlink).
7  *
8  * It inherits from ThermoPhase, but is built on top of the tpx package.
9  */
10 
11 // Copyright 2003 California Institute of Technology
12 
13 #ifndef CT_EOS_TPX_H
14 #define CT_EOS_TPX_H
15 
16 #include "ThermoPhase.h"
17 #include "mix_defs.h"
18 #include "cantera/tpx/Sub.h"
19 
20 namespace Cantera
21 {
22 //! This phase object consists of a single component that can be a
23 //! gas, a liquid, a mixed gas-liquid fluid, or a fluid beyond its
24 //! critical point
25 /*!
26  * The object inherits from ThermoPhase. However, it's built on top
27  * of the tpx package.
28  *
29  * @ingroup thermoprops
30  */
32 {
33 public:
34 
35  //! Empty Base Constructor
37 
38  //! Copy Constructor
39  /*!
40  * @param right Object to be copied
41  */
42  PureFluidPhase(const PureFluidPhase& right);
43 
44  //! Assignment operator
45  /*!
46  * @param right Object to be copied
47  */
49 
50  //! Destructor
51  virtual ~PureFluidPhase();
52 
53  //! Duplication function
54  /*!
55  * This virtual function is used to create a duplicate of the
56  * current phase. It's used to duplicate the phase when given
57  * a ThermoPhase pointer to the phase.
58  *
59  * @return It returns a ThermoPhase pointer.
60  */
62 
63  //! Equation of state type
64  virtual int eosType() const {
65  return cPureFluid;
66  }
67 
68  /// Molar enthalpy. Units: J/kmol.
69  virtual doublereal enthalpy_mole() const;
70 
71  /// Molar internal energy. Units: J/kmol.
72  virtual doublereal intEnergy_mole() const;
73 
74  /// Molar entropy. Units: J/kmol/K.
75  virtual doublereal entropy_mole() const;
76 
77  /// Molar Gibbs function. Units: J/kmol.
78  virtual doublereal gibbs_mole() const;
79 
80  /// Molar heat capacity at constant pressure. Units: J/kmol/K.
81  virtual doublereal cp_mole() const;
82 
83  /// Molar heat capacity at constant volume. Units: J/kmol/K.
84  virtual doublereal cv_mole() const;
85 
86  //! Return the thermodynamic pressure (Pa).
87  /*!
88  * This method calculates the current pressure consistent with the
89  * independent variables, T, rho.
90  */
91  virtual doublereal pressure() const;
92 
93  //! sets the thermodynamic pressure (Pa).
94  /*!
95  * This method calculates the density that is consistent with the
96  * desired pressure, given the temperature.
97  *
98  * @param p Pressure (Pa)
99  */
100  virtual void setPressure(doublereal p);
101 
102  //! Get the species chemical potentials. Units: J/kmol.
103  /*!
104  * This function returns a vector of chemical potentials of the
105  * species in solution at the current temperature, pressure
106  * and mole fraction of the solution.
107  *
108  * @param mu Output vector of species chemical
109  * potentials. Length: m_kk. Units: J/kmol
110  */
111  virtual void getChemPotentials(doublereal* mu) const {
112  mu[0] = gibbs_mole();
113  }
114 
115  //! Returns an array of partial molar enthalpies for the species
116  //! in the mixture. Units (J/kmol)
117  /*!
118  * @param hbar Output vector of species partial molar enthalpies.
119  * Length: m_kk. units are J/kmol.
120  */
121  virtual void getPartialMolarEnthalpies(doublereal* hbar) const;
122 
123  //! Returns an array of partial molar entropies of the species in the
124  //! solution. Units: J/kmol/K.
125  /*!
126  * @param sbar Output vector of species partial molar entropies.
127  * Length = m_kk. units are J/kmol/K.
128  */
129  virtual void getPartialMolarEntropies(doublereal* sbar) const;
130 
131  //! Return an array of partial molar internal energies for the
132  //! species in the mixture. Units: J/kmol.
133  /*!
134  * @param ubar Output vector of species partial molar internal energies.
135  * Length = m_kk. units are J/kmol.
136  */
137  virtual void getPartialMolarIntEnergies(doublereal* ubar) const;
138 
139  //! Return an array of partial molar heat capacities for the
140  //! species in the mixture. Units: J/kmol/K
141  /*!
142  * @param cpbar Output vector of species partial molar heat
143  * capacities at constant pressure.
144  * Length = m_kk. units are J/kmol/K.
145  */
146  virtual void getPartialMolarCp(doublereal* cpbar) const;
147 
148  //! Return an array of partial molar volumes for the
149  //! species in the mixture. Units: m^3/kmol.
150  /*!
151  * @param vbar Output vector of species partial molar volumes.
152  * Length = m_kk. units are m^3/kmol.
153  */
154  virtual void getPartialMolarVolumes(doublereal* vbar) const;
155 
156  //! This method returns an array of generalized concentrations
157  /*!
158  * \f$ C^a_k\f$ are defined such that \f$ a_k = C^a_k /
159  * C^0_k, \f$ where \f$ C^0_k \f$ is a standard concentration
160  * defined below and \f$ a_k \f$ are activities used in the
161  * thermodynamic functions. These activity (or generalized)
162  * concentrations are used
163  * by kinetics manager classes to compute the forward and
164  * reverse rates of elementary reactions. Note that they may
165  * or may not have units of concentration --- they might be
166  * partial pressures, mole fractions, or surface coverages,
167  * for example.
168  *
169  * @param c Output array of generalized concentrations. The
170  * units depend upon the implementation of the
171  * reaction rate expressions within the phase.
172  */
173  virtual void getActivityConcentrations(doublereal* c) const;
174 
175  //! Return the standard concentration for the kth species
176  /*!
177  * The standard concentration \f$ C^0_k \f$ used to normalize
178  * the activity (i.e., generalized) concentration. In many cases, this quantity
179  * will be the same for all species in a phase - for example,
180  * for an ideal gas \f$ C^0_k = P/\hat R T \f$. For this
181  * reason, this method returns a single value, instead of an
182  * array. However, for phases in which the standard
183  * concentration is species-specific (e.g. surface species of
184  * different sizes), this method may be called with an
185  * optional parameter indicating the species.
186  *
187  * @param k Optional parameter indicating the species. The default
188  * is to assume this refers to species 0.
189  * @return
190  * Returns the standard concentration. The units are by definition
191  * dependent on the ThermoPhase and kinetics manager representation.
192  */
193  virtual doublereal standardConcentration(size_t k=0) const;
194 
195  //! Get the array of non-dimensional activities at
196  //! the current solution temperature, pressure, and solution concentration.
197  /*!
198  * Note, for molality based formulations, this returns the
199  * molality based activities.
200  *
201  * We resolve this function at this level by calling
202  * on the activityConcentration function. However,
203  * derived classes may want to override this default
204  * implementation.
205  *
206  * @param a Output vector of activities. Length: m_kk.
207  */
208  virtual void getActivities(doublereal* a) const;
209 
210  //! Returns the isothermal compressibility. Units: 1/Pa.
211  /*!
212  * The isothermal compressibility is defined as
213  * \f[
214  * \kappa_T = -\frac{1}{v}\left(\frac{\partial v}{\partial P}\right)_T
215  * \f]
216  */
217  virtual doublereal isothermalCompressibility() const;
218 
219  //! Return the volumetric thermal expansion coefficient. Units: 1/K.
220  /*!
221  * The thermal expansion coefficient is defined as
222  * \f[
223  * \beta = \frac{1}{v}\left(\frac{\partial v}{\partial T}\right)_P
224  * \f]
225  */
226  virtual doublereal thermalExpansionCoeff() const;
227 
228  //! Returns a reference to the substance object
230 
231  //@}
232  /// @name Properties of the Standard State of the Species in the Solution
233  /*!
234  * The standard state of the pure fluid is defined as the real properties
235  * of the pure fluid at the most stable state of the fluid at the current
236  * temperature and pressure of the solution. With this definition, the
237  * activity of the fluid is always then defined to be equal to one.
238  */
239  //@{
240 
241  //! Get the array of chemical potentials at unit activity for the species
242  //! at their standard states at the current <I>T</I> and <I>P</I> of the solution.
243  /*!
244  * These are the standard state chemical potentials \f$ \mu^0_k(T,P)
245  * \f$. The values are evaluated at the current
246  * temperature and pressure of the solution
247  *
248  * @param mu Output vector of chemical potentials.
249  * Length: m_kk.
250  */
251  virtual void getStandardChemPotentials(doublereal* mu) const;
252 
253  //! Get the nondimensional Enthalpy functions for the species
254  //! at their standard states at the current <I>T</I> and <I>P</I> of the solution.
255  /*!
256  * @param hrt Output vector of nondimensional standard state enthalpies.
257  * Length: m_kk.
258  */
259  virtual void getEnthalpy_RT(doublereal* hrt) const;
260 
261  //! Get the array of nondimensional Entropy functions for the
262  //! standard state species at the current <I>T</I> and <I>P</I> of the solution.
263  /*!
264  * @param sr Output vector of nondimensional standard state entropies.
265  * Length: m_kk.
266  */
267  virtual void getEntropy_R(doublereal* sr) const;
268 
269  //! Get the nondimensional Gibbs functions for the species
270  //! in their standard states at the current <I>T</I> and <I>P</I> of the solution.
271  /*!
272  * @param grt Output vector of nondimensional standard state Gibbs free energies
273  * Length: m_kk.
274  */
275  virtual void getGibbs_RT(doublereal* grt) const;
276 
277  //@}
278 
279  /// @name Thermodynamic Values for the Species Reference States
280  /*!
281  * The species reference state for pure fluids is defined as an ideal gas at the
282  * reference pressure and current temperature of the fluid.
283  */
284  //@{
285 
286  //! Returns the vector of nondimensional enthalpies of the reference state at the current temperature
287  //! of the solution and the reference pressure for the species.
288  /*!
289  * @param hrt Output vector containing the nondimensional reference state enthalpies
290  * Length: m_kk.
291  */
292  virtual void getEnthalpy_RT_ref(doublereal* hrt) const;
293 
294  //! Returns the vector of nondimensional Gibbs Free Energies of the reference state at the current temperature
295  //! of the solution and the reference pressure for the species.
296  /*!
297  * @param grt Output vector containing the nondimensional reference state
298  * Gibbs Free energies. Length: m_kk.
299  */
300  virtual void getGibbs_RT_ref(doublereal* grt) const;
301 
302  //! Returns the vector of the Gibbs function of the reference state at the current temperature
303  //! of the solution and the reference pressure for the species.
304  /*!
305  * units = J/kmol
306  *
307  * @param g Output vector containing the reference state
308  * Gibbs Free energies. Length: m_kk. Units: J/kmol.
309  */
310  virtual void getGibbs_ref(doublereal* g) const;
311 
312  //! Returns the vector of nondimensional entropies of the reference state at the current temperature
313  //! of the solution and the reference pressure for each species.
314  /*!
315  * @param er Output vector containing the nondimensional reference state
316  * entropies. Length: m_kk.
317  */
318  virtual void getEntropy_R_ref(doublereal* er) const;
319 
320  /**
321  * @name Setting the State
322  *
323  * These methods set all or part of the thermodynamic state.
324  * @{
325  */
326 
327  //! Set the internally stored specific enthalpy (J/kg) and pressure (Pa) of the phase.
328  /*!
329  * @param h Specific enthalpy (J/kg)
330  * @param p Pressure (Pa)
331  * @param tol Optional parameter setting the tolerance of the
332  * calculation.
333  */
334  virtual void setState_HP(doublereal h, doublereal p,
335  doublereal tol = 1.e-8);
336 
337  //! Set the specific internal energy (J/kg) and specific volume (m^3/kg).
338  /*!
339  * This function fixes the internal state of the phase so that
340  * the specific internal energy and specific volume have the value of the input parameters.
341  *
342  * @param u specific internal energy (J/kg)
343  * @param v specific volume (m^3/kg).
344  * @param tol Optional parameter setting the tolerance of the
345  * calculation.
346  */
347  virtual void setState_UV(doublereal u, doublereal v,
348  doublereal tol = 1.e-8);
349 
350  //! Set the specific entropy (J/kg/K) and specific volume (m^3/kg).
351  /*!
352  * This function fixes the internal state of the phase so that
353  * the specific entropy and specific volume have the value of the input parameters.
354  *
355  * @param s specific entropy (J/kg/K)
356  * @param v specific volume (m^3/kg).
357  * @param tol Optional parameter setting the tolerance of the
358  * calculation.
359  */
360  virtual void setState_SV(doublereal s, doublereal v,
361  doublereal tol = 1.e-8);
362 
363  //! Set the specific entropy (J/kg/K) and pressure (Pa).
364  /*!
365  * This function fixes the internal state of the phase so that
366  * the specific entropy and the pressure have the value of the input parameters.
367  *
368  * @param s specific entropy (J/kg/K)
369  * @param p specific pressure (Pa).
370  * @param tol Optional parameter setting the tolerance of the
371  * calculation.
372  */
373  virtual void setState_SP(doublereal s, doublereal p,
374  doublereal tol = 1.e-8);
375  //@}
376 
377  //! @name Critical State Properties
378  //@{
379 
380  //! critical temperature
381  virtual doublereal critTemperature() const;
382 
383  //! critical pressure
384  virtual doublereal critPressure() const;
385 
386  //! critical density
387  virtual doublereal critDensity() const;
388 
389  //@}
390 
391  //! @name Saturation properties.
392  //@{
393 
394  //! saturation temperature
395  /*!
396  * @param p Pressure (Pa)
397  */
398  virtual doublereal satTemperature(doublereal p) const;
399 
400  //! Return the saturation pressure given the temperature
401  /*!
402  * @param t Temperature (Kelvin)
403  */
404  virtual doublereal satPressure(doublereal t);
405 
406  //! Return the fraction of vapor at the current conditions
407  virtual doublereal vaporFraction() const;
408 
409  //! Set the state to a saturated system at a particular temperature
410  /*!
411  * @param t Temperature (kelvin)
412  * @param x Fraction of vapor
413  */
414  virtual void setState_Tsat(doublereal t, doublereal x);
415 
416  //! Set the state to a saturated system at a particular pressure
417  /*!
418  * @param p Pressure (Pa)
419  * @param x Fraction of vapor
420  */
421  virtual void setState_Psat(doublereal p, doublereal x);
422  //@}
423 
424  //! Initialize the ThermoPhase object after all species have been set up
425  /*!
426  * @internal Initialize.
427  *
428  * This method is provided to allow
429  * subclasses to perform any initialization required after all
430  * species have been added. For example, it might be used to
431  * resize internal work arrays that must have an entry for
432  * each species. The base class implementation does nothing,
433  * and subclasses that do not require initialization do not
434  * need to overload this method. When importing a CTML phase
435  * description, this method is called from ThermoPhase::initThermoXML(),
436  * which is called from importPhase(),
437  * just prior to returning from function importPhase().
438  */
439  virtual void initThermo();
440 
441  //! Set equation of state parameter values from XML entries.
442  /*!
443  * This method is called by function importPhase() when processing a phase
444  * definition in an input file. It should be overloaded in subclasses to set
445  * any parameters that are specific to that particular phase
446  * model. Note, this method is called before the phase is
447  * initialized with elements and/or species.
448  *
449  * @param eosdata An XML_Node object corresponding to
450  * the "thermo" entry for this phase in the input file.
451  */
452  virtual void setParametersFromXML(const XML_Node& eosdata);
453 
454  //! returns a summary of the state of the phase as a string
455  /*!
456  * @param show_thermo If true, extra information is printed out
457  * about the thermodynamic state of the system.
458  * @param threshold Unused in this subclass
459  */
460  virtual std::string report(bool show_thermo=true,
461  doublereal threshold=1e-14) const;
462 
463 protected:
464 
465  //! Main call to the tpx level to set the state of the system
466  /*!
467  * @param n Integer indicating which 2 thermo components are held constant
468  * @param x Value of the first component
469  * @param y Value of the second component
470  */
471  void Set(tpx::PropertyPair::type n, double x, double y) const;
472 
473  //! Sets the state using a TPX::TV call
474  void setTPXState() const;
475 
476 private:
477  //! Pointer to the underlying tpx object Substance that does the work
479 
480  //! Int indicating the type of the fluid
481  /*!
482  * The tpx package uses an int to indicate what fluid is being sought.
483  */
485 
486  //! Molecular weight of the substance (kg kmol-1)
487  doublereal m_mw;
488 
489  //! flag to turn on some printing.
490  bool m_verbose;
491 };
492 
493 }
494 
495 #endif
PureFluidPhase()
Empty Base Constructor.
virtual ~PureFluidPhase()
Destructor.
virtual void setPressure(doublereal p)
sets the thermodynamic pressure (Pa).
virtual doublereal gibbs_mole() const
Molar Gibbs function. Units: J/kmol.
bool m_verbose
flag to turn on some printing.
This phase object consists of a single component that can be a gas, a liquid, a mixed gas-liquid flui...
virtual void getEnthalpy_RT_ref(doublereal *hrt) const
Returns the vector of nondimensional enthalpies of the reference state at the current temperature of ...
virtual void getActivities(doublereal *a) const
Get the array of non-dimensional activities at the current solution temperature, pressure, and solution concentration.
virtual doublereal intEnergy_mole() const
Molar internal energy. Units: J/kmol.
virtual void setState_UV(doublereal u, doublereal v, doublereal tol=1.e-8)
Set the specific internal energy (J/kg) and specific volume (m^3/kg).
tpx::Substance & TPX_Substance()
Returns a reference to the substance object.
virtual void getEntropy_R(doublereal *sr) const
Get the array of nondimensional Entropy functions for the standard state species at the current T and...
Class XML_Node is a tree-based representation of the contents of an XML file.
Definition: xml.h:100
int m_subflag
Int indicating the type of the fluid.
virtual void getPartialMolarVolumes(doublereal *vbar) const
Return an array of partial molar volumes for the species in the mixture.
doublereal m_mw
Molecular weight of the substance (kg kmol-1)
virtual void getPartialMolarCp(doublereal *cpbar) const
Return an array of partial molar heat capacities for the species in the mixture.
virtual void setState_Tsat(doublereal t, doublereal x)
Set the state to a saturated system at a particular temperature.
virtual doublereal isothermalCompressibility() const
Returns the isothermal compressibility. Units: 1/Pa.
virtual void getPartialMolarEntropies(doublereal *sbar) const
Returns an array of partial molar entropies of the species in the solution.
virtual void setState_SV(doublereal s, doublereal v, doublereal tol=1.e-8)
Set the specific entropy (J/kg/K) and specific volume (m^3/kg).
virtual doublereal critPressure() const
critical pressure
virtual void getGibbs_RT(doublereal *grt) const
Get the nondimensional Gibbs functions for the species in their standard states at the current T and ...
Base class for a phase with thermodynamic properties.
Definition: ThermoPhase.h:97
tpx::Substance * m_sub
Pointer to the underlying tpx object Substance that does the work.
virtual void setState_Psat(doublereal p, doublereal x)
Set the state to a saturated system at a particular pressure.
virtual doublereal satTemperature(doublereal p) const
saturation temperature
virtual doublereal critTemperature() const
critical temperature
virtual doublereal thermalExpansionCoeff() const
Return the volumetric thermal expansion coefficient. Units: 1/K.
ThermoPhase * duplMyselfAsThermoPhase() const
Duplication function.
void setTPXState() const
Sets the state using a TPX::TV call.
virtual void getChemPotentials(doublereal *mu) const
Get the species chemical potentials. Units: J/kmol.
virtual doublereal standardConcentration(size_t k=0) const
Return the standard concentration for the kth species.
virtual doublereal pressure() const
Return the thermodynamic pressure (Pa).
virtual std::string report(bool show_thermo=true, doublereal threshold=1e-14) const
returns a summary of the state of the phase as a string
virtual doublereal critDensity() const
critical density
virtual void getEnthalpy_RT(doublereal *hrt) const
Get the nondimensional Enthalpy functions for the species at their standard states at the current T a...
virtual doublereal enthalpy_mole() const
Molar enthalpy. Units: J/kmol.
virtual void getStandardChemPotentials(doublereal *mu) const
Get the array of chemical potentials at unit activity for the species at their standard states at the...
virtual void setParametersFromXML(const XML_Node &eosdata)
Set equation of state parameter values from XML entries.
virtual void setState_HP(doublereal h, doublereal p, doublereal tol=1.e-8)
Set the internally stored specific enthalpy (J/kg) and pressure (Pa) of the phase.
virtual doublereal entropy_mole() const
Molar entropy. Units: J/kmol/K.
PureFluidPhase & operator=(const PureFluidPhase &right)
Assignment operator.
virtual void getGibbs_ref(doublereal *g) const
Returns the vector of the Gibbs function of the reference state at the current temperature of the sol...
virtual void initThermo()
Initialize the ThermoPhase object after all species have been set up.
virtual doublereal cv_mole() const
Molar heat capacity at constant volume. Units: J/kmol/K.
virtual void getPartialMolarIntEnergies(doublereal *ubar) const
Return an array of partial molar internal energies for the species in the mixture.
virtual void getPartialMolarEnthalpies(doublereal *hbar) const
Returns an array of partial molar enthalpies for the species in the mixture.
virtual void getActivityConcentrations(doublereal *c) const
This method returns an array of generalized concentrations.
virtual doublereal cp_mole() const
Molar heat capacity at constant pressure. Units: J/kmol/K.
void Set(tpx::PropertyPair::type n, double x, double y) const
Main call to the tpx level to set the state of the system.
Header file for class ThermoPhase, the base class for phases with thermodynamic properties, and the text for the Module thermoprops (see Thermodynamic Properties and class ThermoPhase).
virtual doublereal satPressure(doublereal t)
Return the saturation pressure given the temperature.
virtual void setState_SP(doublereal s, doublereal p, doublereal tol=1.e-8)
Set the specific entropy (J/kg/K) and pressure (Pa).
virtual void getGibbs_RT_ref(doublereal *grt) const
Returns the vector of nondimensional Gibbs Free Energies of the reference state at the current temper...
virtual doublereal vaporFraction() const
Return the fraction of vapor at the current conditions.
virtual void getEntropy_R_ref(doublereal *er) const
Returns the vector of nondimensional entropies of the reference state at the current temperature of t...
virtual int eosType() const
Equation of state type.