Cantera  3.1.0a1
IdealSolidSolnPhase.h
Go to the documentation of this file.
1 /**
2  * @file IdealSolidSolnPhase.h Header file for an ideal solid solution model
3  * with incompressible thermodynamics (see @ref thermoprops and
4  * @link Cantera::IdealSolidSolnPhase IdealSolidSolnPhase@endlink).
5  *
6  * This class inherits from the %Cantera class ThermoPhase and implements an
7  * ideal solid solution model with incompressible thermodynamics.
8  */
9 
10 // This file is part of Cantera. See License.txt in the top-level directory or
11 // at https://cantera.org/license.txt for license and copyright information.
12 
13 #ifndef CT_IDEALSOLIDSOLNPHASE_H
14 #define CT_IDEALSOLIDSOLNPHASE_H
15 
16 #include "ThermoPhase.h"
17 
18 namespace Cantera
19 {
20 
21 /**
22  * Class IdealSolidSolnPhase represents a condensed phase ideal solution
23  * compound. The phase and the pure species phases which comprise the standard
24  * states of the species are assumed to have zero volume expansivity and zero
25  * isothermal compressibility. Each species does, however, have constant but
26  * distinct partial molar volumes equal to their pure species molar volumes. The
27  * class derives from class ThermoPhase, and overloads the virtual methods
28  * defined there with ones that use expressions appropriate for ideal solution
29  * mixtures.
30  *
31  * The generalized concentrations can have three different forms depending on
32  * the value of the member attribute #m_formGC, which is supplied in the
33  * constructor and in the input file. The value and form of the generalized
34  * concentration will affect reaction rate constants involving species in this
35  * phase.
36  *
37  * @ingroup thermoprops
38  */
40 {
41 public:
42  //! Construct and initialize an IdealSolidSolnPhase ThermoPhase object
43  //! directly from an input file
44  /*!
45  * This constructor will also fully initialize the object.
46  * The generalized concentrations can have three different forms
47  * depending on the value of the member attribute #m_formGC, which
48  * is supplied in the constructor or read from the input file.
49  *
50  * @param infile File name for the input file containing information
51  * for this phase. If blank, an empty phase will be
52  * created.
53  * @param id The name of this phase. This is used to look up
54  * the phase in the input file.
55  */
56  explicit IdealSolidSolnPhase(const string& infile="", const string& id="");
57 
58  string type() const override {
59  return "ideal-condensed";
60  }
61 
62  bool isIdeal() const override {
63  return true;
64  }
65 
66  bool isCompressible() const override {
67  return false;
68  }
69 
70  //! @name Molar Thermodynamic Properties of the Solution
71  //! @{
72 
73  /**
74  * Molar enthalpy of the solution. Units: J/kmol. For an ideal, constant
75  * partial molar volume solution mixture with pure species phases which
76  * exhibit zero volume expansivity and zero isothermal compressibility:
77  * @f[
78  * \hat h(T,P) = \sum_k X_k \hat h^0_k(T) + (P - P_{ref}) (\sum_k X_k \hat V^0_k)
79  * @f]
80  * The reference-state pure-species enthalpies at the reference pressure Pref
81  * @f$ \hat h^0_k(T) @f$, are computed by the species thermodynamic
82  * property manager. They are polynomial functions of temperature.
83  * @see MultiSpeciesThermo
84  */
85  double enthalpy_mole() const override;
86 
87  /**
88  * Molar entropy of the solution. Units: J/kmol/K. For an ideal, constant
89  * partial molar volume solution mixture with pure species phases which
90  * exhibit zero volume expansivity:
91  * @f[
92  * \hat s(T, P, X_k) = \sum_k X_k \hat s^0_k(T) - \hat R \sum_k X_k \ln(X_k)
93  * @f]
94  * The reference-state pure-species entropies
95  * @f$ \hat s^0_k(T,p_{ref}) @f$ are computed by the species thermodynamic
96  * property manager. The pure species entropies are independent of
97  * pressure since the volume expansivities are equal to zero.
98  * @see MultiSpeciesThermo
99  */
100  double entropy_mole() const override;
101 
102  /**
103  * Molar Gibbs free energy of the solution. Units: J/kmol. For an ideal,
104  * constant partial molar volume solution mixture with pure species phases
105  * which exhibit zero volume expansivity:
106  * @f[
107  * \hat g(T, P) = \sum_k X_k \hat g^0_k(T,P) + \hat R T \sum_k X_k \ln(X_k)
108  * @f]
109  * The reference-state pure-species Gibbs free energies
110  * @f$ \hat g^0_k(T) @f$ are computed by the species thermodynamic
111  * property manager, while the standard state Gibbs free energies
112  * @f$ \hat g^0_k(T,P) @f$ are computed by the member function, gibbs_RT().
113  * @see MultiSpeciesThermo
114  */
115  double gibbs_mole() const override;
116 
117  /**
118  * Molar heat capacity at constant pressure of the solution.
119  * Units: J/kmol/K.
120  * For an ideal, constant partial molar volume solution mixture with
121  * pure species phases which exhibit zero volume expansivity:
122  * @f[
123  * \hat c_p(T,P) = \sum_k X_k \hat c^0_{p,k}(T) .
124  * @f]
125  * The heat capacity is independent of pressure. The reference-state pure-
126  * species heat capacities @f$ \hat c^0_{p,k}(T) @f$ are computed by the
127  * species thermodynamic property manager.
128  * @see MultiSpeciesThermo
129  */
130  double cp_mole() const override;
131 
132  /**
133  * Molar heat capacity at constant volume of the solution. Units: J/kmol/K.
134  * For an ideal, constant partial molar volume solution mixture with pure
135  * species phases which exhibit zero volume expansivity:
136  * @f[ \hat c_v(T,P) = \hat c_p(T,P) @f]
137  * The two heat capacities are equal.
138  */
139  double cv_mole() const override {
140  return cp_mole();
141  }
142 
143  //! @}
144  //! @name Mechanical Equation of State Properties
145  //!
146  //! In this equation of state implementation, the density is a function only
147  //! of the mole fractions. Therefore, it can't be an independent variable.
148  //! Instead, the pressure is used as the independent variable. Functions
149  //! which try to set the thermodynamic state by calling setDensity() will
150  //! cause an exception to be thrown.
151  //! @{
152 
153  /**
154  * Pressure. Units: Pa. For this incompressible system, we return the
155  * internally stored independent value of the pressure.
156  */
157  double pressure() const override {
158  return m_Pcurrent;
159  }
160 
161  /**
162  * Set the pressure at constant temperature. Units: Pa. This method sets a
163  * constant within the object. The mass density is not a function of
164  * pressure.
165  *
166  * @param p Input Pressure (Pa)
167  */
168  void setPressure(double p) override;
169 
170  /**
171  * Calculate the density of the mixture using the partial molar volumes and
172  * mole fractions as input
173  *
174  * The formula for this is
175  *
176  * @f[
177  * \rho = \frac{\sum_k{X_k W_k}}{\sum_k{X_k V_k}}
178  * @f]
179  *
180  * where @f$ X_k @f$ are the mole fractions, @f$ W_k @f$ are the molecular
181  * weights, and @f$ V_k @f$ are the pure species molar volumes.
182  *
183  * Note, the basis behind this formula is that in an ideal solution the
184  * partial molar volumes are equal to the pure species molar volumes. We
185  * have additionally specified in this class that the pure species molar
186  * volumes are independent of temperature and pressure.
187  */
188  virtual void calcDensity();
189 
190  //! @}
191  //! @name Chemical Potentials and Activities
192  //!
193  //! The activity @f$ a_k @f$ of a species in solution is related to the
194  //! chemical potential by
195  //! @f[
196  //! \mu_k(T,P,X_k) = \mu_k^0(T,P) + \hat R T \ln a_k.
197  //! @f]
198  //! The quantity @f$ \mu_k^0(T,P) @f$ is the standard state chemical potential
199  //! at unit activity. It may depend on the pressure and the temperature.
200  //! However, it may not depend on the mole fractions of the species in the
201  //! solid solution.
202  //!
203  //! The activities are related to the generalized concentrations, @f$ \tilde
204  //! C_k @f$, and standard concentrations, @f$ C^0_k @f$, by the following
205  //! formula:
206  //!
207  //! @f[
208  //! a_k = \frac{\tilde C_k}{C^0_k}
209  //! @f]
210  //! The generalized concentrations are used in the kinetics classes to
211  //! describe the rates of progress of reactions involving the species. Their
212  //! formulation depends upon the specification of the rate constants for
213  //! reaction, especially the units used in specifying the rate constants. The
214  //! bridge between the thermodynamic equilibrium expressions that use a_k and
215  //! the kinetics expressions which use the generalized concentrations is
216  //! provided by the multiplicative factor of the standard concentrations.
217  //! @{
218 
219  Units standardConcentrationUnits() const override;
220 
221  /**
222  * This method returns the array of generalized concentrations. The
223  * generalized concentrations are used in the evaluation of the rates of
224  * progress for reactions involving species in this phase. The generalized
225  * concentration divided by the standard concentration is also equal to the
226  * activity of species.
227  *
228  * For this implementation the activity is defined to be the mole fraction
229  * of the species. The generalized concentration is defined to be equal to
230  * the mole fraction divided by the partial molar volume. The generalized
231  * concentrations for species in this phase therefore have units of
232  * kmol/m^3. Rate constants must reflect this fact.
233  *
234  * On a general note, the following must be true. For an ideal solution, the
235  * generalized concentration must consist of the mole fraction multiplied by
236  * a constant. The constant may be fairly arbitrarily chosen, with
237  * differences adsorbed into the reaction rate expression. 1/V_N, 1/V_k, or
238  * 1 are equally good, as long as the standard concentration is adjusted
239  * accordingly. However, it must be a constant (and not the concentration,
240  * btw, which is a function of the mole fractions) in order for the ideal
241  * solution properties to hold at the same time having the standard
242  * concentration to be independent of the mole fractions.
243  *
244  * In this implementation the form of the generalized concentrations
245  * depend upon the member attribute, #m_formGC.
246  *
247  * HKM Note: We have absorbed the pressure dependence of the pure species
248  * state into the thermodynamics functions. Therefore the standard
249  * state on which the activities are based depend on both temperature
250  * and pressure. If we hadn't, it would have appeared in this
251  * function in a very awkward exp[] format.
252  *
253  * @param c Pointer to array of doubles of length m_kk, which on exit
254  * will contain the generalized concentrations.
255  */
256  void getActivityConcentrations(double* c) const override;
257 
258  /**
259  * The standard concentration @f$ C^0_k @f$ used to normalize the
260  * generalized concentration. In many cases, this quantity will be the
261  * same for all species in a phase. However, for this case, we will return
262  * a distinct concentration for each species. This is the inverse of the
263  * species molar volume. Units for the standard concentration are kmol/m^3.
264  *
265  * @param k Species number: this is a require parameter, a change from the
266  * ThermoPhase base class, where it was an optional parameter.
267  */
268  double standardConcentration(size_t k) const override;
269 
270  //! Get the array of species activity coefficients
271  /*!
272  * @param ac output vector of activity coefficients. Length: m_kk
273  */
274  void getActivityCoefficients(double* ac) const override;
275 
276  /**
277  * Get the species chemical potentials. Units: J/kmol.
278  *
279  * This function returns a vector of chemical potentials of the
280  * species in solution.
281  * @f[
282  * \mu_k = \mu^{ref}_k(T) + V_k * (p - p_o) + R T \ln(X_k)
283  * @f]
284  * or another way to phrase this is
285  * @f[
286  * \mu_k = \mu^o_k(T,p) + R T \ln(X_k)
287  * @f]
288  * where @f$ \mu^o_k(T,p) = \mu^{ref}_k(T) + V_k * (p - p_o) @f$
289  *
290  * @param mu Output vector of chemical potentials.
291  */
292  void getChemPotentials(double* mu) const override;
293 
294  //! @}
295  //! @name Partial Molar Properties of the Solution
296  //! @{
297 
298  //! Returns an array of partial molar enthalpies for the species in the
299  //! mixture.
300  /*!
301  * Units (J/kmol). For this phase, the partial molar enthalpies are equal to
302  * the pure species enthalpies
303  * @f[
304  * \bar h_k(T,P) = \hat h^{ref}_k(T) + (P - P_{ref}) \hat V^0_k
305  * @f]
306  * The reference-state pure-species enthalpies, @f$ \hat h^{ref}_k(T) @f$,
307  * at the reference pressure,@f$ P_{ref} @f$, are computed by the species
308  * thermodynamic property manager. They are polynomial functions of
309  * temperature.
310  * @see MultiSpeciesThermo
311  *
312  * @param hbar Output vector containing partial molar enthalpies.
313  * Length: m_kk.
314  */
315  void getPartialMolarEnthalpies(double* hbar) const override;
316 
317  /**
318  * Returns an array of partial molar entropies of the species in the
319  * solution. Units: J/kmol/K. For this phase, the partial molar entropies
320  * are equal to the pure species entropies plus the ideal solution
321  * contribution.
322  * @f[
323  * \bar s_k(T,P) = \hat s^0_k(T) - R \ln(X_k)
324  * @f]
325  * The reference-state pure-species entropies,@f$ \hat s^{ref}_k(T) @f$, at
326  * the reference pressure, @f$ P_{ref} @f$, are computed by the species
327  * thermodynamic property manager. They are polynomial functions of
328  * temperature.
329  * @see MultiSpeciesThermo
330  *
331  * @param sbar Output vector containing partial molar entropies.
332  * Length: m_kk.
333  */
334  void getPartialMolarEntropies(double* sbar) const override;
335 
336  /**
337  * Returns an array of partial molar Heat Capacities at constant pressure of
338  * the species in the solution. Units: J/kmol/K. For this phase, the partial
339  * molar heat capacities are equal to the standard state heat capacities.
340  *
341  * @param cpbar Output vector of partial heat capacities. Length: m_kk.
342  */
343  void getPartialMolarCp(double* cpbar) const override;
344 
345  /**
346  * returns an array of partial molar volumes of the species
347  * in the solution. Units: m^3 kmol-1.
348  *
349  * For this solution, the partial molar volumes are equal to the
350  * constant species molar volumes.
351  *
352  * @param vbar Output vector of partial molar volumes. Length: m_kk.
353  */
354  void getPartialMolarVolumes(double* vbar) const override;
355 
356  //! @}
357  //! @name Properties of the Standard State of the Species in the Solution
358  //! @{
359 
360  /**
361  * Get the standard state chemical potentials of the species. This is the
362  * array of chemical potentials at unit activity @f$ \mu^0_k(T,P) @f$. We
363  * define these here as the chemical potentials of the pure species at the
364  * temperature and pressure of the solution. This function is used in the
365  * evaluation of the equilibrium constant Kc. Therefore, Kc will also depend
366  * on T and P. This is the norm for liquid and solid systems.
367  *
368  * units = J / kmol
369  *
370  * @param mu0 Output vector of standard state chemical potentials.
371  * Length: m_kk.
372  */
373  void getStandardChemPotentials(double* mu0) const override {
374  getPureGibbs(mu0);
375  }
376 
377  //! Get the array of nondimensional Enthalpy functions for the standard
378  //! state species at the current *T* and *P* of the solution.
379  /*!
380  * We assume an incompressible constant partial molar volume here:
381  * @f[
382  * h^0_k(T,P) = h^{ref}_k(T) + (P - P_{ref}) * V_k
383  * @f]
384  * where @f$ V_k @f$ is the molar volume of pure species *k*.
385  * @f$ h^{ref}_k(T) @f$ is the enthalpy of the pure species *k* at the
386  * reference pressure, @f$ P_{ref} @f$.
387  *
388  * @param hrt Vector of length m_kk, which on return hrt[k] will contain the
389  * nondimensional standard state enthalpy of species k.
390  */
391  void getEnthalpy_RT(double* hrt) const override;
392 
393  //! Get the nondimensional Entropies for the species standard states at the
394  //! current T and P of the solution.
395  /*!
396  * Note, this is equal to the reference state entropies due to the zero
397  * volume expansivity: that is, (dS/dP)_T = (dV/dT)_P = 0.0
398  *
399  * @param sr Vector of length m_kk, which on return sr[k] will contain the
400  * nondimensional standard state entropy for species k.
401  */
402  void getEntropy_R(double* sr) const override;
403 
404  /**
405  * Get the nondimensional Gibbs function for the species standard states at
406  * the current T and P of the solution.
407  *
408  * @f[
409  * \mu^0_k(T,P) = \mu^{ref}_k(T) + (P - P_{ref}) * V_k
410  * @f]
411  * where @f$ V_k @f$ is the molar volume of pure species *k*.
412  * @f$ \mu^{ref}_k(T) @f$ is the chemical potential of pure species *k*
413  * at the reference pressure, @f$ P_{ref} @f$.
414  *
415  * @param grt Vector of length m_kk, which on return sr[k] will contain the
416  * nondimensional standard state Gibbs function for species k.
417  */
418  void getGibbs_RT(double* grt) const override;
419 
420  /**
421  * Get the Gibbs functions for the pure species at the current *T* and *P*
422  * of the solution. We assume an incompressible constant partial molar
423  * volume here:
424  * @f[
425  * \mu^0_k(T,P) = \mu^{ref}_k(T) + (P - P_{ref}) * V_k
426  * @f]
427  * where @f$ V_k @f$ is the molar volume of pure species *k*.
428  * @f$ \mu^{ref}_k(T) @f$ is the chemical potential of pure species *k* at
429  * the reference pressure, @f$ P_{ref} @f$.
430  *
431  * @param gpure Output vector of Gibbs functions for species. Length: m_kk.
432  */
433  void getPureGibbs(double* gpure) const override;
434 
435  void getIntEnergy_RT(double* urt) const override;
436 
437  /**
438  * Get the nondimensional heat capacity at constant pressure function for
439  * the species standard states at the current T and P of the solution.
440  * @f[
441  * Cp^0_k(T,P) = Cp^{ref}_k(T)
442  * @f]
443  * where @f$ V_k @f$ is the molar volume of pure species *k*.
444  * @f$ Cp^{ref}_k(T) @f$ is the constant pressure heat capacity of species
445  * *k* at the reference pressure, @f$ p_{ref} @f$.
446  *
447  * @param cpr Vector of length m_kk, which on return cpr[k] will contain the
448  * nondimensional constant pressure heat capacity for species k.
449  */
450  void getCp_R(double* cpr) const override;
451 
452  void getStandardVolumes(double* vol) const override;
453 
454  //! @}
455  //! @name Thermodynamic Values for the Species Reference States
456  //! @{
457 
458  void getEnthalpy_RT_ref(double* hrt) const override;
459  void getGibbs_RT_ref(double* grt) const override;
460  void getGibbs_ref(double* g) const override;
461  void getEntropy_R_ref(double* er) const override;
462  void getIntEnergy_RT_ref(double* urt) const override;
463  void getCp_R_ref(double* cprt) const override;
464 
465  /**
466  * Returns a reference to the vector of nondimensional enthalpies of the
467  * reference state at the current temperature. Real reason for its existence
468  * is that it also checks to see if a recalculation of the reference
469  * thermodynamics functions needs to be done.
470  */
471  const vector<double>& enthalpy_RT_ref() const;
472 
473  /**
474  * Returns a reference to the vector of nondimensional enthalpies of the
475  * reference state at the current temperature. Real reason for its existence
476  * is that it also checks to see if a recalculation of the reference
477  * thermodynamics functions needs to be done.
478  */
479  const vector<double>& gibbs_RT_ref() const {
480  _updateThermo();
481  return m_g0_RT;
482  }
483 
484  /**
485  * Returns a reference to the vector of nondimensional enthalpies of the
486  * reference state at the current temperature. Real reason for its existence
487  * is that it also checks to see if a recalculation of the reference
488  * thermodynamics functions needs to be done.
489  */
490  const vector<double>& entropy_R_ref() const;
491 
492  /**
493  * Returns a reference to the vector of nondimensional enthalpies of the
494  * reference state at the current temperature. Real reason for its existence
495  * is that it also checks to see if a recalculation of the reference
496  * thermodynamics functions needs to be done.
497  */
498  const vector<double>& cp_R_ref() const {
499  _updateThermo();
500  return m_cp0_R;
501  }
502 
503  //! @}
504  //! @name Utility Functions
505  //! @{
506 
507  bool addSpecies(shared_ptr<Species> spec) override;
508  void initThermo() override;
509  void getParameters(AnyMap& phaseNode) const override;
510  void getSpeciesParameters(const string& name, AnyMap& speciesNode) const override;
511  void setToEquilState(const double* mu_RT) override;
512 
513  //! Set the form for the standard and generalized concentrations
514  /*!
515  * Must be one of 'unity', 'species-molar-volume', or
516  * 'solvent-molar-volume'. The default is 'unity'.
517  *
518  * | m_formGC | GeneralizedConc | StandardConc |
519  * | -------------------- | --------------- | ------------ |
520  * | unity | X_k | 1.0 |
521  * | species-molar-volume | X_k / V_k | 1.0 / V_k |
522  * | solvent-molar-volume | X_k / V_N | 1.0 / V_N |
523  *
524  * The value and form of the generalized concentration will affect
525  * reaction rate constants involving species in this phase.
526  */
527  void setStandardConcentrationModel(const string& model);
528 
529  /**
530  * Report the molar volume of species k
531  *
532  * units - @f$ m^3 kmol^-1 @f$
533  *
534  * @param k species index
535  */
536  double speciesMolarVolume(int k) const;
537 
538  /**
539  * Fill in a return vector containing the species molar volumes.
540  *
541  * units - @f$ m^3 kmol^-1 @f$
542  *
543  * @param smv output vector containing species molar volumes.
544  * Length: m_kk.
545  */
546  void getSpeciesMolarVolumes(double* smv) const;
547 
548  //! @}
549 
550 protected:
551  void compositionChanged() override;
552 
553  /**
554  * The standard concentrations can have one of three different forms:
555  * 0 = 'unity', 1 = 'species-molar-volume', 2 = 'solvent-molar-volume'. See
556  * setStandardConcentrationModel().
557  */
558  int m_formGC = 0;
559 
560  /**
561  * Value of the reference pressure for all species in this phase. The T
562  * dependent polynomials are evaluated at the reference pressure. Note,
563  * because this is a single value, all species are required to have the same
564  * reference pressure.
565  */
566  double m_Pref = OneAtm;
567 
568  /**
569  * m_Pcurrent = The current pressure
570  * Since the density isn't a function of pressure, but only of the
571  * mole fractions, we need to independently specify the pressure.
572  * The density variable which is inherited as part of the State class,
573  * m_dens, is always kept current whenever T, P, or X[] change.
574  */
575  double m_Pcurrent = OneAtm;
576 
577  //! Vector of molar volumes for each species in the solution
578  /**
579  * Species molar volumes (@f$ m^3 kmol^-1 @f$) at the current mixture state.
580  * For the IdealSolidSolnPhase class, these are constant.
581  */
582  mutable vector<double> m_speciesMolarVolume;
583 
584  //! Vector containing the species reference enthalpies at T = m_tlast
585  mutable vector<double> m_h0_RT;
586 
587  //! Vector containing the species reference constant pressure heat
588  //! capacities at T = m_tlast
589  mutable vector<double> m_cp0_R;
590 
591  //! Vector containing the species reference Gibbs functions at T = m_tlast
592  mutable vector<double> m_g0_RT;
593 
594  //! Vector containing the species reference entropies at T = m_tlast
595  mutable vector<double> m_s0_R;
596 
597  //! Vector containing the species reference exp(-G/RT) functions at
598  //! T = m_tlast
599  mutable vector<double> m_expg0_RT;
600 
601  //! Temporary array used in equilibrium calculations
602  mutable vector<double> m_pp;
603 
604 private:
605  //! @name Utility Functions
606  //! @{
607  /**
608  * This function gets called for every call to functions in this class. It
609  * checks to see whether the temperature has changed and thus the reference
610  * thermodynamics functions for all of the species must be recalculated. If
611  * the temperature has changed, the species thermo manager is called to
612  * recalculate G, Cp, H, and S at the current temperature.
613  */
614  virtual void _updateThermo() const;
615 
616  //! @}
617 };
618 }
619 
620 #endif
Header file for class ThermoPhase, the base class for phases with thermodynamic properties,...
A map of string keys to values whose type can vary at runtime.
Definition: AnyMap.h:427
Class IdealSolidSolnPhase represents a condensed phase ideal solution compound.
const vector< double > & entropy_R_ref() const
Returns a reference to the vector of nondimensional enthalpies of the reference state at the current ...
void getStandardChemPotentials(double *mu0) const override
Get the standard state chemical potentials of the species.
double enthalpy_mole() const override
Molar enthalpy of the solution.
bool isIdeal() const override
Boolean indicating whether phase is ideal.
void getPartialMolarEnthalpies(double *hbar) const override
Returns an array of partial molar enthalpies for the species in the mixture.
void getChemPotentials(double *mu) const override
Get the species chemical potentials.
double pressure() const override
Pressure.
vector< double > m_g0_RT
Vector containing the species reference Gibbs functions at T = m_tlast.
bool isCompressible() const override
Return whether phase represents a compressible substance.
const vector< double > & cp_R_ref() const
Returns a reference to the vector of nondimensional enthalpies of the reference state at the current ...
void getSpeciesParameters(const string &name, AnyMap &speciesNode) const override
Get phase-specific parameters of a Species object such that an identical one could be reconstructed a...
const vector< double > & gibbs_RT_ref() const
Returns a reference to the vector of nondimensional enthalpies of the reference state at the current ...
void getEntropy_R(double *sr) const override
Get the nondimensional Entropies for the species standard states at the current T and P of the soluti...
vector< double > m_h0_RT
Vector containing the species reference enthalpies at T = m_tlast.
void getGibbs_ref(double *g) const override
Returns the vector of the Gibbs function of the reference state at the current temperature of the sol...
double speciesMolarVolume(int k) const
Report the molar volume of species k.
vector< double > m_pp
Temporary array used in equilibrium calculations.
double m_Pref
Value of the reference pressure for all species in this phase.
void getCp_R(double *cpr) const override
Get the nondimensional heat capacity at constant pressure function for the species standard states at...
void getParameters(AnyMap &phaseNode) const override
Store the parameters of a ThermoPhase object such that an identical one could be reconstructed using ...
string type() const override
String indicating the thermodynamic model implemented.
void initThermo() override
Initialize the ThermoPhase object after all species have been set up.
void getActivityConcentrations(double *c) const override
This method returns the array of generalized concentrations.
void getSpeciesMolarVolumes(double *smv) const
Fill in a return vector containing the species molar volumes.
void setPressure(double p) override
Set the pressure at constant temperature.
void getPartialMolarVolumes(double *vbar) const override
returns an array of partial molar volumes of the species in the solution.
double cv_mole() const override
Molar heat capacity at constant volume of the solution.
void getPureGibbs(double *gpure) const override
Get the Gibbs functions for the pure species at the current T and P of the solution.
double standardConcentration(size_t k) const override
The standard concentration used to normalize the generalized concentration.
void setStandardConcentrationModel(const string &model)
Set the form for the standard and generalized concentrations.
void getIntEnergy_RT_ref(double *urt) const override
Returns the vector of nondimensional internal Energies of the reference state at the current temperat...
void getEnthalpy_RT(double *hrt) const override
Get the array of nondimensional Enthalpy functions for the standard state species at the current T an...
void getEntropy_R_ref(double *er) const override
Returns the vector of nondimensional entropies of the reference state at the current temperature of t...
vector< double > m_s0_R
Vector containing the species reference entropies at T = m_tlast.
void getGibbs_RT(double *grt) const override
Get the nondimensional Gibbs function for the species standard states at the current T and P of the s...
double entropy_mole() const override
Molar entropy of the solution.
int m_formGC
The standard concentrations can have one of three different forms: 0 = 'unity', 1 = 'species-molar-vo...
vector< double > m_speciesMolarVolume
Vector of molar volumes for each species in the solution.
void getCp_R_ref(double *cprt) const override
Returns the vector of nondimensional constant pressure heat capacities of the reference state at the ...
void getStandardVolumes(double *vol) const override
Get the molar volumes of the species standard states at the current T and P of the solution.
double cp_mole() const override
Molar heat capacity at constant pressure of the solution.
void getIntEnergy_RT(double *urt) const override
Returns the vector of nondimensional Internal Energies of the standard state species at the current T...
Units standardConcentrationUnits() const override
Returns the units of the "standard concentration" for this phase.
IdealSolidSolnPhase(const string &infile="", const string &id="")
Construct and initialize an IdealSolidSolnPhase ThermoPhase object directly from an input file.
void getPartialMolarCp(double *cpbar) const override
Returns an array of partial molar Heat Capacities at constant pressure of the species in the solution...
void compositionChanged() override
Apply changes to the state which are needed after the composition changes.
double gibbs_mole() const override
Molar Gibbs free energy of the solution.
vector< double > m_cp0_R
Vector containing the species reference constant pressure heat capacities at T = m_tlast.
bool addSpecies(shared_ptr< Species > spec) override
Add a Species to this Phase.
virtual void _updateThermo() const
This function gets called for every call to functions in this class.
void setToEquilState(const double *mu_RT) override
This method is used by the ChemEquil equilibrium solver.
void getGibbs_RT_ref(double *grt) const override
Returns the vector of nondimensional Gibbs Free Energies of the reference state at the current temper...
void getActivityCoefficients(double *ac) const override
Get the array of species activity coefficients.
vector< double > m_expg0_RT
Vector containing the species reference exp(-G/RT) functions at T = m_tlast.
void getPartialMolarEntropies(double *sbar) const override
Returns an array of partial molar entropies of the species in the solution.
double m_Pcurrent
m_Pcurrent = The current pressure Since the density isn't a function of pressure, but only of the mol...
void getEnthalpy_RT_ref(double *hrt) const override
Returns the vector of nondimensional enthalpies of the reference state at the current temperature of ...
virtual void calcDensity()
Calculate the density of the mixture using the partial molar volumes and mole fractions as input.
const vector< double > & enthalpy_RT_ref() const
Returns a reference to the vector of nondimensional enthalpies of the reference state at the current ...
string name() const
Return the name of the phase.
Definition: Phase.cpp:20
Base class for a phase with thermodynamic properties.
Definition: ThermoPhase.h:390
A representation of the units associated with a dimensional quantity.
Definition: Units.h:35
const double OneAtm
One atmosphere [Pa].
Definition: ct_defs.h:96
Namespace for the Cantera kernel.
Definition: AnyMap.cpp:564