Cantera  2.1.2
IdealSolnGasVPSS.h
Go to the documentation of this file.
1 /**
2  * @file IdealSolnGasVPSS.h
3  * Definition file for a derived class of ThermoPhase that assumes either
4  * an ideal gas or ideal solution approximation and handles
5  * variable pressure standard state methods for calculating
6  * thermodynamic properties (see \ref thermoprops and
7  * class \link Cantera::IdealSolnGasVPSS IdealSolnGasVPSS\endlink).
8  */
9 /*
10  * Copyright (2005) Sandia Corporation. Under the terms of
11  * Contract DE-AC04-94AL85000 with Sandia Corporation, the
12  * U.S. Government retains certain rights in this software.
13  */
14 #ifndef CT_IDEALSOLNGASVPSS_H
15 #define CT_IDEALSOLNGASVPSS_H
16 
17 #include "VPStandardStateTP.h"
18 #include "VPSSMgr.h"
19 
20 namespace Cantera
21 {
22 class XML_Node;
23 class PDSS;
24 
25 /*!
26  * @name CONSTANTS
27  * Models for the Standard State of an IdealSolnPhase
28  */
29 //@{
30 const int cIdealSolnGasPhaseG = 6009;
31 const int cIdealSolnGasPhase0 = 6010;
32 const int cIdealSolnGasPhase1 = 6011;
33 const int cIdealSolnGasPhase2 = 6012;
34 
35 /**
36  * @ingroup thermoprops
37  *
38  * An ideal solution or an ideal gas approximation of a phase. Uses variable
39  * pressure standard state methods for calculating thermodynamic properties.
40  *
41  * @nosubgrouping
42  */
44 {
45 public:
46  /*!
47  * @name Constructors and Duplicators for %IdealSolnGasVPSS
48  */
49  //! @{
50 
51  /// Constructor.
53 
54  /// Create an object from an XML input file
55  IdealSolnGasVPSS(const std::string& infile, std::string id="");
56 
57  /// Copy Constructor.
59 
60  /// Assignment operator
62 
63  //! Duplication routine
64  virtual ThermoPhase* duplMyselfAsThermoPhase() const;
65 
66  //@}
67  //! @name Utilities (IdealSolnGasVPSS)
68  //@{
69  /**
70  * Equation of state type flag. The base class returns
71  * zero. Subclasses should define this to return a unique
72  * non-zero value. Constants defined for this purpose are
73  * listed in mix_defs.h.
74  */
75  virtual int eosType() const;
76 
77  //! @}
78  //! @name Molar Thermodynamic Properties
79  //! @{
80 
81  /// Molar enthalpy. Units: J/kmol.
82  doublereal enthalpy_mole() const;
83 
84  /// Molar internal energy. Units: J/kmol.
85  doublereal intEnergy_mole() const;
86 
87  /// Molar entropy. Units: J/kmol/K.
88  doublereal entropy_mole() const;
89 
90  /// Molar Gibbs function. Units: J/kmol.
91  doublereal gibbs_mole() const;
92 
93  /// Molar heat capacity at constant pressure. Units: J/kmol/K.
94  doublereal cp_mole() const;
95 
96  /// Molar heat capacity at constant volume. Units: J/kmol/K.
97  doublereal cv_mole() const;
98 
99  //! @}
100  //! @name Mechanical Properties
101  //! @{
102 
103  //! Set the pressure in the fluid
104  /*!
105  * @param p pressure in pascals.
106  */
107  void setPressure(doublereal p);
108 
109  //! Returns the isothermal compressibility. Units: 1/Pa.
110  /*!
111  * The isothermal compressibility is defined as
112  * \f[
113  * \kappa_T = -\frac{1}{v}\left(\frac{\partial v}{\partial P}\right)_T
114  * \f]
115  */
116  virtual doublereal isothermalCompressibility() const;
117 
118 protected:
119  /**
120  * Calculate the density of the mixture using the partial
121  * molar volumes and mole fractions as input
122  *
123  * The formula for this is
124  *
125  * \f[
126  * \rho = \frac{\sum_k{X_k W_k}}{\sum_k{X_k V_k}}
127  * \f]
128  *
129  * where \f$X_k\f$ are the mole fractions, \f$W_k\f$ are
130  * the molecular weights, and \f$V_k\f$ are the pure species
131  * molar volumes.
132  *
133  * Note, the basis behind this formula is that in an ideal
134  * solution the partial molar volumes are equal to the
135  * species standard state molar volumes.
136  * The species molar volumes may be functions
137  * of temperature and pressure.
138  */
139  virtual void calcDensity();
140  //! @}
141 
142 public:
143  //! This method returns an array of generalized concentrations
144  /*!
145  * \f$ C^a_k\f$ are defined such that \f$ a_k = C^a_k /
146  * C^0_k, \f$ where \f$ C^0_k \f$ is a standard concentration
147  * defined below and \f$ a_k \f$ are activities used in the
148  * thermodynamic functions. These activity (or generalized)
149  * concentrations are used
150  * by kinetics manager classes to compute the forward and
151  * reverse rates of elementary reactions. Note that they may
152  * or may not have units of concentration --- they might be
153  * partial pressures, mole fractions, or surface coverages,
154  * for example.
155  *
156  * @param c Output array of generalized concentrations. The
157  * units depend upon the implementation of the
158  * reaction rate expressions within the phase.
159  */
160  virtual void getActivityConcentrations(doublereal* c) const;
161 
162  //! Returns the standard concentration \f$ C^0_k \f$, which is used to normalize
163  //! the generalized concentration.
164  /*!
165  * This is defined as the concentration by which the generalized
166  * concentration is normalized to produce the activity.
167  * In many cases, this quantity will be the same for all species in a phase.
168  * Since the activity for an ideal gas mixture is
169  * simply the mole fraction, for an ideal gas \f$ C^0_k = P/\hat R T \f$.
170  *
171  * @param k Optional parameter indicating the species. The default
172  * is to assume this refers to species 0.
173  * @return
174  * Returns the standard Concentration in units of m3 kmol-1.
175  */
176  virtual doublereal standardConcentration(size_t k=0) const;
177 
178  //! Returns the natural logarithm of the standard
179  //! concentration of the kth species
180  /*!
181  * @param k index of the species. (defaults to zero)
182  */
183  virtual doublereal logStandardConc(size_t k=0) const;
184 
185  //! Returns the units of the standard and generalized concentrations.
186  /*!
187  * Note they have the same units, as their
188  * ratio is defined to be equal to the activity of the kth
189  * species in the solution, which is unitless.
190  *
191  * This routine is used in print out applications where the
192  * units are needed. Usually, MKS units are assumed throughout
193  * the program and in the XML input files.
194  *
195  * The base %ThermoPhase class assigns the default quantities
196  * of (kmol/m3) for all species.
197  * Inherited classes are responsible for overriding the default
198  * values if necessary.
199  *
200  * @param uA Output vector containing the units
201  * uA[0] = kmol units - default = 1
202  * uA[1] = m units - default = -nDim(), the number of spatial
203  * dimensions in the Phase class.
204  * uA[2] = kg units - default = 0;
205  * uA[3] = Pa(pressure) units - default = 0;
206  * uA[4] = Temperature units - default = 0;
207  * uA[5] = time units - default = 0
208  * @param k species index. Defaults to 0.
209  * @param sizeUA output int containing the size of the vector.
210  * Currently, this is equal to 6.
211  * @deprecated
212  */
213  virtual void getUnitsStandardConc(double* uA, int k = 0,
214  int sizeUA = 6) const;
215 
216  //! Get the array of non-dimensional activity coefficients at
217  //! the current solution temperature, pressure, and solution concentration.
218  /*!
219  * For ideal gases, the activity coefficients are all equal to one.
220  *
221  * @param ac Output vector of activity coefficients. Length: m_kk.
222  */
223  virtual void getActivityCoefficients(doublereal* ac) const;
224 
225 
226  /// @name Partial Molar Properties of the Solution
227  //@{
228 
229  //! Get the array of non-dimensional species chemical potentials
230  //! These are partial molar Gibbs free energies.
231  /*!
232  * \f$ \mu_k / \hat R T \f$.
233  * Units: unitless
234  *
235  * We close the loop on this function, here, calling
236  * getChemPotentials() and then dividing by RT. No need for child
237  * classes to handle.
238  *
239  * @param mu Output vector of non-dimensional species chemical potentials
240  * Length: m_kk.
241  */
242  void getChemPotentials_RT(doublereal* mu) const;
243 
244  //! Get the species chemical potentials. Units: J/kmol.
245  /*!
246  * This function returns a vector of chemical potentials of the
247  * species in solution at the current temperature, pressure
248  * and mole fraction of the solution.
249  *
250  * @param mu Output vector of species chemical
251  * potentials. Length: m_kk. Units: J/kmol
252  */
253  virtual void getChemPotentials(doublereal* mu) const;
254 
255  //! Get the species partial molar enthalpies. Units: J/kmol.
256  /*!
257  * @param hbar Output vector of species partial molar enthalpies.
258  * Length: m_kk. units are J/kmol.
259  */
260  virtual void getPartialMolarEnthalpies(doublereal* hbar) const;
261 
262  //! Get the species partial molar entropies. Units: J/kmol/K.
263  /*!
264  * @param sbar Output vector of species partial molar entropies.
265  * Length = m_kk. units are J/kmol/K.
266  */
267  virtual void getPartialMolarEntropies(doublereal* sbar) const;
268 
269  //! Get the species partial molar enthalpies. Units: J/kmol.
270  /*!
271  * @param ubar Output vector of species partial molar internal energies.
272  * Length = m_kk. units are J/kmol.
273  */
274  virtual void getPartialMolarIntEnergies(doublereal* ubar) const;
275 
276  //! Get the partial molar heat capacities Units: J/kmol/K
277  /*!
278  * @param cpbar Output vector of species partial molar heat capacities
279  * at constant pressure.
280  * Length = m_kk. units are J/kmol/K.
281  */
282  virtual void getPartialMolarCp(doublereal* cpbar) const;
283 
284  //! Get the species partial molar volumes. Units: m^3/kmol.
285  /*!
286  * @param vbar Output vector of species partial molar volumes.
287  * Length = m_kk. units are m^3/kmol.
288  */
289  virtual void getPartialMolarVolumes(doublereal* vbar) const;
290  //@}
291 
292 public:
293  //! @name Initialization Methods - For Internal use
294  /*!
295  * The following methods are used in the process of constructing
296  * the phase and setting its parameters from a specification in an
297  * input file. They are not normally used in application programs.
298  * To see how they are used, see files importCTML.cpp and
299  * ThermoFactory.cpp.
300  */
301  //@{
302 
303  //! Set equation of state parameter values from XML entries.
304  /*!
305  * This method is called by function importPhase in
306  * file importCTML.cpp when processing a phase definition in
307  * an input file. It should be overloaded in subclasses to set
308  * any parameters that are specific to that particular phase
309  * model.
310  *
311  * @param thermoNode An XML_Node object corresponding to
312  * the "thermo" entry for this phase in the input file.
313  */
314  virtual void setParametersFromXML(const XML_Node& thermoNode);
315 
316  //! @internal Initialize the object
317  /*!
318  * This method is provided to allow
319  * subclasses to perform any initialization required after all
320  * species have been added. For example, it might be used to
321  * resize internal work arrays that must have an entry for
322  * each species. The base class implementation does nothing,
323  * and subclasses that do not require initialization do not
324  * need to overload this method. When importing a CTML phase
325  * description, this method is called just prior to returning
326  * from function importPhase().
327  *
328  * @see importCTML.cpp
329  */
330  virtual void initThermo();
331 
332  //!This method is used by the ChemEquil equilibrium solver.
333  /*!
334  * It sets the state such that the chemical potentials satisfy
335  * \f[ \frac{\mu_k}{\hat R T} = \sum_m A_{k,m}
336  * \left(\frac{\lambda_m} {\hat R T}\right) \f] where
337  * \f$ \lambda_m \f$ is the element potential of element m. The
338  * temperature is unchanged. Any phase (ideal or not) that
339  * implements this method can be equilibrated by ChemEquil.
340  *
341  * @param lambda_RT Input vector of dimensionless element potentials
342  * The length is equal to nElements().
343  */
344  void setToEquilState(const doublereal* lambda_RT);
345 
346  //! Initialize a ThermoPhase object, potentially reading activity
347  //! coefficient information from an XML database.
348  /*!
349  * This routine initializes the lengths in the current object and
350  * then calls the parent routine.
351  * This method is provided to allow
352  * subclasses to perform any initialization required after all
353  * species have been added. For example, it might be used to
354  * resize internal work arrays that must have an entry for
355  * each species. The base class implementation does nothing,
356  * and subclasses that do not require initialization do not
357  * need to overload this method. When importing a CTML phase
358  * description, this method is called just prior to returning
359  * from function importPhase().
360  *
361  * @param phaseNode This object must be the phase node of a
362  * complete XML tree
363  * description of the phase, including all of the
364  * species data. In other words while "phase" must
365  * point to an XML phase object, it must have
366  * sibling nodes "speciesData" that describe
367  * the species in the phase.
368  * @param id ID of the phase. If nonnull, a check is done
369  * to see if phaseNode is pointing to the phase
370  * with the correct id.
371  */
372  virtual void initThermoXML(XML_Node& phaseNode, const std::string& id);
373 
374 private:
375  //! @internal Initialize the internal lengths in this object.
376  /*!
377  * Note this is not a virtual function and only handles this object
378  */
379  void initLengths();
380 
381  //@}
382 
383 protected:
384  //! boolean indicating what ideal solution this is
385  /*!
386  * - 1 = ideal gas
387  * - 0 = ideal soln
388  */
390 
391  //! form of the generalized concentrations
392  /*!
393  * - 0 unity
394  * - 1 1/V_k
395  * - 2 1/V_0
396  */
397  int m_formGC;
398 
399  //! Temporary storage - length = m_kk.
401 };
402 }
403 
404 #endif
virtual void getPartialMolarEnthalpies(doublereal *hbar) const
Get the species partial molar enthalpies. Units: J/kmol.
virtual void getPartialMolarEntropies(doublereal *sbar) const
Get the species partial molar entropies. Units: J/kmol/K.
Declaration file for a virtual base class that manages the calculation of standard state properties f...
virtual doublereal isothermalCompressibility() const
Returns the isothermal compressibility. Units: 1/Pa.
Class XML_Node is a tree-based representation of the contents of an XML file.
Definition: xml.h:100
doublereal cv_mole() const
Molar heat capacity at constant volume. Units: J/kmol/K.
virtual ThermoPhase * duplMyselfAsThermoPhase() const
Duplication routine.
vector_fp m_pp
Temporary storage - length = m_kk.
doublereal cp_mole() const
Molar heat capacity at constant pressure. Units: J/kmol/K.
Base class for a phase with thermodynamic properties.
Definition: ThermoPhase.h:101
virtual void getPartialMolarCp(doublereal *cpbar) const
Get the partial molar heat capacities Units: J/kmol/K.
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.
doublereal enthalpy_mole() const
Molar enthalpy. Units: J/kmol.
int m_formGC
form of the generalized concentrations
This is a filter class for ThermoPhase that implements some prepatory steps for efficiently handling ...
virtual void getActivityConcentrations(doublereal *c) const
This method returns an array of generalized concentrations.
IdealSolnGasVPSS & operator=(const IdealSolnGasVPSS &)
Assignment operator.
virtual void getPartialMolarVolumes(doublereal *vbar) const
Get the species partial molar volumes. Units: m^3/kmol.
int m_idealGas
boolean indicating what ideal solution this is
void setPressure(doublereal p)
Set the pressure in the fluid.
virtual void setParametersFromXML(const XML_Node &thermoNode)
Set equation of state parameter values from XML entries.
Header file for a derived class of ThermoPhase that handles variable pressure standard state methods ...
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:165
virtual int eosType() const
Equation of state type flag.
virtual void getActivityCoefficients(doublereal *ac) const
Get the array of non-dimensional activity coefficients at the current solution temperature, pressure, and solution concentration.
virtual void getChemPotentials(doublereal *mu) const
Get the species chemical potentials. Units: J/kmol.
virtual doublereal logStandardConc(size_t k=0) const
Returns the natural logarithm of the standard concentration of the kth species.
virtual void getUnitsStandardConc(double *uA, int k=0, int sizeUA=6) const
Returns the units of the standard and generalized concentrations.
doublereal intEnergy_mole() const
Molar internal energy. Units: J/kmol.
doublereal entropy_mole() const
Molar entropy. Units: J/kmol/K.
virtual void initThermoXML(XML_Node &phaseNode, const std::string &id)
Initialize a ThermoPhase object, potentially reading activity coefficient information from an XML dat...
doublereal gibbs_mole() const
Molar Gibbs function. Units: J/kmol.
virtual void getPartialMolarIntEnergies(doublereal *ubar) const
Get the species partial molar enthalpies. Units: J/kmol.
void setToEquilState(const doublereal *lambda_RT)
This method is used by the ChemEquil equilibrium solver.
An ideal solution or an ideal gas approximation of a phase.
void getChemPotentials_RT(doublereal *mu) const
Get the array of non-dimensional species chemical potentials These are partial molar Gibbs free energ...