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