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