Cantera  2.1.2
SingleSpeciesTP.h
Go to the documentation of this file.
1 /**
2  * @file SingleSpeciesTP.h
3  * Header for the %SingleSpeciesTP class, which is a filter class for %ThermoPhase,
4  * that eases the construction of single species phases
5  * ( see \ref thermoprops and class \link Cantera::SingleSpeciesTP SingleSpeciesTP\endlink).
6  *
7  */
8 /*
9  * Copyright (2005) Sandia Corporation. Under the terms of
10  * Contract DE-AC04-94AL85000 with Sandia Corporation, the
11  * U.S. Government retains certain rights in this software.
12  */
13 #ifndef CT_SINGLESPECIESTP_H
14 #define CT_SINGLESPECIESTP_H
15 
16 #include "ThermoPhase.h"
17 
18 
19 namespace Cantera
20 {
21 
22 /**
23  * @ingroup thermoprops
24  *
25  * The SingleSpeciesTP class is a filter class for ThermoPhase.
26  * What it does is to simplify the construction of ThermoPhase
27  * objects by assuming that the phase consists of one and
28  * only one type of species. In other words, it's a stoichiometric
29  * phase. However, no assumptions are made concerning the
30  * thermodynamic functions or the equation of state of the
31  * phase. Therefore it's an incomplete description of
32  * the thermodynamics. The complete description must be
33  * made in a derived class of %SingleSpeciesTP.
34  *
35  * Several different groups of thermodynamic functions are resolved
36  * at this level by this class. For example, All partial molar property
37  * routines call their single species standard state equivalents.
38  * All molar solution thermodynamic routines call the single species
39  * standard state equivalents.
40  * Activities routines are resolved at this level, as there is only
41  * one species.
42  *
43  * It is assumed that the reference state thermodynamics may be
44  * obtained by a pointer to a populated species thermodynamic property
45  * manager class (see ThermoPhase::m_spthermo). How to relate pressure
46  * changes to the reference state thermodynamics is again left open
47  * to implementation.
48  *
49  * Mole fraction and Mass fraction vectors are assumed to be equal
50  * to x[0] = 1 y[0] = 1, respectively. Simplifications to the interface
51  * of setState_TPY() and setState_TPX() functions result and are made
52  * within the class.
53  *
54  * Note, this class can handle the thermodynamic description of one
55  * phase of one species. It can not handle the description of phase
56  * equilibrium between two phases of a stoichiometric compound
57  * (e.g. water liquid and water vapor, below the critical point).
58  * However, it may be used to describe the thermodynamics of one phase
59  * of such a compound even past the phase equilibrium point, up to the
60  * point where the phase itself ceases to be a stable phase.
61  *
62  * This class doesn't do much at the initialization level.
63  * Its SingleSpeciesTP::initThermo()
64  * member does check that one and only one species has been defined
65  * to occupy the phase.
66  *
67  * \nosubgrouping
68  */
70 {
71 public:
72  //! Base empty constructor.
74 
75  //! Copy constructor
76  /*!
77  * @param right Object to be copied
78  */
79  SingleSpeciesTP(const SingleSpeciesTP& right);
80 
81  //! Assignment operator
82  /*!
83  * @param right Object to be copied
84  */
86 
87  //! Duplication function
88  /*!
89  * This virtual function is used to create a duplicate of the
90  * current phase. It's used to duplicate the phase when given
91  * a ThermoPhase pointer to the phase.
92  *
93  * @return It returns a ThermoPhase pointer.
94  */
96 
97  /**
98  * Returns the equation of state type flag.
99  * This is a modified base class.
100  * Therefore, if not overridden in derivied classes,
101  * this call will throw an exception.
102  */
103  virtual int eosType() const;
104 
105  /**
106  * @name Molar Thermodynamic Properties of the Solution
107  *
108  * These functions are resolved at this level, by reference
109  * to the partial molar functions and standard state
110  * functions for species 0. Derived classes don't need
111  * to supply entries for these functions.
112  * @{
113  */
114 
115  /// Molar enthalpy. Units: J/kmol.
116  /*!
117  * This function is resolved here by calling the standard state
118  * thermo function.
119  */
120  doublereal enthalpy_mole() const;
121 
122  /// Molar internal energy. Units: J/kmol.
123  /*!
124  * This function is resolved here by calling the standard state
125  * thermo function.
126  */
127  doublereal intEnergy_mole() const;
128 
129  /// Molar entropy. Units: J/kmol/K.
130  /*!
131  * This function is resolved here by calling the standard state
132  * thermo function.
133  */
134  doublereal entropy_mole() const;
135 
136  /// Molar Gibbs function. Units: J/kmol.
137  /*!
138  * This function is resolved here by calling the standard state
139  * thermo function.
140  */
141  doublereal gibbs_mole() const;
142 
143  /// Molar heat capacity at constant pressure. Units: J/kmol/K.
144  /*!
145  * This function is resolved here by calling the standard state
146  * thermo function.
147  */
148  doublereal cp_mole() const;
149 
150  /// Molar heat capacity at constant volume. Units: J/kmol/K.
151  /*!
152  * This function is resolved here by calling the standard state
153  * thermo function.
154  */
155  doublereal cv_mole() const;
156 
157  /**
158  * @}
159  * @name Activities, Standard State, and Activity Concentrations
160  *
161  * The activity \f$a_k\f$ of a species in solution is
162  * related to the chemical potential by \f[ \mu_k = \mu_k^0(T)
163  * + \hat R T \log a_k. \f] The quantity \f$\mu_k^0(T)\f$ is
164  * the chemical potential at unit activity, which depends only
165  * on temperature.
166  * @{
167  */
168 
169  /**
170  * Get the array of non-dimensional activities at
171  * the current solution temperature, pressure, and
172  * solution concentration.
173  *
174  * We redefine this function to just return 1.0 here.
175  *
176  * @param a Output vector of activities. Length: 1.
177  */
178  virtual void getActivities(doublereal* a) const {
179  a[0] = 1.0;
180  }
181 
182  /**
183  * Get the array of non-dimensional activity coefficients at
184  * the current solution temperature, pressure, and
185  * solution concentration.
186  *
187  * @param ac Output vector of activity coefficients. Length: 1.
188  */
189  virtual void getActivityCoefficients(doublereal* ac) const {
190  if (m_kk == 1) {
191  ac[0] = 1.0;
192  } else {
193  err("getActivityCoefficients");
194  }
195  }
196 
197  //@}
198  /// @name Partial Molar Properties of the Solution
199  ///
200  /// These functions are resolved at this level, by reference
201  /// to the partial molar functions and standard state
202  /// functions for species 0. Derived classes don't need
203  /// to supply entries for these functions.
204  //@{
205 
206  //! Get the array of non-dimensional species chemical potentials
207  //! These are partial molar Gibbs free energies.
208  /*!
209  * These are the phase, partial molar, and the standard state
210  * dimensionless chemical potentials.
211  * \f$ \mu_k / \hat R T \f$.
212  *
213  * Units: unitless
214  *
215  * @param murt On return, Contains the chemical potential / RT of the single species
216  * and the phase. Units are unitless. Length = 1
217  */
218  void getChemPotentials_RT(doublereal* murt) const;
219 
220  //! Get the array of chemical potentials
221  /*!
222  * These are the phase, partial molar, and the standard state chemical potentials.
223  * \f$ \mu(T,P) = \mu^0_k(T,P) \f$.
224  *
225  * @param mu On return, Contains the chemical potential of the single species
226  * and the phase. Units are J / kmol . Length = 1
227  */
228  void getChemPotentials(doublereal* mu) const;
229 
230  //! Get the species electrochemical potentials. Units: J/kmol.
231  /*!
232  * This method adds a term \f$ Fz_k \phi_k \f$ to
233  * each chemical potential.
234  *
235  * This is resolved here. A single species phase
236  * is not allowed to have anything other than a zero charge.
237  *
238  * @param mu On return, Contains the electrochemical potential of the single species
239  * and the phase. Units J/kmol . Length = 1
240  */
241  void getElectrochemPotentials(doublereal* mu) const;
242 
243  //! Get the species partial molar enthalpies. Units: J/kmol.
244  /*!
245  * These are the phase enthalpies. \f$ h_k \f$.
246  *
247  * @param hbar Output vector of species partial molar enthalpies.
248  * Length: 1. units are J/kmol.
249  */
250  void getPartialMolarEnthalpies(doublereal* hbar) const;
251 
252  //! Get the species partial molar internal energies. Units: J/kmol.
253  /*!
254  * These are the phase internal energies. \f$ u_k \f$.
255  *
256  * @param ubar On return, Contains the internal energy of the single species
257  * and the phase. Units are J / kmol . Length = 1
258  */
259  virtual void getPartialMolarIntEnergies(doublereal* ubar) const;
260 
261  //! Get the species partial molar entropy. Units: J/kmol K.
262  /*!
263  * This is the phase entropy. \f$ s(T,P) = s_o(T,P) \f$.
264  *
265  * @param sbar On return, Contains the entropy of the single species
266  * and the phase. Units are J / kmol / K . Length = 1
267  */
268  void getPartialMolarEntropies(doublereal* sbar) const;
269 
270  //! Get the species partial molar Heat Capacities. Units: J/ kmol /K.
271  /*!
272  * This is the phase heat capacity. \f$ Cp(T,P) = Cp_o(T,P) \f$.
273  *
274  * @param cpbar On return, Contains the heat capacity of the single species
275  * and the phase. Units are J / kmol / K . Length = 1
276  */
277  void getPartialMolarCp(doublereal* cpbar) const;
278 
279  //! Get the species partial molar volumes. Units: m^3/kmol.
280  /*!
281  * This is the phase molar volume. \f$ V(T,P) = V_o(T,P) \f$.
282  *
283  * @param vbar On return, Contains the molar volume of the single species
284  * and the phase. Units are m^3 / kmol. Length = 1
285  */
286  void getPartialMolarVolumes(doublereal* vbar) const;
287 
288  //@}
289  /// @name Properties of the Standard State of the Species in the Solution
290  /// These functions are the primary way real properties are
291  /// supplied to derived thermodynamics classes of SingleSpeciesTP.
292  /// These functions must be supplied in derived classes. They
293  /// are not resolved at the SingleSpeciesTP level.
294  //@{
295 
296  /**
297  * Get the dimensional Gibbs functions for the standard
298  * state of the species at the current T and P.
299  *
300  * @param gpure returns a vector of size 1, containing the Gibbs function
301  * Units: J/kmol.
302  */
303  void getPureGibbs(doublereal* gpure) const;
304 
305  //! Get the molar volumes of each species in their standard
306  //! states at the current <I>T</I> and <I>P</I> of the solution.
307  /*!
308  * units = m^3 / kmol
309  *
310  * We resolve this function at this level, by assigning
311  * the molecular weight divided by the phase density
312  *
313  * @param vbar On output this contains the standard volume of the species
314  * and phase (m^3/kmol). Vector of length 1
315  */
316  void getStandardVolumes(doublereal* vbar) const;
317 
318  //@}
319  /// @name Thermodynamic Values for the Species Reference State
320  ///
321  /// Almost all functions in this group are resolved by this
322  /// class. It is assumed that the m_spthermo species thermo
323  /// pointer is populated and yields the reference state thermodynamics
324  /// The internal energy function is not given by this
325  /// class, since it would involve a specification of the
326  /// equation of state.
327  //@{
328 
329 #ifdef H298MODIFY_CAPABILITY
330 
331  //! Modify the value of the 298 K Heat of Formation of one species in the phase (J kmol-1)
332  /*!
333  * The 298K heat of formation is defined as the enthalpy change to create the standard state
334  * of the species from its constituent elements in their standard states at 298 K and 1 bar.
335  *
336  * @param k Species k
337  * @param Hf298New Specify the new value of the Heat of Formation at 298K and 1 bar
338  */
339  virtual void modifyOneHf298SS(const size_t& k, const doublereal Hf298New) {
340  m_spthermo->modifyOneHf298(k, Hf298New);
341  m_tlast += 0.0001234;
342  }
343 #endif
344 
345  /*!
346  * Returns the vector of nondimensional
347  * enthalpies of the reference state at the current temperature
348  * of the solution and the reference pressure for the species.
349  *
350  * This function is resolved in this class. It is assumed that the m_spthermo species thermo
351  * pointer is populated and yields the reference state.
352  *
353  * @param hrt Output vector containing the nondimensional reference state enthalpies
354  * Length: m_kk.
355  */
356  virtual void getEnthalpy_RT_ref(doublereal* hrt) const;
357 
358  /*!
359  * Returns the vector of nondimensional
360  * enthalpies of the reference state at the current temperature
361  * of the solution and the reference pressure for the species.
362  *
363  * This function is resolved in this class. It is assumed that the m_spthermo species thermo
364  * pointer is populated and yields the reference state.
365  *
366  * @param grt Output vector containing the nondimensional reference state
367  * Gibbs Free energies. Length: m_kk.
368  */
369  virtual void getGibbs_RT_ref(doublereal* grt) const;
370 
371  /*!
372  * Returns the vector of the
373  * gibbs function of the reference state at the current temperature
374  * of the solution and the reference pressure for the species.
375  * units = J/kmol
376  *
377  * This function is resolved in this class. It is assumed that the m_spthermo species thermo
378  * pointer is populated and yields the reference state.
379  *
380  * @param g Output vector containing the reference state
381  * Gibbs Free energies. Length: m_kk. Units: J/kmol.
382  */
383  virtual void getGibbs_ref(doublereal* g) const;
384 
385  /*!
386  * Returns the vector of nondimensional
387  * entropies of the reference state at the current temperature
388  * of the solution and the reference pressure for each species.
389  *
390  * This function is resolved in this class. It is assumed that the m_spthermo species thermo
391  * pointer is populated and yields the reference state.
392  *
393  * @param er Output vector containing the nondimensional reference state
394  * entropies. Length: m_kk.
395  */
396  virtual void getEntropy_R_ref(doublereal* er) const;
397 
398  /*!
399  * Returns the vector of nondimensional
400  * constant pressure heat capacities of the reference state
401  * at the current temperature of the solution
402  * and reference pressure for each species.
403  *
404  * This function is resolved in this class. It is assumed that the m_spthermo species thermo
405  * pointer is populated and yields the reference state.
406  *
407  * @param cprt Output vector of nondimensional reference state
408  * heat capacities at constant pressure for the species.
409  * Length: m_kk
410  */
411  virtual void getCp_R_ref(doublereal* cprt) const;
412 
413  /**
414  * @name Setting the State
415  *
416  * These methods set all or part of the thermodynamic state.
417  * @{
418  */
419 
420  //! Set the temperature (K), pressure (Pa), and mole fractions.
421  /*!
422  * Note, the mole fractions are set to X[0] = 1.0.
423  * Setting the pressure may involve the solution of a nonlinear equation.
424  *
425  * @param t Temperature (K)
426  * @param p Pressure (Pa)
427  * @param x Vector of mole fractions.
428  * Length is equal to m_kk.
429  */
430  void setState_TPX(doublereal t, doublereal p, const doublereal* x);
431 
432  //! Set the temperature (K), pressure (Pa), and mole fractions.
433  /*!
434  * Note, the mole fractions are set to X[0] = 1.0.
435  * Setting the pressure may involve the solution of a nonlinear equation.
436  *
437  * @param t Temperature (K)
438  * @param p Pressure (Pa)
439  * @param x String containing a composition map of the mole fractions. Species not in
440  * the composition map are assumed to have zero mole fraction
441  */
442  void setState_TPX(doublereal t, doublereal p, compositionMap& x);
443 
444  //! Set the temperature (K), pressure (Pa), and mole fractions.
445  /*!
446  * Note, the mole fractions are set to X[0] = 1.0.
447  * Setting the pressure may involve the solution of a nonlinear equation.
448  *
449  * @param t Temperature (K)
450  * @param p Pressure (Pa)
451  * @param x String containing a composition map of the mole fractions. Species not in
452  * the composition map are assumed to have zero mole fraction
453  */
454  void setState_TPX(doublereal t, doublereal p, const std::string& x);
455 
456  //! Set the internally stored temperature (K), pressure (Pa), and mass fractions of the phase.
457  /*!
458  * Note, the mass fractions are set to Y[0] = 1.0.
459  * Setting the pressure may involve the solution of a nonlinear equation.
460  *
461  * @param t Temperature (K)
462  * @param p Pressure (Pa)
463  * @param y Vector of mass fractions.
464  * Length is equal to m_kk.
465  */
466  void setState_TPY(doublereal t, doublereal p, const doublereal* y);
467 
468  //! Set the internally stored temperature (K), pressure (Pa), and mass fractions of the phase
469  /*!
470  * Note, the mass fractions are set to Y[0] = 1.0.
471  * Setting the pressure may involve the solution of a nonlinear equation.
472  *
473  * @param t Temperature (K)
474  * @param p Pressure (Pa)
475  * @param y Composition map of mass fractions. Species not in
476  * the composition map are assumed to have zero mass fraction
477  */
478  void setState_TPY(doublereal t, doublereal p, compositionMap& y);
479 
480  //! Set the internally stored temperature (K), pressure (Pa), and mass fractions of the phase
481  /*!
482  * Note, the mass fractions are set to Y[0] = 1.0.
483  * Setting the pressure may involve the solution of a nonlinear equation.
484  *
485  * @param t Temperature (K)
486  * @param p Pressure (Pa)
487  * @param y String containing a composition map of the mass fractions. Species not in
488  * the composition map are assumed to have zero mass fraction
489  */
490  void setState_TPY(doublereal t, doublereal p, const std::string& y);
491 
492  //! Set the pressure (Pa) and mole fractions.
493  /*!
494  * Note, the mole fractions are set to X[0] = 1.0.
495  * Setting the pressure may involve the solution of a nonlinear equation.
496  *
497  * @param p Pressure (Pa)
498  * @param x Vector of mole fractions.
499  * Length is equal to m_kk.
500  */
501  void setState_PX(doublereal p, doublereal* x);
502 
503  //! Set the internally stored pressure (Pa) and mass fractions.
504  /*!
505  * Note, the mass fractions are set to Y[0] = 1.0.
506  * Note, the temperature is held constant during this operation.
507  * Setting the pressure may involve the solution of a nonlinear equation.
508  *
509  * @param p Pressure (Pa)
510  * @param y Vector of mass fractions.
511  * Length is equal to m_kk.
512  */
513  void setState_PY(doublereal p, doublereal* y);
514 
515  //! Set the internally stored specific enthalpy (J/kg) and pressure (Pa) of the phase.
516  /*!
517  * @param h Specific enthalpy (J/kg)
518  * @param p Pressure (Pa)
519  * @param tol Optional parameter setting the tolerance of the
520  * calculation.
521  */
522  virtual void setState_HP(doublereal h, doublereal p,
523  doublereal tol = 1.e-8);
524 
525  //! Set the specific internal energy (J/kg) and specific volume (m^3/kg).
526  /*!
527  * This function fixes the internal state of the phase so that
528  * the specific internal energy and specific volume have the value of the input parameters.
529  *
530  * @param u specific internal energy (J/kg)
531  * @param v specific volume (m^3/kg).
532  * @param tol Optional parameter setting the tolerance of the
533  * calculation.
534  */
535  virtual void setState_UV(doublereal u, doublereal v,
536  doublereal tol = 1.e-8);
537 
538  //! Set the specific entropy (J/kg/K) and pressure (Pa).
539  /*!
540  * This function fixes the internal state of the phase so that
541  * the specific entropy and the pressure have the value of the input parameters.
542  *
543  * @param s specific entropy (J/kg/K)
544  * @param p specific pressure (Pa).
545  * @param tol Optional parameter setting the tolerance of the
546  * calculation.
547  */
548  virtual void setState_SP(doublereal s, doublereal p,
549  doublereal tol = 1.e-8);
550 
551  //! Set the specific entropy (J/kg/K) and specific volume (m^3/kg).
552  /*!
553  * This function fixes the internal state of the phase so that
554  * the specific entropy and specific volume have the value of the input parameters.
555  *
556  * @param s specific entropy (J/kg/K)
557  * @param v specific volume (m^3/kg).
558  * @param tol Optional parameter setting the tolerance of the
559  * calculation.
560  */
561  virtual void setState_SV(doublereal s, doublereal v,
562  doublereal tol = 1.e-8);
563 
564  /**
565  * @internal
566  * Set equation of state parameters. The number and meaning of
567  * these depends on the subclass.
568  * @param n number of parameters
569  * @param c array of n coefficients
570  * @deprecated Unimplemented
571  */
572  virtual void setParameters(int n, doublereal* const c) {
573  warn_deprecated("SingleSpeciesTP::setParameters");
574  }
575 
576  //! @deprecated Unimplemented
577  virtual void getParameters(int& n, doublereal* const c) const {
578  warn_deprecated("SingleSpeciesTP::getParameters");
579  }
580 
581  /**
582  * Set equation of state parameter values from XML
583  * entries. This method is called by function importPhase in
584  * file importCTML.cpp when processing a phase definition in
585  * an input file. It should be overloaded in subclasses to set
586  * any parameters that are specific to that particular phase
587  * model.
588  *
589  * @param eosdata An XML_Node object corresponding to
590  * the "thermo" entry for this phase in the input file.
591  */
592  virtual void setParametersFromXML(const XML_Node& eosdata) {}
593 
594  //@}
595  /// @name Saturation properties.
596  /// These methods are only implemented by subclasses that
597  /// implement full liquid-vapor equations of state.
598  ///
599  virtual doublereal satTemperature(doublereal p) const {
600  err("satTemperature");
601  return -1.0;
602  }
603 
604  virtual doublereal satPressure(doublereal t) {
605  err("satPressure");
606  return -1.0;
607  }
608 
609  virtual doublereal vaporFraction() const {
610  err("vaprFraction");
611  return -1.0;
612  }
613 
614  virtual void setState_Tsat(doublereal t, doublereal x) {
615  err("setState_sat");
616  }
617 
618  virtual void setState_Psat(doublereal p, doublereal x) {
619  err("setState_sat");
620  }
621 
622  //@}
623 
624  /**
625  * @internal Initialize.
626  *
627  * This method is provided to allow
628  * subclasses to perform any initialization required after all
629  * species have been added. For example, it might be used to
630  * resize internal work arrays that must have an entry for
631  * each species. When importing a CTML phase
632  * description, this method is called just prior to returning
633  * from function importPhase().
634  *
635  * Inheriting objects should call this function
636  *
637  * This version sets the mole fraction vector to x[0] = 1.0, and then
638  * calls the ThermoPhase::initThermo() function.
639  *
640  * @see importCTML.cpp
641  */
642  virtual void initThermo();
643 
644 protected:
645  //! The current pressure of the solution (Pa)
646  /*!
647  * It gets initialized to 1 atm.
648  */
649  doublereal m_press;
650 
651  /*!
652  * Reference pressure (Pa) must be the same for all species
653  * - defaults to 1 atm.
654  */
655  doublereal m_p0;
656 
657  //! Last temperature used to evaluate the thermodynamic polynomial.
658  mutable doublereal m_tlast;
659 
660  //! Dimensionless enthalpy at the (mtlast, m_p0)
662  //! Dimensionless heat capacity at the (mtlast, m_p0)
664  //! Dimensionless entropy at the (mtlast, m_p0)
665  mutable vector_fp m_s0_R;
666 
667 protected:
668  /**
669  * @internal
670  * This crucial internal routine calls the species thermo
671  * update program to calculate new species Cp0, H0, and
672  * S0 whenever the temperature has changed.
673  */
674  void _updateThermo() const;
675 
676 private:
677 
678  //! Error return for unhandled cases.
679  /*!
680  * It's used when this class doesn't have an answer for the question given
681  * to it, because the derived class isn't overriding a function.
682  *
683  * @param msg String message
684  */
685  doublereal err(const std::string& msg) const;
686 };
687 
688 }
689 
690 #endif
virtual doublereal satPressure(doublereal t)
Return the saturation pressure given the temperature.
std::map< std::string, doublereal > compositionMap
Map connecting a string name with a double.
Definition: ct_defs.h:162
virtual void setParametersFromXML(const XML_Node &eosdata)
Set equation of state parameter values from XML entries.
doublereal cv_mole() const
Molar heat capacity at constant volume. Units: J/kmol/K.
ThermoPhase * duplMyselfAsThermoPhase() const
Duplication function.
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 satTemperature(doublereal p) const
Return the saturation temperature given the pressure.
virtual void getEnthalpy_RT_ref(doublereal *hrt) const
doublereal err(const std::string &msg) const
Error return for unhandled cases.
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).
virtual void modifyOneHf298SS(const int k, const doublereal Hf298New)
Modify the value of the 298 K Heat of Formation of one species in the phase (J kmol-1) ...
Definition: ThermoPhase.h:227
virtual void setParameters(int n, doublereal *const c)
Class XML_Node is a tree-based representation of the contents of an XML file.
Definition: xml.h:100
void warn_deprecated(const std::string &method, const std::string &extra)
Print a warning indicating that method is deprecated.
Definition: global.cpp:76
virtual void setState_Psat(doublereal p, doublereal x)
Set the state to a saturated system at a particular pressure.
virtual void getEntropy_R_ref(doublereal *er) const
void setState_TPY(doublereal t, doublereal p, const doublereal *y)
Set the internally stored temperature (K), pressure (Pa), and mass fractions of the phase...
void getPartialMolarCp(doublereal *cpbar) const
Get the species partial molar Heat Capacities. Units: J/ kmol /K.
vector_fp m_s0_R
Dimensionless entropy at the (mtlast, m_p0)
Base class for a phase with thermodynamic properties.
Definition: ThermoPhase.h:101
void setState_PY(doublereal p, doublereal *y)
Set the internally stored pressure (Pa) and mass fractions.
void getChemPotentials_RT(doublereal *murt) const
Get the array of non-dimensional species chemical potentials These are partial molar Gibbs free energ...
void getPartialMolarEnthalpies(doublereal *hbar) const
Get the species partial molar enthalpies. Units: J/kmol.
virtual void getActivities(doublereal *a) const
Get the array of non-dimensional activities at the current solution temperature, pressure, and solution concentration.
virtual void getGibbs_RT_ref(doublereal *grt) const
virtual void getCp_R_ref(doublereal *cprt) const
void getPartialMolarEntropies(doublereal *sbar) const
Get the species partial molar entropy. Units: J/kmol K.
virtual int eosType() const
Returns the equation of state type flag.
doublereal enthalpy_mole() const
Molar enthalpy. Units: J/kmol.
void getChemPotentials(doublereal *mu) const
Get the array of chemical potentials.
SingleSpeciesTP & operator=(const SingleSpeciesTP &right)
Assignment operator.
void getElectrochemPotentials(doublereal *mu) const
Get the species electrochemical potentials. Units: J/kmol.
doublereal entropy_mole() const
Molar entropy. Units: J/kmol/K.
vector_fp m_h0_RT
Dimensionless enthalpy at the (mtlast, m_p0)
vector_fp m_cp0_R
Dimensionless heat capacity at the (mtlast, m_p0)
virtual void getGibbs_ref(doublereal *g) const
void getStandardVolumes(doublereal *vbar) const
Get the molar volumes of each species in their standard states at the current T and P of the solution...
virtual void setState_SP(doublereal s, doublereal p, doublereal tol=1.e-8)
Set the specific entropy (J/kg/K) and pressure (Pa).
doublereal gibbs_mole() const
Molar Gibbs function. Units: J/kmol.
doublereal intEnergy_mole() const
Molar internal energy. Units: J/kmol.
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 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).
void setState_TPX(doublereal t, doublereal p, const doublereal *x)
Set the temperature (K), pressure (Pa), and mole fractions.
virtual void getPartialMolarIntEnergies(doublereal *ubar) const
Get the species partial molar internal energies. Units: J/kmol.
doublereal m_press
The current pressure of the solution (Pa)
doublereal cp_mole() const
Molar heat capacity at constant pressure. Units: J/kmol/K.
doublereal m_tlast
Last temperature used to evaluate the thermodynamic polynomial.
void getPartialMolarVolumes(doublereal *vbar) const
Get the species partial molar volumes. Units: m^3/kmol.
size_t m_kk
Number of species in the phase.
Definition: Phase.h:716
void getPureGibbs(doublereal *gpure) const
Get the dimensional Gibbs functions for the standard state of the species at the current T and P...
virtual void setState_Tsat(doublereal t, doublereal x)
Set the state to a saturated system at a particular temperature.
SingleSpeciesTP()
Base empty constructor.
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 void getActivityCoefficients(doublereal *ac) const
Get the array of non-dimensional activity coefficients at the current solution temperature, pressure, and solution concentration.
SpeciesThermo * m_spthermo
Pointer to the calculation manager for species reference-state thermodynamic properties.
Definition: ThermoPhase.h:1625
void setState_PX(doublereal p, doublereal *x)
Set the pressure (Pa) and mole fractions.
The SingleSpeciesTP class is a filter class for ThermoPhase.
virtual doublereal vaporFraction() const
Return the fraction of vapor at the current conditions.
virtual void getParameters(int &n, doublereal *const c) const