Cantera  2.0
RedlichKisterVPSSTP.h
Go to the documentation of this file.
1 /**
2  * @file RedlichKisterVPSSTP.h
3  * Header for intermediate ThermoPhase object for phases which
4  * employ gibbs excess free energy based formulations
5  * (see \ref thermoprops
6  * and class \link Cantera::RedlichKisterVPSSTP RedlichKisterVPSSTP\endlink).
7  *
8  * Header file for a derived class of ThermoPhase that handles
9  * variable pressure standard state methods for calculating
10  * thermodynamic properties that are further based upon activities
11  * based on the molality scale. These include most of the methods for
12  * calculating liquid electrolyte thermodynamics.
13  */
14 /*
15  * Copyright (2006) Sandia Corporation. Under the terms of
16  * Contract DE-AC04-94AL85000 with Sandia Corporation, the
17  * U.S. Government retains certain rights in this software.
18  */
19 
20 #ifndef CT_REDLICHKISTERVPSSTP_H
21 #define CT_REDLICHKISTERVPSSTP_H
22 
24 #include "cantera/base/Array.h"
25 
26 namespace Cantera
27 {
28 
29 /**
30  * @ingroup thermoprops
31  */
32 
33 //! RedlichKisterVPSSTP is a derived class of GibbsExcessVPSSTP that employs
34 //! the Redlich-Kister approximation for the excess gibbs free energy
35 /*!
36  *
37  * %RedlichKisterVPSSTP derives from class GibbsExcessVPSSTP which is derived
38  * from VPStandardStateTP, and overloads the virtual methods defined there with ones that
39  * use expressions appropriate for the Redlich Kister Excess gibbs free energy approximation.
40  *
41  * The independent unknowns are pressure, temperature, and mass fraction.
42  *
43  * Several concepts are introduced. The first concept is there are temporary
44  * variables for holding the species standard state values of Cp, H, S, G, and V at the
45  * last temperature and pressure called. These functions are not recalculated
46  * if a new call is made using the previous temperature and pressure. Currently,
47  * these variables and the calculation method are handled by the VPSSMgr class,
48  * for which VPStandardStateTP owns a pointer to.
49  *
50  * To support the above functionality, pressure and temperature variables,
51  * m_plast_ss and m_tlast_ss, are kept which store the last pressure and temperature
52  * used in the evaluation of standard state properties.
53  *
54  * This class is usually used for nearly incompressible phases. For those phases, it
55  * makes sense to change the equation of state independent variable from
56  * density to pressure. The variable m_Pcurrent contains the current value of the
57  * pressure within the phase.
58  *
59  *
60  * <HR>
61  * <H2> Specification of Species Standard %State Properties </H2>
62  * <HR>
63  *
64  * All species are defined to have standard states that depend upon both
65  * the temperature and the pressure. The Redlich-Kister approximation assumes
66  * symmetric standard states, where all of the standard state assume
67  * that the species are in pure component states at the temperature
68  * and pressure of the solution. I don't think it prevents, however,
69  * some species from being dilute in the solution.
70  *
71  *
72  * <HR>
73  * <H2> Specification of Solution Thermodynamic Properties </H2>
74  * <HR>
75  *
76  * The molar excess Gibbs free energy is given by the following formula which is a sum over interactions <I>i</I>.
77  * Each of the interactions are binary interactions involving two of the species in the phase, denoted, <I>Ai</I>
78  * and <I>Bi</I>.
79  * This is the generalization of the Redlich-Kister formulation for a phase that has more than 2 species.
80  *
81  * \f[
82  * G^E = \sum_{i} G^E_{i}
83  * \f]
84  *
85  * where
86  *
87  * \f[
88  * G^E_{i} = n X_{Ai} X_{Bi} \sum_m \left( A^{i}_m {\left( X_{Ai} - X_{Bi} \right)}^m \right)
89  * \f]
90  *
91  * and where we can break down the gibbs free energy contributions into enthalpy and entropy contributions
92  *
93  * \f[
94  * H^E_i = n X_{Ai} X_{Bi} \sum_m \left( H^{i}_m {\left( X_{Ai} - X_{Bi} \right)}^m \right)
95  * \f]
96  *
97  * \f[
98  * S^E_i = n X_{Ai} X_{Bi} \sum_m \left( S^{i}_m {\left( X_{Ai} - X_{Bi} \right)}^m \right)
99  * \f]
100  *
101  * where n is the total moles in the solution.
102  *
103  * The activity of a species defined in the phase is given by an excess Gibbs free energy formulation.
104  *
105  * \f[
106  * a_k = \gamma_k X_k
107  * \f]
108  *
109  * where
110  *
111  * \f[
112  * R T \ln( \gamma_k )= \frac{d(n G^E)}{d(n_k)}\Bigg|_{n_i}
113  * \f]
114  *
115  * Taking the derivatives results in the following expression
116  * \f[
117  * R T \ln( \gamma_k )= \sum_i \delta_{Ai,k} (1 - X_{Ai}) X_{Bi} \sum_m \left( A^{i}_m {\left( X_{Ai} - X_{Bi} \right)}^m \right)
118  * + \sum_i \delta_{Ai,k} X_{Ai} X_{Bi} \sum_m \left( A^{i}_0 + A^{i}_m {\left( X_{Ai} - X_{Bi} \right)}^{m-1} (1 - X_{Ai} + X_{Bi}) \right)
119  * \f]
120  * where
121  *
122  *
123  * This object inherits from the class VPStandardStateTP. Therefore, the specification and
124  * calculation of all standard state and reference state values are handled at that level. Various functional
125  * forms for the standard state are permissible.
126  * The chemical potential for species <I>k</I> is equal to
127  *
128  * \f[
129  * \mu_k(T,P) = \mu^o_k(T, P) + R T \ln(\gamma_k X_k)
130  * \f]
131  *
132  * The partial molar entropy for species <I>k</I> is given by the following relation,
133  *
134  * \f[
135  * \tilde{s}_k(T,P) = s^o_k(T,P) - R \ln( \gamma_k X_k )
136  * - R T \frac{d \ln(\gamma_k) }{dT}
137  * \f]
138  *
139  * The partial molar enthalpy for species <I>k</I> is given by
140  *
141  * \f[
142  * \tilde{h}_k(T,P) = h^o_k(T,P) - R T^2 \frac{d \ln(\gamma_k)}{dT}
143  * \f]
144  *
145  * The partial molar volume for species <I>k</I> is
146  *
147  * \f[
148  * \tilde V_k(T,P) = V^o_k(T,P) + R T \frac{d \ln(\gamma_k) }{dP}
149  * \f]
150  *
151  * The partial molar Heat Capacity for species <I>k</I> is
152  *
153  * \f[
154  * \tilde{C}_{p,k}(T,P) = C^o_{p,k}(T,P) - 2 R T \frac{d \ln( \gamma_k )}{dT}
155  * - R T^2 \frac{d^2 \ln(\gamma_k) }{{dT}^2}
156  * \f]
157  *
158  * <HR>
159  * <H2> %Application within %Kinetics Managers </H2>
160  * <HR>
161  *
162  * \f$ C^a_k\f$ are defined such that \f$ a_k = C^a_k /
163  * C^s_k, \f$ where \f$ C^s_k \f$ is a standard concentration
164  * defined below and \f$ a_k \f$ are activities used in the
165  * thermodynamic functions. These activity (or generalized)
166  * concentrations are used
167  * by kinetics manager classes to compute the forward and
168  * reverse rates of elementary reactions.
169  * The activity concentration,\f$ C^a_k \f$,is given by the following expression.
170  *
171  * \f[
172  * C^a_k = C^s_k X_k = \frac{P}{R T} X_k
173  * \f]
174  *
175  * The standard concentration for species <I>k</I> is independent of <I>k</I> and equal to
176  *
177  * \f[
178  * C^s_k = C^s = \frac{P}{R T}
179  * \f]
180  *
181  * For example, a bulk-phase binary gas reaction between species j and k, producing
182  * a new gas species l would have the
183  * following equation for its rate of progress variable, \f$ R^1 \f$, which has
184  * units of kmol m-3 s-1.
185  *
186  * \f[
187  * R^1 = k^1 C_j^a C_k^a = k^1 (C^s a_j) (C^s a_k)
188  * \f]
189  * where
190  * \f[
191  * C_j^a = C^s a_j \mbox{\quad and \quad} C_k^a = C^s a_k
192  * \f]
193  *
194  *
195  * \f$ C_j^a \f$ is the activity concentration of species j, and
196  * \f$ C_k^a \f$ is the activity concentration of species k. \f$ C^s \f$
197  * is the standard concentration. \f$ a_j \f$ is
198  * the activity of species j which is equal to the mole fraction of j.
199  *
200  * The reverse rate constant can then be obtained from the law of microscopic reversibility
201  * and the equilibrium expression for the system.
202  *
203  * \f[
204  * \frac{a_j a_k}{ a_l} = K_a^{o,1} = \exp(\frac{\mu^o_l - \mu^o_j - \mu^o_k}{R T} )
205  * \f]
206  *
207  * \f$ K_a^{o,1} \f$ is the dimensionless form of the equilibrium constant, associated with
208  * the pressure dependent standard states \f$ \mu^o_l(T,P) \f$ and their associated activities,
209  * \f$ a_l \f$, repeated here:
210  *
211  * \f[
212  * \mu_l(T,P) = \mu^o_l(T, P) + R T \log(a_l)
213  * \f]
214  *
215  * We can switch over to expressing the equilibrium constant in terms of the reference
216  * state chemical potentials
217  *
218  * \f[
219  * K_a^{o,1} = \exp(\frac{\mu^{ref}_l - \mu^{ref}_j - \mu^{ref}_k}{R T} ) * \frac{P_{ref}}{P}
220  * \f]
221  *
222  * The concentration equilibrium constant, \f$ K_c \f$, may be obtained by changing over
223  * to activity concentrations. When this is done:
224  *
225  * \f[
226  * \frac{C^a_j C^a_k}{ C^a_l} = C^o K_a^{o,1} = K_c^1 =
227  * \exp(\frac{\mu^{ref}_l - \mu^{ref}_j - \mu^{ref}_k}{R T} ) * \frac{P_{ref}}{RT}
228  * \f]
229  *
230  * %Kinetics managers will calculate the concentration equilibrium constant, \f$ K_c \f$,
231  * using the second and third part of the above expression as a definition for the concentration
232  * equilibrium constant.
233  *
234  * For completeness, the pressure equilibrium constant may be obtained as well
235  *
236  * \f[
237  * \frac{P_j P_k}{ P_l P_{ref}} = K_p^1 = \exp(\frac{\mu^{ref}_l - \mu^{ref}_j - \mu^{ref}_k}{R T} )
238  * \f]
239  *
240  * \f$ K_p \f$ is the simplest form of the equilibrium constant for ideal gases. However, it isn't
241  * necessarily the simplest form of the equilibrium constant for other types of phases; \f$ K_c \f$ is
242  * used instead because it is completely general.
243  *
244  * The reverse rate of progress may be written down as
245  * \f[
246  * R^{-1} = k^{-1} C_l^a = k^{-1} (C^o a_l)
247  * \f]
248  *
249  * where we can use the concept of microscopic reversibility to
250  * write the reverse rate constant in terms of the
251  * forward reate constant and the concentration equilibrium
252  * constant, \f$ K_c \f$.
253  *
254  * \f[
255  * k^{-1} = k^1 K^1_c
256  * \f]
257  *
258  * \f$k^{-1} \f$ has units of s-1.
259  *
260  *
261  * <HR>
262  * <H2> Instantiation of the Class </H2>
263  * <HR>
264  *
265  *
266  * The constructor for this phase is located in the default ThermoFactory
267  * for %Cantera. A new %IdealGasPhase may be created by the following code
268  * snippet:
269  *
270  * @code
271  * XML_Node *xc = get_XML_File("silane.xml");
272  * XML_Node * const xs = xc->findNameID("phase", "silane");
273  * ThermoPhase *silane_tp = newPhase(*xs);
274  * IdealGasPhase *silaneGas = dynamic_cast <IdealGasPhase *>(silane_tp);
275  * @endcode
276  *
277  * or by the following constructor:
278  *
279  * @code
280  * XML_Node *xc = get_XML_File("silane.xml");
281  * XML_Node * const xs = xc->findNameID("phase", "silane");
282  * IdealGasPhase *silaneGas = new IdealGasPhase(*xs);
283  * @endcode
284  *
285  * <HR>
286  * <H2> XML Example </H2>
287  * <HR>
288  * An example of an XML Element named phase setting up a IdealGasPhase
289  * object named silane is given below.
290  *
291  *
292  * @verbatim
293  <!-- phase silane -->
294  <phase dim="3" id="silane">
295  <elementArray datasrc="elements.xml"> Si H He </elementArray>
296  <speciesArray datasrc="#species_data">
297  H2 H HE SIH4 SI SIH SIH2 SIH3 H3SISIH SI2H6
298  H2SISIH2 SI3H8 SI2 SI3
299  </speciesArray>
300  <reactionArray datasrc="#reaction_data"/>
301  <thermo model="IdealGas"/>
302  <kinetics model="GasKinetics"/>
303  <transport model="None"/>
304  </phase>
305  @endverbatim
306  *
307  * The model attribute "IdealGas" of the thermo XML element identifies the phase as
308  * being of the type handled by the IdealGasPhase object.
309  *
310  * @ingroup thermoprops
311  *
312 
313 */
315 {
316 
317 public:
318 
319  //! Constructor
320  /*!
321  * This doesn't do much more than initialize constants with
322  * default values.
323  */
325 
326  //! Construct and initialize a RedlichKisterVPSSTP ThermoPhase object
327  //! directly from an xml input file
328  /*!
329  * Working constructors
330  *
331  * The two constructors below are the normal way the phase initializes itself. They are shells that call
332  * the routine initThermo(), with a reference to the XML database to get the info for the phase.
333  *
334  * @param inputFile Name of the input file containing the phase XML data
335  * to set up the object
336  * @param id ID of the phase in the input file. Defaults to the
337  * empty string.
338  */
339  RedlichKisterVPSSTP(std::string inputFile, std::string id = "");
340 
341  //! Construct and initialize a RedlichKisterVPSSTP ThermoPhase object
342  //! directly from an XML database
343  /*!
344  * @param phaseRef XML phase node containing the description of the phase
345  * @param id id attribute containing the name of the phase.
346  * (default is the empty string)
347  */
348  RedlichKisterVPSSTP(XML_Node& phaseRef, std::string id = "");
349 
350 
351  //! Special constructor for a hard-coded problem
352  /*!
353  *
354  * @param testProb Hard-coded value. Only the value of 1 is
355  * used. It's for
356  * a LiKCl system
357  * -> test to predict the eutectic and liquidus correctly.
358  */
359  RedlichKisterVPSSTP(int testProb);
360 
361  //! Copy constructor
362  /*!
363  * Note this stuff will not work until the underlying phase
364  * has a working copy constructor
365  *
366  * @param b class to be copied
367  */
369 
370  //! Assignment operator
371  /*!
372  *
373  * @param b class to be copied.
374  */
376 
377  //! Destructor
378  virtual ~RedlichKisterVPSSTP();
379 
380  //! Duplication routine for objects which inherit from ThermoPhase.
381  /*!
382  * This virtual routine can be used to duplicate thermophase objects
383  * inherited from ThermoPhase even if the application only has
384  * a pointer to ThermoPhase to work with.
385  */
386  virtual ThermoPhase* duplMyselfAsThermoPhase() const;
387 
388  /**
389  *
390  * @name Utilities
391  * @{
392  */
393 
394 
395  //! Equation of state type flag.
396  /*!
397  * The ThermoPhase base class returns
398  * zero. Subclasses should define this to return a unique
399  * non-zero value. Known constants defined for this purpose are
400  * listed in mix_defs.h. The MolalityVPSSTP class also returns
401  * zero, as it is a non-complete class.
402  */
403  virtual int eosType() const;
404 
405  //! Initialization of a phase using an xml file
406  /*!
407  * This routine is a precursor to
408  * routine, which does most of the work.
409  *
410  * @param inputFile XML file containing the description of the
411  * phase
412  *
413  * @param id Optional parameter identifying the name of the
414  * phase. If none is given, the first XML
415  * phase element will be used.
416  */
417  void constructPhaseFile(std::string inputFile, std::string id);
418 
419  //! Import and initialize a phase
420  //! specification in an XML tree into the current object.
421  /*!
422  * Here we read an XML description of the phase.
423  * We import descriptions of the elements that make up the
424  * species in a phase.
425  * We import information about the species, including their
426  * reference state thermodynamic polynomials. We then freeze
427  * the state of the species.
428  *
429  * Then, we read the species molar volumes from the xml
430  * tree to finish the initialization.
431  *
432  * @param phaseNode This object must be the phase node of a
433  * complete XML tree
434  * description of the phase, including all of the
435  * species data. In other words while "phase" must
436  * point to an XML phase object, it must have
437  * sibling nodes "speciesData" that describe
438  * the species in the phase.
439  *
440  * @param id ID of the phase. If nonnull, a check is done
441  * to see if phaseNode is pointing to the phase
442  * with the correct id.
443  */
444  void constructPhaseXML(XML_Node& phaseNode, std::string id);
445 
446  /**
447  * @}
448  * @name Molar Thermodynamic Properties
449  * @{
450  */
451 
452 
453  /**
454  * @}
455  * @name Utilities for Solvent ID and Molality
456  * @{
457  */
458 
459 
460 
461 
462  /**
463  * @}
464  * @name Mechanical Properties
465  * @{
466  */
467 
468  /**
469  * @}
470  * @name Potential Energy
471  *
472  * Species may have an additional potential energy due to the
473  * presence of external gravitation or electric fields. These
474  * methods allow specifying a potential energy for individual
475  * species.
476  * @{
477  */
478 
479  /**
480  * @}
481  * @name Activities, Standard States, and Activity Concentrations
482  *
483  * The activity \f$a_k\f$ of a species in solution is
484  * related to the chemical potential by \f[ \mu_k = \mu_k^0(T)
485  * + \hat R T \log a_k. \f] The quantity \f$\mu_k^0(T,P)\f$ is
486  * the chemical potential at unit activity, which depends only
487  * on temperature and pressure.
488  * @{
489  */
490 
491  //! Get the array of non-dimensional molar-based ln activity coefficients at
492  //! the current solution temperature, pressure, and solution concentration.
493  /*!
494  * @param lnac Output vector of ln activity coefficients. Length: m_kk.
495  */
496  virtual void getLnActivityCoefficients(doublereal* lnac) const;
497 
498  //@}
499  /// @name Partial Molar Properties of the Solution
500  //@{
501 
502  //! Get the species chemical potentials. Units: J/kmol.
503  /*!
504  * This function returns a vector of chemical potentials of the
505  * species in solution at the current temperature, pressure
506  * and mole fraction of the solution.
507  *
508  * @param mu Output vector of species chemical
509  * potentials. Length: m_kk. Units: J/kmol
510  */
511  virtual void getChemPotentials(doublereal* mu) const;
512 
513  /// Molar enthalpy. Units: J/kmol.
514  virtual doublereal enthalpy_mole() const;
515 
516  /// Molar entropy. Units: J/kmol.
517  virtual doublereal entropy_mole() const;
518 
519  /// Molar heat capacity at constant pressure. Units: J/kmol/K.
520  virtual doublereal cp_mole() const;
521 
522  /// Molar heat capacity at constant volume. Units: J/kmol/K.
523  virtual doublereal cv_mole() const;
524 
525  //! Returns an array of partial molar enthalpies for the species
526  //! in the mixture.
527  /*!
528  * Units (J/kmol)
529  *
530  * For this phase, the partial molar enthalpies are equal to the
531  * standard state enthalpies modified by the derivative of the
532  * molality-based activity coefficient wrt temperature
533  *
534  * \f[
535  * \bar h_k(T,P) = h^o_k(T,P) - R T^2 \frac{d \ln(\gamma_k)}{dT}
536  * \f]
537  *
538  * @param hbar Vector of returned partial molar enthalpies
539  * (length m_kk, units = J/kmol)
540  */
541  virtual void getPartialMolarEnthalpies(doublereal* hbar) const;
542 
543  //! Returns an array of partial molar entropies for the species
544  //! in the mixture.
545  /*!
546  * Units (J/kmol)
547  *
548  * For this phase, the partial molar enthalpies are equal to the
549  * standard state enthalpies modified by the derivative of the
550  * activity coefficient wrt temperature
551  *
552  * \f[
553  * \bar s_k(T,P) = s^o_k(T,P) - R T^2 \frac{d \ln(\gamma_k)}{dT}
554  * - R \ln( \gamma_k X_k)
555  * - R T \frac{d \ln(\gamma_k) }{dT}
556  * \f]
557  *
558  * @param sbar Vector of returned partial molar entropies
559  * (length m_kk, units = J/kmol/K)
560  */
561  virtual void getPartialMolarEntropies(doublereal* sbar) const;
562 
563  //! Returns an array of partial molar entropies for the species
564  //! in the mixture.
565  /*!
566  * Units (J/kmol)
567  *
568  * For this phase, the partial molar enthalpies are equal to the
569  * standard state enthalpies modified by the derivative of the
570  * activity coefficient wrt temperature
571  *
572  * \f[
573  * ???????????????
574  * \bar s_k(T,P) = s^o_k(T,P) - R T^2 \frac{d \ln(\gamma_k)}{dT}
575  * - R \ln( \gamma_k X_k)
576  * - R T \frac{d \ln(\gamma_k) }{dT}
577  * ???????????????
578  * \f]
579  *
580  * @param cpbar Vector of returned partial molar heat capacities
581  * (length m_kk, units = J/kmol/K)
582  */
583  virtual void getPartialMolarCp(doublereal* cpbar) const;
584 
585 
586  //! Return an array of partial molar volumes for the
587  //! species in the mixture. Units: m^3/kmol.
588  /*!
589  * Frequently, for this class of thermodynamics representations,
590  * the excess Volume due to mixing is zero. Here, we set it as
591  * a default. It may be overridden in derived classes.
592  *
593  * @param vbar Output vector of species partial molar volumes.
594  * Length = m_kk. units are m^3/kmol.
595  */
596  virtual void getPartialMolarVolumes(doublereal* vbar) const;
597 
598  //! Get the species electrochemical potentials.
599  /*!
600  * These are partial molar quantities.
601  * This method adds a term \f$ Fz_k \phi_k \f$ to the
602  * to each chemical potential.
603  *
604  * Units: J/kmol
605  *
606  * @param mu output vector containing the species electrochemical potentials.
607  * Length: m_kk., units = J/kmol
608  */
609  void getElectrochemPotentials(doublereal* mu) const;
610 
611  //! Get the array of temperature second derivatives of the log activity coefficients
612  /*!
613  * This function is a virtual class, but it first appears in GibbsExcessVPSSTP
614  * class and derived classes from GibbsExcessVPSSTP.
615  *
616  * units = 1/Kelvin
617  *
618  * @param d2lnActCoeffdT2 Output vector of temperature 2nd derivatives of the
619  * log Activity Coefficients. length = m_kk
620  *
621  */
622  virtual void getd2lnActCoeffdT2(doublereal* d2lnActCoeffdT2) const;
623 
624  //! Get the array of temperature derivatives of the log activity coefficients
625  /*!
626  * This function is a virtual class, but it first appears in GibbsExcessVPSSTP
627  * class and derived classes from GibbsExcessVPSSTP.
628  *
629  * units = 1/Kelvin
630  *
631  * @param dlnActCoeffdT Output vector of temperature derivatives of the
632  * log Activity Coefficients. length = m_kk
633  *
634  */
635  virtual void getdlnActCoeffdT(doublereal* dlnActCoeffdT) const;
636 
637 
638 
639  //@}
640  /// @name Properties of the Standard State of the Species in the Solution
641  //@{
642 
643 
644 
645  //@}
646  /// @name Thermodynamic Values for the Species Reference States
647  //@{
648 
649 
650  ///////////////////////////////////////////////////////
651  //
652  // The methods below are not virtual, and should not
653  // be overloaded.
654  //
655  //////////////////////////////////////////////////////
656 
657  /**
658  * @name Specific Properties
659  * @{
660  */
661 
662 
663  /**
664  * @name Setting the State
665  *
666  * These methods set all or part of the thermodynamic
667  * state.
668  * @{
669  */
670 
671 
672 
673  //@}
674 
675  /**
676  * @name Chemical Equilibrium
677  * Routines that implement the Chemical equilibrium capability
678  * for a single phase, based on the element-potential method.
679  * @{
680  */
681 
682 
683 
684  //@}
685 
686 
687 
688  /// The following methods are used in the process of constructing
689  /// the phase and setting its parameters from a specification in an
690  /// input file. They are not normally used in application programs.
691  /// To see how they are used, see files importCTML.cpp and
692  /// ThermoFactory.cpp.
693 
694 
695  /*!
696  * @internal Initialize. This method is provided to allow
697  * subclasses to perform any initialization required after all
698  * species have been added. For example, it might be used to
699  * resize internal work arrays that must have an entry for
700  * each species. The base class implementation does nothing,
701  * and subclasses that do not require initialization do not
702  * need to overload this method. When importing a CTML phase
703  * description, this method is called just prior to returning
704  * from function importPhase.
705  *
706  * @see importCTML.cpp
707  */
708  virtual void initThermo();
709 
710 
711  /**
712  * Import and initialize a ThermoPhase object
713  *
714  * @param phaseNode This object must be the phase node of a
715  * complete XML tree
716  * description of the phase, including all of the
717  * species data. In other words while "phase" must
718  * point to an XML phase object, it must have
719  * sibling nodes "speciesData" that describe
720  * the species in the phase.
721  * @param id ID of the phase. If nonnull, a check is done
722  * to see if phaseNode is pointing to the phase
723  * with the correct id.
724  */
725  void initThermoXML(XML_Node& phaseNode, std::string id);
726 
727  /**
728  * @}
729  * @name Derivatives of Thermodynamic Variables needed for Applications
730  * @{
731  */
732 
733  //! Get the change in activity coefficients w.r.t. change in state (temp, mole fraction, etc.) along
734  //! a line in parameter space or along a line in physical space
735  /*!
736  *
737  * @param dTds Input of temperature change along the path
738  * @param dXds Input vector of changes in mole fraction along the path. length = m_kk
739  * Along the path length it must be the case that the mole fractions sum to one.
740  * @param dlnActCoeffds Output vector of the directional derivatives of the
741  * log Activity Coefficients along the path. length = m_kk
742  * units are 1/units(s). if s is a physical coordinate then the units are 1/m.
743  */
744  virtual void getdlnActCoeffds(const doublereal dTds, const doublereal* const dXds, doublereal* dlnActCoeffds) const;
745 
746  //! Get the array of log concentration-like derivatives of the
747  //! log activity coefficients - diagonal component
748  /*!
749  * This function is a virtual method. For ideal mixtures
750  * (unity activity coefficients), this can return zero.
751  * Implementations should take the derivative of the
752  * logarithm of the activity coefficient with respect to the
753  * logarithm of the mole fraction.
754  *
755  * units = dimensionless
756  *
757  * @param dlnActCoeffdlnX_diag Output vector of the diagonal component of the log(mole fraction)
758  * derivatives of the log Activity Coefficients.
759  * length = m_kk
760  */
761  virtual void getdlnActCoeffdlnX_diag(doublereal* dlnActCoeffdlnX_diag) const;
762 
763  //! Get the array of derivatives of the log activity coefficients wrt mole numbers - diagonal only
764  /*!
765  * This function is a virtual method. For ideal mixtures
766  * (unity activity coefficients), this can return zero.
767  * Implementations should take the derivative of the
768  * logarithm of the activity coefficient with respect to the
769  * logarithm of the concentration-like variable (i.e. mole fraction,
770  * molality, etc.) that represents the standard state.
771  *
772  * units = dimensionless
773  *
774  * @param dlnActCoeffdlnN_diag Output vector of the diagonal entries for the log(mole fraction)
775  * derivatives of the log Activity Coefficients.
776  * length = m_kk
777  */
778  virtual void getdlnActCoeffdlnN_diag(doublereal* dlnActCoeffdlnN_diag) const;
779 
780 
781  //! Get the array of derivatives of the ln activity coefficients with respect to the ln species mole numbers
782  /*!
783  * Implementations should take the derivative of the logarithm of the activity coefficient with respect to a
784  * log of a species mole number (with all other species mole numbers held constant)
785  *
786  * units = 1 / kmol
787  *
788  * dlnActCoeffdlnN[ ld * k + m] will contain the derivative of log act_coeff for the <I>m</I><SUP>th</SUP>
789  * species with respect to the number of moles of the <I>k</I><SUP>th</SUP> species.
790  *
791  * \f[
792  * \frac{d \ln(\gamma_m) }{d \ln( n_k ) }\Bigg|_{n_i}
793  * \f]
794  *
795  * @param ld Number of rows in the matrix
796  * @param dlnActCoeffdlnN Output vector of derivatives of the
797  * log Activity Coefficients. length = m_kk * m_kk
798  */
799  virtual void getdlnActCoeffdlnN(const size_t ld, doublereal* const dlnActCoeffdlnN) ;
800 
801  //@}
802 
803 private:
804 
805  //! Process an XML node called "binaryNeutralSpeciesParameters"
806  /*!
807  * This node contains all of the parameters necessary to describe
808  * the Redlich-Kister model for a particular binary interaction.
809  * This function reads the XML file and writes the coefficients
810  * it finds to an internal data structures.
811  *
812  * @param xmlBinarySpecies Reference to the XML_Node named "binaryNeutralSpeciesParameters"
813  * containing the binary interaction
814  */
815  void readXMLBinarySpecies(XML_Node& xmlBinarySpecies);
816 
817  //! Resize internal arrays within the object that depend upon the number
818  //! of binary Redlich-Kister interaction terms
819  /*!
820  * @param num Number of binary Redlich-Kister interaction terms
821  */
822  void resizeNumInteractions(const size_t num);
823 
824 
825  //! Initialize lengths of local variables after all species have
826  //! been identified.
827  void initLengths();
828 
829  //! Update the activity coefficients
830  /*!
831  * This function will be called to update the internally stored
832  * natural logarithm of the activity coefficients
833  */
834  void s_update_lnActCoeff() const;
835 
836  //! Update the derivative of the log of the activity coefficients wrt T
837  /*!
838  * This function will be called to update the internally stored
839  * derivative of the natural logarithm of the activity coefficients
840  * wrt temperature.
841  */
842  void s_update_dlnActCoeff_dT() const;
843 
844  //! Internal routine that calculates the derivative of the activity coefficients wrt
845  //! the mole fractions.
846  /*!
847  * This routine calculates the the derivative of the activity coefficients wrt to mole fraction
848  * with all other mole fractions held constant. This is strictly not permitted. However, if the
849  * resulting matrix is multiplied by a permissible deltaX vector then everything is ok.
850  *
851  * This is the natural way to handle concentration derivatives in this routine.
852  */
853  void s_update_dlnActCoeff_dX_() const;
854 
855 #ifdef DEBUG_MODE
856 public:
857  //! Utility routine that calculates a literature expression
858  /*!
859  * @param VintOut Output contribution to the voltage corresponding to nonideal term
860  * @param voltsOut Output contribution to the voltage corresponding to nonideal term and mf term
861  */
862  void Vint(double& VintOut, double& voltsOut) ;
863 #endif
864 
865 private:
866  //! Error function
867  /*!
868  * Print an error string and exit
869  *
870  * @param msg Message to be printed
871  */
872  doublereal err(std::string msg) const;
873 
874 protected:
875 
876  //! number of binary interaction expressions
878 
879  //! vector of species indices representing species A in the interaction
880  /*!
881  * Each Redlich-Kister excess Gibbs free energy term involves two species, A and B.
882  * This vector identifies species A.
883  */
884  std::vector<size_t> m_pSpecies_A_ij;
885 
886  //! vector of species indices representing species B in the interaction
887  /*!
888  * Each Redlich-Kister excess Gibbs free energy term involves two species, A and B.
889  * This vector identifies species B.
890  */
891  std::vector<size_t> m_pSpecies_B_ij;
892 
893 
894  //! Vector of the length of the polynomial for the interaction.
895  std::vector<size_t> m_N_ij;
896 
897 
898  //! Enthalpy term for the binary mole fraction interaction of the
899  //! excess gibbs free energy expression
900  mutable std::vector< vector_fp> m_HE_m_ij;
901 
902 
903  //! Entropy term for the binary mole fraction interaction of the
904  //! excess gibbs free energy expression
905  mutable std::vector< vector_fp> m_SE_m_ij;
906 
907  //! form of the RedlichKister interaction expression
908  /*!
909  * Currently there is only one form.
910  */
912 
913  //! form of the temperature dependence of the Redlich-Kister interaction expression
914  /*!
915  * Currently there is only one form -> constant wrt temperature.
916  */
918 
919 
920  //! Two dimensional array of derivatives of activity coefficients wrt mole fractions
922 
923 
924 };
925 
926 
927 
928 }
929 
930 #endif
931 
932 
933 
934 
935