Cantera  2.0
RedlichKwongMFTP.h
Go to the documentation of this file.
1 /**
2  * @file RedlichKwongMFTP.h
3  * Definition file for a derived class of ThermoPhase that assumes either
4  * an ideal gas or ideal solution approximation and handles
5  * variable pressure standard state methods for calculating
6  * thermodynamic properties (see \ref thermoprops and
7  * class \link Cantera::RedlichKwongMFTP RedlichKwongMFTP\endlink).
8  */
9 /*
10  * Copyright (2005) Sandia Corporation. Under the terms of
11  * Contract DE-AC04-94AL85000 with Sandia Corporation, the
12  * U.S. Government retains certain rights in this software.
13  */
14 
15 #ifndef CT_REDLICHKWONGMFTP_H
16 #define CT_REDLICHKWONGMFTP_H
17 
18 #include "MixtureFugacityTP.h"
19 
20 namespace Cantera
21 {
22 
23 class XML_Node;
24 class PDSS;
25 
26 /*!
27  * @name CONSTANTS - Models for the Standard State of IdealSolnPhase's
28  */
29 //@{
30 
31 /**
32  * @ingroup thermoprops
33  *
34  * This class can handle either an ideal solution or an ideal gas approximation
35  * of a phase.
36  *
37  *
38  * @nosubgrouping
39  */
41 {
42 
43 public:
44 
45  /*!
46  *
47  * @name Constructors and Duplicators for %RedlichKwongMFTP
48  *
49  */
50  //! Base constructor.
52 
53  //! Construct and initialize a RedlichKwongMFTP ThermoPhase object
54  //! directly from an ASCII input file
55  /*!
56  * Working constructors
57  *
58  * The two constructors below are the normal way the phase initializes itself. They are shells that call
59  * the routine initThermo(), with a reference to the
60  * XML database to get the info for the phase.
61  *
62  * @param infile Name of the input file containing the phase XML data
63  * to set up the object
64  * @param id ID of the phase in the input file. Defaults to the empty string.
65  */
66  RedlichKwongMFTP(std::string infile, std::string id="");
67 
68  //! Construct and initialize a RedlichKwongMFTP ThermoPhase object
69  //! directly from an XML database
70  /*!
71  * @param phaseRef XML phase node containing the description of the phase
72  * @param id id attribute containing the name of the phase. (default is the empty string)
73  */
74  RedlichKwongMFTP(XML_Node& phaseRef, std::string id = "");
75 
76  //! This is a special constructor, used to replicate test problems
77  //! during the initial verification of the object
78  /*!
79  *
80  * test problems:
81  * 1: Pure CO2 problem
82  * input file = CO2_RedlickKwongMFTP.xml
83  *
84  * @param testProb Hard -coded test problem to instantiate.
85  * Current valid values are 1.
86  */
87  RedlichKwongMFTP(int testProb);
88 
89  //! Copy Constructor
90  /*!
91  * Copy constructor for the object. Constructed object will be a clone of this object, but will
92  * also own all of its data. This is a wrapper around the assignment operator
93  *
94  * @param right Object to be copied.
95  */
96  RedlichKwongMFTP(const RedlichKwongMFTP& right);
97 
98  //! Assignment operator
99  /*!
100  * Assignment operator for the object. Constructed object will be a clone of this object, but will
101  * also own all of its data.
102  *
103  * @param right Object to be copied.
104  */
106 
107  //! Destructor.
108  virtual ~RedlichKwongMFTP();
109 
110 
111  //! Duplicator from the ThermoPhase parent class
112  /*!
113  * Given a pointer to a ThermoPhase object, this function will
114  * duplicate the ThermoPhase object and all underlying structures.
115  * This is basically a wrapper around the copy constructor.
116  *
117  * @return returns a pointer to a ThermoPhase
118  */
119  virtual ThermoPhase* duplMyselfAsThermoPhase() const;
120 
121  //@}
122 
123  /**
124  * @name Utilities (RedlichKwongMFTP)
125  */
126  //@{
127  /**
128  * Equation of state type flag. The base class returns
129  * zero. Subclasses should define this to return a unique
130  * non-zero value. Constants defined for this purpose are
131  * listed in mix_defs.h.
132  */
133  virtual int eosType() const;
134 
135  //@}
136 
137  /// Molar enthalpy. Units: J/kmol.
138  virtual doublereal enthalpy_mole() const;
139 
140  /// Molar internal energy. Units: J/kmol.
141  virtual doublereal intEnergy_mole() const;
142 
143  /// Molar entropy. Units: J/kmol/K.
144  virtual doublereal entropy_mole() const;
145 
146  /// Molar Gibbs function. Units: J/kmol.
147  virtual doublereal gibbs_mole() const;
148 
149  /// Molar heat capacity at constant pressure. Units: J/kmol/K.
150  virtual doublereal cp_mole() const;
151 
152  /// Molar heat capacity at constant volume. Units: J/kmol/K.
153  virtual doublereal cv_mole() const;
154 
155  /**
156  * @}
157  * @name Mechanical Properties
158  * @{
159  */
160 
161  //! Return the thermodynamic pressure (Pa).
162  /*!
163  * Since the mass density, temperature, and mass fractions are stored,
164  * this method uses these values to implement the
165  * mechanical equation of state \f$ P(T, \rho, Y_1, \dots, Y_K) \f$.
166  *
167  * \f[
168  * P = \frac{RT}{v-b_{mix}} - \frac{a_{mix}}{T^{0.5} v \left( v + b_{mix} \right) }
169  * \f]
170  *
171  */
172  virtual doublereal pressure() const;
173 
174  //! Returns the isothermal compressibility. Units: 1/Pa.
175  /*!
176  * The isothermal compressibility is defined as
177  * \f[
178  * \kappa_T = -\frac{1}{v}\left(\frac{\partial v}{\partial P}\right)_T
179  * \f]
180  */
181  virtual doublereal isothermalCompressibility() const;
182 
183 protected:
184  /**
185  * Calculate the density of the mixture using the partial
186  * molar volumes and mole fractions as input
187  *
188  * The formula for this is
189  *
190  * \f[
191  * \rho = \frac{\sum_k{X_k W_k}}{\sum_k{X_k V_k}}
192  * \f]
193  *
194  * where \f$X_k\f$ are the mole fractions, \f$W_k\f$ are
195  * the molecular weights, and \f$V_k\f$ are the pure species
196  * molar volumes.
197  *
198  * Note, the basis behind this formula is that in an ideal
199  * solution the partial molar volumes are equal to the
200  * species standard state molar volumes.
201  * The species molar volumes may be functions
202  * of temperature and pressure.
203  *
204  * NOTE: This is a non-virtual function, which is not a
205  * member of the ThermoPhase base class.
206  */
207  virtual void calcDensity();
208 
209 protected:
210  //! Set the temperature (K)
211  /*!
212  * Overwritten setTemperature(double) from State.h. This
213  * function sets the temperature, and makes sure that
214  * the value propagates to underlying objects
215  *
216  * @todo Make Phase::setTemperature a virtual function
217  *
218  * @param temp Temperature in kelvin
219  */
220  virtual void setTemperature(const doublereal temp);
221 
222  //! Set the mass fractions to the specified values, and then
223  //! normalize them so that they sum to 1.0.
224  /*!
225  * @param y Array of unnormalized mass fraction values (input).
226  * Must have a length greater than or equal to the number of species.
227  */
228  virtual void setMassFractions(const doublereal* const y);
229 
230  //!Set the mass fractions to the specified values without normalizing.
231  /*!
232  * This is useful when the normalization
233  * condition is being handled by some other means, for example
234  * by a constraint equation as part of a larger set of
235  * equations.
236  *
237  * @param y Input vector of mass fractions.
238  * Length is m_kk.
239  */
240  virtual void setMassFractions_NoNorm(const doublereal* const y);
241 
242  //! Set the mole fractions to the specified values, and then
243  //! normalize them so that they sum to 1.0.
244  /*!
245  * @param x Array of unnormalized mole fraction values (input).
246  * Must have a length greater than or equal to the number of species.
247  */
248  virtual void setMoleFractions(const doublereal* const x);
249 
250  //! Set the mole fractions to the specified values without normalizing.
251  /*!
252  * This is useful when the normalization
253  * condition is being handled by some other means, for example
254  * by a constraint equation as part of a larger set ofequations.
255  *
256  * @param x Input vector of mole fractions.
257  * Length is m_kk.
258  */
259  virtual void setMoleFractions_NoNorm(const doublereal* const x);
260 
261 
262  //! Set the concentrations to the specified values within the phase.
263  /*!
264  * @param c The input vector to this routine is in dimensional
265  * units. For volumetric phases c[k] is the
266  * concentration of the kth species in kmol/m3.
267  * For surface phases, c[k] is the concentration
268  * in kmol/m2. The length of the vector is the number
269  * of species in the phase.
270  */
271  virtual void setConcentrations(const doublereal* const c);
272 
273 
274 public:
275 
276  //! This method returns an array of generalized concentrations
277  /*!
278  * \f$ C^a_k\f$ are defined such that \f$ a_k = C^a_k /
279  * C^0_k, \f$ where \f$ C^0_k \f$ is a standard concentration
280  * defined below and \f$ a_k \f$ are activities used in the
281  * thermodynamic functions. These activity (or generalized)
282  * concentrations are used
283  * by kinetics manager classes to compute the forward and
284  * reverse rates of elementary reactions. Note that they may
285  * or may not have units of concentration --- they might be
286  * partial pressures, mole fractions, or surface coverages,
287  * for example.
288  *
289  * @param c Output array of generalized concentrations. The
290  * units depend upon the implementation of the
291  * reaction rate expressions within the phase.
292  */
293  virtual void getActivityConcentrations(doublereal* c) const;
294 
295  //! Returns the standard concentration \f$ C^0_k \f$, which is used to normalize
296  //! the generalized concentration.
297  /*!
298  * This is defined as the concentration by which the generalized
299  * concentration is normalized to produce the activity.
300  * In many cases, this quantity will be the same for all species in a phase.
301  * Since the activity for an ideal gas mixture is
302  * simply the mole fraction, for an ideal gas \f$ C^0_k = P/\hat R T \f$.
303  *
304  * @param k Optional parameter indicating the species. The default
305  * is to assume this refers to species 0.
306  * @return
307  * Returns the standard Concentration in units of m3 kmol-1.
308  */
309  virtual doublereal standardConcentration(size_t k=0) const;
310 
311  //! Returns the natural logarithm of the standard
312  //! concentration of the kth species
313  /*!
314  * @param k index of the species. (defaults to zero)
315  */
316  virtual doublereal logStandardConc(size_t k=0) const;
317 
318  //! Returns the units of the standard and generalized concentrations.
319  /*!
320  * Note they have the same units, as their
321  * ratio is defined to be equal to the activity of the kth
322  * species in the solution, which is unitless.
323  *
324  * This routine is used in print out applications where the
325  * units are needed. Usually, MKS units are assumed throughout
326  * the program and in the XML input files.
327  *
328  * The base %ThermoPhase class assigns the default quantities
329  * of (kmol/m3) for all species.
330  * Inherited classes are responsible for overriding the default
331  * values if necessary.
332  *
333  * @param uA Output vector containing the units
334  * uA[0] = kmol units - default = 1
335  * uA[1] = m units - default = -nDim(), the number of spatial
336  * dimensions in the Phase class.
337  * uA[2] = kg units - default = 0;
338  * uA[3] = Pa(pressure) units - default = 0;
339  * uA[4] = Temperature units - default = 0;
340  * uA[5] = time units - default = 0
341  * @param k species index. Defaults to 0.
342  * @param sizeUA output int containing the size of the vector.
343  * Currently, this is equal to 6.
344  */
345  virtual void getUnitsStandardConc(double* uA, int k = 0, int sizeUA = 6) const;
346 
347  //! Get the array of non-dimensional activity coefficients at
348  //! the current solution temperature, pressure, and solution concentration.
349  /*!
350  * For all objects with the Mixture Fugacity approximation, we define the
351  * standard state as an ideal gas at the current temperature and pressure
352  * of the solution. The activities are based on this standard state.
353  *
354  * @param ac Output vector of activity coefficients. Length: m_kk.
355  */
356  virtual void getActivityCoefficients(doublereal* ac) const;
357 
358 
359  /// @name Partial Molar Properties of the Solution (RedlichKwongMFTP)
360  //@{
361 
362  //! Get the array of non-dimensional species chemical potentials.
363  //! These are partial molar Gibbs free energies.
364  /*!
365  * \f$ \mu_k / \hat R T \f$.
366  * Units: unitless
367  *
368  * We close the loop on this function, here, calling
369  * getChemPotentials() and then dividing by RT. No need for child
370  * classes to handle.
371  *
372  * @param mu Output vector of non-dimensional species chemical potentials
373  * Length: m_kk.
374  */
375  void getChemPotentials_RT(doublereal* mu) const;
376 
377  //! Get the species chemical potentials. Units: J/kmol.
378  /*!
379  * This function returns a vector of chemical potentials of the
380  * species in solution at the current temperature, pressure
381  * and mole fraction of the solution.
382  *
383  * @param mu Output vector of species chemical
384  * potentials. Length: m_kk. Units: J/kmol
385  */
386  virtual void getChemPotentials(doublereal* mu) const;
387 
388  //! Get the species partial molar enthalpies. Units: J/kmol.
389  /*!
390  * @param hbar Output vector of species partial molar enthalpies.
391  * Length: m_kk. units are J/kmol.
392  */
393  virtual void getPartialMolarEnthalpies(doublereal* hbar) const;
394 
395  //! Get the species partial molar entropies. Units: J/kmol/K.
396  /*!
397  * @param sbar Output vector of species partial molar entropies.
398  * Length = m_kk. units are J/kmol/K.
399  */
400  virtual void getPartialMolarEntropies(doublereal* sbar) const;
401 
402  //! Get the species partial molar enthalpies. Units: J/kmol.
403  /*!
404  * @param ubar Output vector of species partial molar internal energies.
405  * Length = m_kk. units are J/kmol.
406  */
407  virtual void getPartialMolarIntEnergies(doublereal* ubar) const;
408 
409  //! Get the partial molar heat capacities Units: J/kmol/K
410  /*!
411  * @param cpbar Output vector of species partial molar heat capacities
412  * at constant pressure.
413  * Length = m_kk. units are J/kmol/K.
414  */
415  virtual void getPartialMolarCp(doublereal* cpbar) const;
416 
417  //! Get the species partial molar volumes. Units: m^3/kmol.
418  /*!
419  * @param vbar Output vector of species partial molar volumes.
420  * Length = m_kk. units are m^3/kmol.
421  */
422  virtual void getPartialMolarVolumes(doublereal* vbar) const;
423 
424  //@}
425 
426  /*!
427  * @name Properties of the Standard State of the Species in the Solution
428  *
429  * Properties of the standard states are delegated to the VPSSMgr object.
430  * The values are cached within this object, and are not recalculated unless
431  * the temperature or pressure changes.
432  */
433  //@{
434 
435  //@}
436 
437  /// @name Thermodynamic Values for the Species Reference States (RedlichKwongMFTP)
438  /*!
439  * Properties of the reference states are delegated to the VPSSMgr object.
440  * The values are cached within this object, and are not recalculated unless
441  * the temperature or pressure changes.
442  */
443  //@{
444 
445  //@}
446 
447 
448 
449  //---------------------------------------------------------
450  /// @name Critical State Properties.
451  /// These methods are only implemented by some subclasses, and may
452  /// be moved out of ThermoPhase at a later date.
453 
454  //@{
455 
456  /// Critical temperature (K).
457  virtual doublereal critTemperature() const;
458 
459  /// Critical pressure (Pa).
460  virtual doublereal critPressure() const;
461 
462  /// Critical density (kg/m3).
463  virtual doublereal critDensity() const;
464  //@}
465 
466 
467 
468 
469 public:
470 
471  //! @name Initialization Methods - For Internal use (VPStandardState)
472  /*!
473  * The following methods are used in the process of constructing
474  * the phase and setting its parameters from a specification in an
475  * input file. They are not normally used in application programs.
476  * To see how they are used, see files importCTML.cpp and
477  * ThermoFactory.cpp.
478  */
479  //@{
480 
481 
482  //! Set equation of state parameter values from XML
483  //! entries.
484  /*!
485  * This method is called by function importPhase in
486  * file importCTML.cpp when processing a phase definition in
487  * an input file. It should be overloaded in subclasses to set
488  * any parameters that are specific to that particular phase
489  * model.
490  *
491  * @param thermoNode An XML_Node object corresponding to
492  * the "thermo" entry for this phase in the input file.
493  */
494  virtual void setParametersFromXML(const XML_Node& thermoNode);
495 
496  //! @internal Initialize the object
497  /*!
498  * This method is provided to allow
499  * subclasses to perform any initialization required after all
500  * species have been added. For example, it might be used to
501  * resize internal work arrays that must have an entry for
502  * each species. The base class implementation does nothing,
503  * and subclasses that do not require initialization do not
504  * need to overload this method. When importing a CTML phase
505  * description, this method is called just prior to returning
506  * from function importPhase().
507  *
508  * @see importCTML.cpp
509  */
510  virtual void initThermo();
511 
512 
513  //!This method is used by the ChemEquil equilibrium solver.
514  /*!
515  * It sets the state such that the chemical potentials satisfy
516  * \f[ \frac{\mu_k}{\hat R T} = \sum_m A_{k,m}
517  * \left(\frac{\lambda_m} {\hat R T}\right) \f] where
518  * \f$ \lambda_m \f$ is the element potential of element m. The
519  * temperature is unchanged. Any phase (ideal or not) that
520  * implements this method can be equilibrated by ChemEquil.
521  *
522  * @param lambda_RT Input vector of dimensionless element potentials
523  * The length is equal to nElements().
524  */
525  void setToEquilState(const doublereal* lambda_RT);
526 
527  //! Initialize a ThermoPhase object, potentially reading activity
528  //! coefficient information from an XML database.
529  /*!
530  *
531  * This routine initializes the lengths in the current object and
532  * then calls the parent routine.
533  * This method is provided to allow
534  * subclasses to perform any initialization required after all
535  * species have been added. For example, it might be used to
536  * resize internal work arrays that must have an entry for
537  * each species. The base class implementation does nothing,
538  * and subclasses that do not require initialization do not
539  * need to overload this method. When importing a CTML phase
540  * description, this method is called just prior to returning
541  * from function importPhase().
542  *
543  * @param phaseNode This object must be the phase node of a
544  * complete XML tree
545  * description of the phase, including all of the
546  * species data. In other words while "phase" must
547  * point to an XML phase object, it must have
548  * sibling nodes "speciesData" that describe
549  * the species in the phase.
550  * @param id ID of the phase. If nonnull, a check is done
551  * to see if phaseNode is pointing to the phase
552  * with the correct id.
553  */
554  virtual void initThermoXML(XML_Node& phaseNode, std::string id);
555 
556 private:
557  //! Read the pure species RedlichKwong input parameters
558  /*!
559  * @param pureFluidParam XML_Node for the pure fluid parameters
560  */
561  void readXMLPureFluid(XML_Node& pureFluidParam);
562 
563 
564  //! Apply mixing rules for a coefficients
566 
567 
568  //! Read the cross species RedlichKwong input parameters
569  /*!
570  * @param pureFluidParam XML_Node for the cross fluid parameters
571  */
572  void readXMLCrossFluid(XML_Node& pureFluidParam);
573 
574 
575 
576  //==============================================================================
577 private:
578  //! @internal Initialize the internal lengths in this object.
579  /*!
580  * Note this is not a virtual function and only handles
581  * this object
582  */
583  void initLengths();
584 
585  //==============================================================================
586  // Special functions inherited from MixtureFugacityTP
587 
588 protected:
589 
590  //! Calculate the deviation terms for the total entropy of the mixture from the
591  //! ideal gas mixture
592  /*!
593  * Here we use the current state conditions
594  *
595  * @return Returns the change in entropy in units of J kmol-1 K-1.
596  */
597  virtual doublereal sresid() const;
598 
599  // Calculate the deviation terms for the total enthalpy of the mixture from the
600  // ideal gas mixture
601  /*
602  * Here we use the current state conditions
603  *
604  * @return Returns the change in enthalpy in units of J kmol-1.
605  */
606  virtual doublereal hresid() const;
607 public:
608  //! Estimate for the molar volume of the liquid
609  /*!
610  * Note: this is only used as a starting guess for later routines that actually calculate an
611  * accurate value for the liquid molar volume.
612  * This routine doesn't change the state of the system.
613  *
614  * @param TKelvin temperature in kelvin
615  * @param pres Pressure in Pa. This is used as an initial guess. If the routine
616  * needs to change the pressure to find a stable liquid state, the
617  * new pressure is returned in this variable.
618  *
619  * @return Returns the estimate of the liquid volume.
620  */
621  virtual doublereal liquidVolEst(doublereal TKelvin, doublereal& pres) const;
622 
623 public:
624  //! Calculates the density given the temperature and the pressure and a guess at the density.
625  /*!
626  * Note, below T_c, this is a multivalued function. We do not cross the vapor dome in this.
627  * This is protected because it is called during setState_TP() routines. Infinite loops would result
628  * if it were not protected.
629  *
630  * -> why is this not const?
631  *
632  * parameters:
633  * @param TKelvin Temperature in Kelvin
634  * @param pressure Pressure in Pascals (Newton/m**2)
635  * @param phase int representing the phase whose density we are requesting. If we put
636  * a gas or liquid phase here, we will attempt to find a volume in that
637  * part of the volume space, only, in this routine. A value of FLUID_UNDEFINED
638  * means that we will accept anything.
639  *
640  * @param rhoguess Guessed density of the fluid. A value of -1.0 indicates that there
641  * is no guessed density
642  *
643  *
644  * @return We return the density of the fluid at the requested phase. If we have not found any
645  * acceptable density we return a -1. If we have found an accectable density at a
646  * different phase, we return a -2.
647  */
648  virtual doublereal densityCalc(doublereal TKelvin, doublereal pressure, int phase, doublereal rhoguess);
649 
650 public:
651  //! Return the value of the density at the liquid spinodal point (on the liquid side)
652  //! for the current temperature.
653  /*!
654  * @return returns the density with units of kg m-3
655  */
656  virtual doublereal densSpinodalLiquid() const;
657 
658 
659  //! Return the value of the density at the gas spinodal point (on the gas side)
660  //! for the current temperature.
661  /*!
662  * @return returns the density with units of kg m-3
663  */
664  virtual doublereal densSpinodalGas() const;
665 
666 
667 
668  //! Calculate the pressure given the temperature and the molar volume
669  /*!
670  * Calculate the pressure given the temperature and the molar volume
671  *
672  * @param TKelvin temperature in kelvin
673  * @param molarVol molar volume ( m3/kmol)
674  *
675  * @return Returns the pressure.
676  */
677  virtual doublereal pressureCalc(doublereal TKelvin, doublereal molarVol) const;
678 
679 
680  //! Calculate the pressure and the pressure derivative given the temperature and the molar volume
681  /*!
682  * Temperature and mole number are held constant
683  *
684  * @param TKelvin temperature in kelvin
685  * @param molarVol molar volume ( m3/kmol)
686  *
687  * @param presCalc Returns the pressure.
688  *
689  * @return Returns the derivative of the pressure wrt the molar volume
690  */
691  virtual doublereal dpdVCalc(doublereal TKelvin, doublereal molarVol, doublereal& presCalc) const;
692 
693 
694  //! Calculate dpdV and dpdT at the current conditions
695  /*!
696  * These are stored internally.
697  */
698  void pressureDerivatives() const;
699 
700 
701  virtual void updateMixingExpressions();
702 
703 
704  //! Update the a and b parameters
705  /*!
706  * The a and the b parameters depend on the mole fraction and the temperature.
707  * This function updates the internal numbers based on the state of the object.
708  */
709  void updateAB();
710 
711 
712  //! Calculate the a and the b parameters given the temperature
713  /*!
714  *
715  * This function doesn't change the internal state of the object, so it is a const
716  * function. It does use the stored mole fractions in the object.
717  *
718  * @param temp Temperature (TKelvin)
719  *
720  * @param aCalc (output) Returns the a value
721  * @param bCalc (output) Returns the b value.
722  */
723  void calculateAB(doublereal temp, doublereal& aCalc, doublereal& bCalc) const;
724 
725 
726  //=========================================================================================
727  // Special functions not inherited from MixtureFugacityTP
728 
729  doublereal da_dt() const;
730 
731  void calcCriticalConditions(doublereal a, doublereal b, doublereal a0_coeff, doublereal aT_coeff,
732  doublereal& pc, doublereal& tc, doublereal& vc) const;
733 
734 
735 
736  int NicholsSolve(double TKelvin, double pres, doublereal a, doublereal b,
737  doublereal Vroot[3]) const;
738 
739  //@}
740  //==============================================================================
741 protected:
742 
743  //! boolean indicating whether standard mixing rules are applied
744  /*!
745  * - 1 = Yes, there are standard cross terms in the a coefficient matrices.
746  * - 0 = No, there are nonstaandard cross terms in the a coefficient matrices.
747  */
749 
750  //! Form of the temperature parameterization
751  /*!
752  * - 0 = There is no temperature parameterization of a or b
753  * - 1 = The a_ij parameter is a linear function of the temperature
754  */
756 
757 
758  //! Value of b in the equation of state
759  /*!
760  * m_b is a function of the temperature and the mole fraction.
761  */
762  doublereal m_b_current;
763 
764  //! Value of a in the equation of state
765  /*!
766  * a_b is a function of the temperature and the mole fraction.
767  */
768  doublereal m_a_current;
769 
770 
771  vector_fp a_vec_Curr_;
772  vector_fp b_vec_Curr_;
773 
774  Array2D a_coeff_vec;
775 
776 
777  vector_fp m_pc_Species;
778  vector_fp m_tc_Species;
779  vector_fp m_vc_Species;
780 
781  int NSolns_;
782 
783  doublereal Vroot_[3];
784 
785 
786 
787  //! Temporary storage - length = m_kk.
788  mutable vector_fp m_pp;
789 
790  //! Temporary storage - length = m_kk.
791  mutable vector_fp m_tmpV;
792 
793  // mutable vector_fp m_tmpV2;
794 
795  // Partial molar volumes of the species
796  mutable vector_fp m_partialMolarVolumes;
797 
798 
799 
800  //! The derivative of the pressure wrt the volume
801  /*!
802  * Calculated at the current conditions
803  * temperature and mole number kept constant
804  */
805  mutable doublereal dpdV_;
806 
807  //! The derivative of the pressure wrt the temperature
808  /*!
809  * Calculated at the current conditions
810  * Total volume and mole number kept constant
811  */
812  mutable doublereal dpdT_;
813 
814  //! Vector of derivatives of pressure wrt mole number
815  /*!
816  * Calculated at the current conditions
817  * Total volume, temperature and other mole number kept constant
818  */
819  mutable vector_fp dpdni_;
820 
821 public:
822  //! Omega constant for a -> value of a in terms of critical properties
823  /*!
824  * this was calculated from a small nonlinear solve
825  */
826  static const doublereal omega_a;
827 
828  //! Omega constant for b
829  static const doublereal omega_b;
830 
831  //! Omega constant for the critical molar volume
832  static const doublereal omega_vc;
833 
834 
835 };
836 }
837 
838 #endif