Cantera  2.1.2
ConstDensityThermo.h
Go to the documentation of this file.
1 /**
2  * @file ConstDensityThermo.h
3  * Header for a Thermo manager for incompressible ThermoPhases
4  * (see \ref thermoprops and \link Cantera::ConstDensityThermo ConstDensityThermo\endlink).
5  */
6 /*
7  * Copyright 2002 California Institute of Technology
8  */
9 
10 #ifndef CT_CONSTRHOTHERMO_H
11 #define CT_CONSTRHOTHERMO_H
12 
13 #include "cantera/base/ct_defs.h"
14 #include "mix_defs.h"
15 #include "ThermoPhase.h"
16 #include "SpeciesThermo.h"
17 #include "cantera/base/utilities.h"
18 
19 namespace Cantera
20 {
21 
22 //! Overloads the virtual methods of class ThermoPhase to implement the
23 //! incompressible equation of state.
24 /**
25  * <b> Specification of Solution Thermodynamic Properties </b>
26  *
27  * The density is assumed to be constant, no matter what the concentration of the solution.
28  *
29  * @ingroup thermoprops
30  */
32 {
33 public:
34  //! Constructor.
36 
37  //! Copy Constructor
38  /*!
39  * @param right Object to be copied
40  */
42 
43  //! Assignment Operator
44  /*!
45  * @param right Object to be copied
46  */
48 
49  //! Duplication routine for objects which inherit from %ThermoPhase
50  /*!
51  * This virtual routine can be used to duplicate objects
52  * derived from %ThermoPhase even if the application only has
53  * a pointer to %ThermoPhase to work with.
54  */
55  virtual ThermoPhase* duplMyselfAsThermoPhase() const;
56 
57  //! Returns a constant corresponding to this class's equation of state
58  virtual int eosType() const;
59 
60  /// Molar enthalpy. Units: J/kmol.
61  virtual doublereal enthalpy_mole() const;
62 
63  /// Molar internal energy. Units: J/kmol.
64  virtual doublereal intEnergy_mole() const;
65 
66  /// Molar entropy. Units: J/kmol/K.
67  virtual doublereal entropy_mole() const;
68 
69  /// Molar Gibbs function. Units: J/kmol.
70  virtual doublereal gibbs_mole() const;
71 
72  /// Molar heat capacity at constant pressure. Units: J/kmol/K.
73  virtual doublereal cp_mole() const;
74 
75  /// Molar heat capacity at constant volume. Units: J/kmol/K.
76  virtual doublereal cv_mole() const;
77 
78  //! Return the thermodynamic pressure (Pa).
79  virtual doublereal pressure() const;
80 
81  //! Set the internally stored pressure (Pa) at constant
82  //! temperature and composition
83  /*!
84  * @param p input Pressure (Pa)
85  */
86  virtual void setPressure(doublereal p);
87 
88  //! This method returns an array of generalized concentrations
89  /*!
90  * \f$ C^a_k\f$ are defined such that \f$ a_k = C^a_k /
91  * C^0_k, \f$ where \f$ C^0_k \f$ is a standard concentration
92  * defined below and \f$ a_k \f$ are activities used in the
93  * thermodynamic functions. These activity (or generalized)
94  * concentrations are used
95  * by kinetics manager classes to compute the forward and
96  * reverse rates of elementary reactions. Note that they may
97  * or may not have units of concentration --- they might be
98  * partial pressures, mole fractions, or surface coverages,
99  * for example.
100  *
101  * @param c Output array of generalized concentrations. The
102  * units depend upon the implementation of the
103  * reaction rate expressions within the phase.
104  */
105  virtual void getActivityConcentrations(doublereal* c) const;
106 
107  //! Get the array of non-dimensional molar-based activity coefficients at
108  //! the current solution temperature, pressure, and solution concentration.
109  /*!
110  * @param ac Output vector of activity coefficients. Length: m_kk.
111  */
112  virtual void getActivityCoefficients(doublereal* ac) const;
113 
114  //! Get the species chemical potentials. Units: J/kmol.
115  /*!
116  * This function returns a vector of chemical potentials of the
117  * species in solution at the current temperature, pressure
118  * and mole fraction of the solution.
119  *
120  * @param mu Output vector of species chemical
121  * potentials. Length: m_kk. Units: J/kmol
122  */
123  virtual void getChemPotentials(doublereal* mu) const;
124 
125  //! Get the array of chemical potentials at unit activity for the species
126  //! at their standard states at the current <I>T</I> and <I>P</I> of the solution.
127  /*!
128  * These are the standard state chemical potentials \f$ \mu^0_k(T,P)
129  * \f$. The values are evaluated at the current
130  * temperature and pressure of the solution
131  *
132  * @param mu0 Output vector of chemical potentials.
133  * Length: m_kk.
134  */
135  virtual void getStandardChemPotentials(doublereal* mu0) const;
136 
137  //! Return the standard concentration for the kth species
138  /*!
139  * The standard concentration \f$ C^0_k \f$ used to normalize
140  * the activity (i.e., generalized) concentration. In many cases, this quantity
141  * will be the same for all species in a phase - for example,
142  * for an ideal gas \f$ C^0_k = P/\hat R T \f$. For this
143  * reason, this method returns a single value, instead of an
144  * array. However, for phases in which the standard
145  * concentration is species-specific (e.g. surface species of
146  * different sizes), this method may be called with an
147  * optional parameter indicating the species.
148  *
149  * @param k Optional parameter indicating the species. The default
150  * is to assume this refers to species 0.
151  * @return
152  * Returns the standard Concentration in units of m3 kmol-1.
153  */
154  virtual doublereal standardConcentration(size_t k=0) const;
155 
156  //! Natural logarithm of the standard concentration of the kth species.
157  /*!
158  * @param k index of the species (defaults to zero)
159  */
160  virtual doublereal logStandardConc(size_t k=0) const;
161 
162  //! Get the Gibbs functions for the standard
163  //! state of the species at the current <I>T</I> and <I>P</I> of the solution
164  /*!
165  * Units are Joules/kmol
166  * @param gpure Output vector of standard state gibbs free energies
167  * Length: m_kk.
168  */
169  virtual void getPureGibbs(doublereal* gpure) const {
170  const vector_fp& gibbsrt = gibbs_RT();
171  scale(gibbsrt.begin(), gibbsrt.end(), gpure, _RT());
172  }
173 
174  //! Get the nondimensional Enthalpy functions for the species
175  //! at their standard states at the current <I>T</I> and <I>P</I> of the solution.
176  /*!
177  * @param hrt Output vector of nondimensional standard state enthalpies.
178  * Length: m_kk.
179  */
180  void getEnthalpy_RT(doublereal* hrt) const {
181  const vector_fp& _h = enthalpy_RT();
182  std::copy(_h.begin(), _h.end(), hrt);
183  }
184 
185  //! Get the array of nondimensional Entropy functions for the
186  //! standard state species at the current <I>T</I> and <I>P</I> of the solution.
187  /*!
188  * @param sr Output vector of nondimensional standard state entropies.
189  * Length: m_kk.
190  */
191  void getEntropy_R(doublereal* sr) const {
192  const vector_fp& _s = entropy_R();
193  std::copy(_s.begin(), _s.end(), sr);
194  }
195 
196  //! Get the nondimensional Gibbs functions for the species
197  //! in their standard states at the current <I>T</I> and <I>P</I> of the solution.
198  /*!
199  * @param grt Output vector of nondimensional standard state gibbs free energies
200  * Length: m_kk.
201  */
202  virtual void getGibbs_RT(doublereal* grt) const {
203  const vector_fp& gibbsrt = gibbs_RT();
204  std::copy(gibbsrt.begin(), gibbsrt.end(), grt);
205  }
206 
207  //! Get the nondimensional Heat Capacities at constant
208  //! pressure for the species standard states
209  //! at the current <I>T</I> and <I>P</I> of the solution
210  /*!
211  * @param cpr Output vector of nondimensional standard state heat capacities
212  * Length: m_kk.
213  */
214  void getCp_R(doublereal* cpr) const {
215  const vector_fp& _cpr = cp_R();
216  std::copy(_cpr.begin(), _cpr.end(), cpr);
217  }
218 
219  //! Returns a reference to the vector of nondimensional
220  //! enthalpies of the reference state at the current temperature
221  //! of the solution and the reference pressure for the species.
222  const vector_fp& enthalpy_RT() const {
223  _updateThermo();
224  return m_h0_RT;
225  }
226 
227  //! Returns a reference to the vector of nondimensional
228  //! Gibbs Free Energies of the reference state at the current temperature
229  //! of the solution and the reference pressure for the species.
230  const vector_fp& gibbs_RT() const {
231  _updateThermo();
232  return m_g0_RT;
233  }
234 
235  //! Returns a reference to the vector of nondimensional
236  //! entropies of the reference state at the current temperature
237  //! of the solution and the reference pressure for each species.
238  const vector_fp& entropy_R() const {
239  _updateThermo();
240  return m_s0_R;
241  }
242 
243  //! Returns a reference to the vector of nondimensional
244  //! constant pressure heat capacities of the reference state
245  //! at the current temperature of the solution
246  //! and reference pressure for each species.
247  const vector_fp& cp_R() const {
248  _updateThermo();
249  return m_cp0_R;
250  }
251 
252  //! Initialize the ThermoPhase object after all species have been set up
253  /*!
254  * @internal Initialize.
255  *
256  * This method is provided to allow
257  * subclasses to perform any initialization required after all
258  * species have been added. For example, it might be used to
259  * resize internal work arrays that must have an entry for
260  * each species. The base class implementation does nothing,
261  * and subclasses that do not require initialization do not
262  * need to overload this method. When importing a CTML phase
263  * description, this method is called from ThermoPhase::initThermoXML(),
264  * which is called from importPhase(),
265  * just prior to returning from function importPhase().
266  *
267  * @see importCTML.cpp
268  */
269  virtual void initThermo();
270 
271  //!This method is used by the ChemEquil equilibrium solver.
272  /*!
273  * It sets the state such that the chemical potentials satisfy
274  * \f[ \frac{\mu_k}{\hat R T} = \sum_m A_{k,m}
275  * \left(\frac{\lambda_m} {\hat R T}\right) \f] where
276  * \f$ \lambda_m \f$ is the element potential of element m. The
277  * temperature is unchanged. Any phase (ideal or not) that
278  * implements this method can be equilibrated by ChemEquil.
279  *
280  * @param lambda_RT Input vector of dimensionless element potentials
281  * The length is equal to nElements().
282  */
283  virtual void setToEquilState(const doublereal* lambda_RT);
284 
285  //! Set the equation of state parameters
286  /*!
287  * @internal
288  * The number and meaning of these depends on the subclass.
289  *
290  * @param n number of parameters
291  * @param c array of \a n coefficients
292  * @deprecated Use setDensity()
293  */
294  virtual void setParameters(int n, doublereal* const c) {
295  warn_deprecated("ConstDensityThermo::setParamters");
296  setDensity(c[0]);
297  }
298 
299  //! Get the equation of state parameters in a vector
300  /*!
301  * @internal
302  * The number and meaning of these depends on the subclass.
303  *
304  * @param n number of parameters
305  * @param c array of \a n coefficients
306  * @deprecated Use density()
307  */
308  virtual void getParameters(int& n, doublereal* const c) const {
309  warn_deprecated("ConstDensityThermo::getParameters");
310  double d = density();
311  c[0] = d;
312  n = 1;
313  }
314 
315  //! Set equation of state parameter values from XML entries.
316  /*!
317  *
318  * This method is called by function importPhase() in
319  * file importCTML.cpp when processing a phase definition in
320  * an input file. It should be overloaded in subclasses to set
321  * any parameters that are specific to that particular phase
322  * model. Note, this method is called before the phase is
323  * initialized with elements and/or species.
324  *
325  * @param eosdata An XML_Node object corresponding to
326  * the "thermo" entry for this phase in the input file.
327  */
328  virtual void setParametersFromXML(const XML_Node& eosdata);
329 
330 protected:
331  //! last value of the temperature processed by reference state
332  mutable doublereal m_tlast;
333 
334  //! Temporary storage for dimensionless reference state enthalpies
336 
337  //! Temporary storage for dimensionless reference state heat capacities
339 
340  //! Temporary storage for dimensionless reference state gibbs energies
342 
343  //! Temporary storage for dimensionless reference state entropies
344  mutable vector_fp m_s0_R;
345 
346  //! Temporary array containing internally calculated partial pressures
347  mutable vector_fp m_pp;
348 
349  //! Current pressure (Pa)
350  doublereal m_press;
351 
352 private:
353  //! Function to update the reference state thermo functions
354  void _updateThermo() const;
355 };
356 }
357 
358 #endif
virtual doublereal density() const
Density (kg/m^3).
Definition: Phase.h:534
void getEntropy_R(doublereal *sr) const
Get the array of nondimensional Entropy functions for the standard state species at the current T and...
vector_fp m_s0_R
Temporary storage for dimensionless reference state entropies.
Various templated functions that carry out common vector operations (see Templated Utility Functions)...
virtual doublereal standardConcentration(size_t k=0) const
Return the standard concentration for the kth species.
virtual void setToEquilState(const doublereal *lambda_RT)
This method is used by the ChemEquil equilibrium solver.
doublereal _RT() const
Return the Gas Constant multiplied by the current temperature.
Definition: ThermoPhase.h:977
virtual void getGibbs_RT(doublereal *grt) const
Get the nondimensional Gibbs functions for the species in their standard states at the current T and ...
virtual void getPureGibbs(doublereal *gpure) const
Get the Gibbs functions for the standard state of the species at the current T and P of the solution...
const vector_fp & cp_R() const
Returns a reference to the vector of nondimensional constant pressure heat capacities of the referenc...
virtual ThermoPhase * duplMyselfAsThermoPhase() const
Duplication routine for objects which inherit from ThermoPhase.
void getEnthalpy_RT(doublereal *hrt) const
Get the nondimensional Enthalpy functions for the species at their standard states at the current T a...
virtual void getActivityCoefficients(doublereal *ac) const
Get the array of non-dimensional molar-based activity coefficients at the current solution temperatur...
Class XML_Node is a tree-based representation of the contents of an XML file.
Definition: xml.h:100
vector_fp m_h0_RT
Temporary storage for dimensionless reference state enthalpies.
This file contains definitions of terms that are used in internal routines and are unlikely to need m...
virtual doublereal cp_mole() const
Molar heat capacity at constant pressure. Units: J/kmol/K.
void warn_deprecated(const std::string &method, const std::string &extra)
Print a warning indicating that method is deprecated.
Definition: global.cpp:76
Virtual base class for the calculation of multiple-species thermodynamic reference-state property man...
doublereal m_press
Current pressure (Pa)
vector_fp m_pp
Temporary array containing internally calculated partial pressures.
Base class for a phase with thermodynamic properties.
Definition: ThermoPhase.h:101
virtual void getChemPotentials(doublereal *mu) const
Get the species chemical potentials. Units: J/kmol.
virtual doublereal gibbs_mole() const
Molar Gibbs function. Units: J/kmol.
virtual void setPressure(doublereal p)
Set the internally stored pressure (Pa) at constant temperature and composition.
virtual void getStandardChemPotentials(doublereal *mu0) const
Get the array of chemical potentials at unit activity for the species at their standard states at the...
void _updateThermo() const
Function to update the reference state thermo functions.
virtual void getParameters(int &n, doublereal *const c) const
Get the equation of state parameters in a vector.
Overloads the virtual methods of class ThermoPhase to implement the incompressible equation of state...
virtual void setParameters(int n, doublereal *const c)
Set the equation of state parameters.
virtual doublereal cv_mole() const
Molar heat capacity at constant volume. Units: J/kmol/K.
const vector_fp & enthalpy_RT() const
Returns a reference to the vector of nondimensional enthalpies of the reference state at the current ...
virtual void setParametersFromXML(const XML_Node &eosdata)
Set equation of state parameter values from XML entries.
const vector_fp & gibbs_RT() const
Returns a reference to the vector of nondimensional Gibbs Free Energies of the reference state at the...
doublereal m_tlast
last value of the temperature processed by reference state
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
void scale(InputIter begin, InputIter end, OutputIter out, S scale_factor)
Multiply elements of an array by a scale factor.
Definition: utilities.h:154
virtual int eosType() const
Returns a constant corresponding to this class's equation of state.
virtual doublereal enthalpy_mole() const
Molar enthalpy. Units: J/kmol.
virtual void getActivityConcentrations(doublereal *c) const
This method returns an array of generalized concentrations.
void getCp_R(doublereal *cpr) const
Get the nondimensional Heat Capacities at constant pressure for the species standard states at the cu...
ConstDensityThermo & operator=(const ConstDensityThermo &right)
Assignment Operator.
const vector_fp & entropy_R() const
Returns a reference to the vector of nondimensional entropies of the reference state at the current t...
virtual doublereal entropy_mole() const
Molar entropy. Units: J/kmol/K.
virtual doublereal pressure() const
Return the thermodynamic pressure (Pa).
vector_fp m_g0_RT
Temporary storage for dimensionless reference state gibbs energies.
virtual doublereal intEnergy_mole() const
Molar internal energy. Units: J/kmol.
vector_fp m_cp0_R
Temporary storage for dimensionless reference state heat capacities.
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 logStandardConc(size_t k=0) const
Natural logarithm of the standard concentration of the kth species.
virtual void setDensity(const doublereal density_)
Set the internally stored density (kg/m^3) of the phase Note the density of a phase is an independent...
Definition: Phase.h:549
virtual void initThermo()
Initialize the ThermoPhase object after all species have been set up.