Cantera  2.3.0
Phase.h
Go to the documentation of this file.
1 /**
2  * @file Phase.h
3  * Header file for class Phase.
4  */
5 
6 // This file is part of Cantera. See License.txt in the top-level directory or
7 // at http://www.cantera.org/license.txt for license and copyright information.
8 
9 #ifndef CT_PHASE_H
10 #define CT_PHASE_H
11 
14 #include "cantera/thermo/Species.h"
16 
17 namespace Cantera
18 {
19 
20 /**
21  * @defgroup phases Models of Phases of Matter
22  *
23  * These classes are used to represent the composition and state of a single
24  * phase of matter. Together these classes form the basis for describing the
25  * species and element compositions of a phase as well as the stoichiometry
26  * of each species, and for describing the current state of the phase. They do
27  * not in themselves contain Thermodynamic equation of state information.
28  * However, they do comprise all of the necessary background functionality to
29  * support thermodynamic calculations (see \ref thermoprops).
30  */
31 
32 //! Class Phase is the base class for phases of matter, managing the species and
33 //! elements in a phase, as well as the independent variables of temperature,
34 //! mass density, species mass/mole fraction, and other generalized forces and
35 //! intrinsic properties (such as electric potential) that define the
36 //! thermodynamic state.
37 /*!
38  *
39  * Class Phase provides information about the elements and species in a
40  * phase - names, index numbers (location in arrays), atomic or molecular
41  * weights, etc. The set of elements must include all those that compose the
42  * species, but may include additional elements.
43  *
44  * It also stores an array of species molecular weights, which are used to
45  * convert between mole and mass representations of the composition. For
46  * efficiency in mass/mole conversion, the vector of mass fractions divided
47  * by molecular weight \f$ Y_k/M_k \f$ is also stored.
48  *
49  * Class Phase is not usually used directly. Its primary use is as a base class
50  * for class ThermoPhase. It is not generally necessary to overloaded any of
51  * class Phase's methods, with the exception of incompressible phases. In that
52  * case, the density must be replaced by the pressure as the independent
53  * variable and functions such as setMassFraction within class Phase must
54  * actually now calculate the density (at constant T and P) instead of leaving
55  * it alone as befits an independent variable. This also applies for nearly-
56  * incompressible phases or phases which utilize standard states based on a
57  * T and P, in which case they need to overload these functions too.
58  *
59  * Class Phase contains a number of utility functions that will set the state
60  * of the phase in its entirety, by first setting the composition, then the
61  * temperature and then the density. An example of this is the function
62  * Phase::setState_TRY(double t, double dens, const double* y).
63  *
64  * Class Phase contains method for saving and restoring the full internal states
65  * of each phase. These are saveState() and restoreState(). These functions
66  * operate on a state vector, which is in general of length (2 + nSpecies()).
67  * The first two entries of the state vector are temperature and density.
68  *
69  * A species name may be referred to via three methods:
70  *
71  * - "speciesName"
72  * - "PhaseId:speciesName"
73  * - "phaseName:speciesName"
74  * .
75  *
76  * The first two methods of naming may not yield a unique species within
77  * complicated assemblies of %Cantera Phases.
78  *
79  * @todo
80  * - Make the concept of saving state vectors more general, so that it can
81  * handle other cases where there are additional internal state variables,
82  * such as the voltage, a potential energy, or a strain field.
83  * - Specify that the input mole, mass, and volume fraction vectors must sum
84  * to one on entry to the set state routines. Non-conforming mole/mass
85  * fraction vectors are not thermodynamically consistent. Moreover, unless
86  * we do this, the calculation of Jacobians will be altered whenever the
87  * treatment of non- conforming mole fractions is changed. Add setState
88  * functions corresponding to specifying mole numbers, which is actually
89  * what is being done (well one of the options, there are many) when non-
90  * conforming mole fractions are input. Note, we realize that most numerical
91  * Jacobian and some analytical Jacobians use non-conforming calculations.
92  * These can easily be changed to the set mole number setState functions.
93  *
94  * @ingroup phases
95  */
96 class Phase
97 {
98 public:
99  Phase(); //!< Default constructor.
100 
101  virtual ~Phase();
102  Phase(const Phase& right);
103  Phase& operator=(const Phase& right);
104 
105  //! Returns a const reference to the XML_Node that describes the phase.
106  /*!
107  * The XML_Node for the phase contains all of the input data used to set up
108  * the model for the phase during its initialization.
109  */
110  XML_Node& xml() const;
111 
112  //! Stores the XML tree information for the current phase
113  /*!
114  * This function now stores the complete XML_Node tree as read into the
115  * code via a file. This is needed to move around within the XML tree
116  * during construction of transport and kinetics mechanisms after copy
117  * construction operations.
118  *
119  * @param xmlPhase Reference to the XML node corresponding to the phase
120  */
121  void setXMLdata(XML_Node& xmlPhase);
122 
123  /*! @name Name and ID
124  * Class Phase contains two strings that identify a phase. The ID is the
125  * value of the ID attribute of the XML phase node that is used to
126  * initialize a phase when it is read. The name field is also initialized
127  * to the value of the ID attribute of the XML phase node.
128  *
129  * However, the name field may be changed to another value during the
130  * course of a calculation. For example, if a phase is located in two
131  * places, but has the same constitutive input, the IDs of the two phases
132  * will be the same, but the names of the two phases may be different.
133  *
134  * It is an error to have two phases in a single problem with the same name
135  * and ID (or the name from one phase being the same as the id of
136  * another phase). Thus, it is expected that there is a 1-1 correspondence
137  * between names and unique phases within a Cantera problem.
138  */
139  //!@{
140 
141  //! Return the string id for the phase.
142  std::string id() const;
143 
144  //! Set the string id for the phase.
145  /*!
146  * @param id String id of the phase
147  */
148  void setID(const std::string& id);
149 
150  //! Return the name of the phase.
151  /*!
152  * Names are unique within a Cantera problem.
153  */
154  std::string name() const;
155 
156  //! Sets the string name for the phase.
157  //! @param nm String name of the phase
158  void setName(const std::string& nm);
159  //!@} end group Name and ID
160 
161  //! @name Element and Species Information
162  //!@{
163 
164  //! Name of the element with index m.
165  //! @param m Element index.
166  std::string elementName(size_t m) const;
167 
168  //! Return the index of element named 'name'. The index is an integer
169  //! assigned to each element in the order it was added. Returns \ref npos
170  //! if the specified element is not found.
171  //! @param name Name of the element
172  size_t elementIndex(const std::string& name) const;
173 
174  //! Return a read-only reference to the vector of element names.
175  const std::vector<std::string>& elementNames() const;
176 
177  //! Atomic weight of element m.
178  //! @param m Element index
179  doublereal atomicWeight(size_t m) const;
180 
181  //! Entropy of the element in its standard state at 298 K and 1 bar
182  //! @param m Element index
183  doublereal entropyElement298(size_t m) const;
184 
185  //! Atomic number of element m.
186  //! @param m Element index
187  int atomicNumber(size_t m) const;
188 
189  //! Return the element constraint type
190  //! Possible types include:
191  //!
192  //! - `CT_ELEM_TYPE_TURNEDOFF -1`
193  //! - `CT_ELEM_TYPE_ABSPOS 0`
194  //! - `CT_ELEM_TYPE_ELECTRONCHARGE 1`
195  //! - `CT_ELEM_TYPE_CHARGENEUTRALITY 2`
196  //! - `CT_ELEM_TYPE_LATTICERATIO 3`
197  //! - `CT_ELEM_TYPE_KINETICFROZEN 4`
198  //! - `CT_ELEM_TYPE_SURFACECONSTRAINT 5`
199  //! - `CT_ELEM_TYPE_OTHERCONSTRAINT 6`
200  //!
201  //! The default is `CT_ELEM_TYPE_ABSPOS`.
202  //! @param m Element index
203  //! @returns the element type
204  int elementType(size_t m) const;
205 
206  //! Change the element type of the mth constraint
207  //! Reassigns an element type.
208  //! @param m Element index
209  //! @param elem_type New elem type to be assigned
210  //! @returns the old element type
211  int changeElementType(int m, int elem_type);
212 
213  //! Return a read-only reference to the vector of atomic weights.
214  const vector_fp& atomicWeights() const;
215 
216  //! Number of elements.
217  size_t nElements() const;
218 
219  //! Check that the specified element index is in range.
220  //! Throws an exception if m is greater than nElements()-1
221  void checkElementIndex(size_t m) const;
222 
223  //! Check that an array size is at least nElements().
224  //! Throws an exception if mm is less than nElements(). Used before calls
225  //! which take an array pointer.
226  void checkElementArraySize(size_t mm) const;
227 
228  //! Number of atoms of element \c m in species \c k.
229  //! @param k species index
230  //! @param m element index
231  doublereal nAtoms(size_t k, size_t m) const;
232 
233  //! Get a vector containing the atomic composition of species k
234  //! @param k species index
235  //! @param atomArray vector containing the atomic number in the species.
236  //! Length: m_mm
237  void getAtoms(size_t k, double* atomArray) const;
238 
239  //! Returns the index of a species named 'name' within the Phase object.
240  //! The first species in the phase will have an index 0, and the last one
241  //! will have an index of nSpecies() - 1.
242  //! @param name String name of the species. It may also be in the form
243  //! phaseName:speciesName
244  //! @return The index of the species. If the name is not found,
245  //! the value \ref npos is returned.
246  size_t speciesIndex(const std::string& name) const;
247 
248  //! Name of the species with index k
249  //! @param k index of the species
250  std::string speciesName(size_t k) const;
251 
252  //! Returns the expanded species name of a species, including the phase name
253  //! This is guaranteed to be unique within a Cantera problem.
254  //! @param k Species index within the phase
255  //! @return The "phaseName:speciesName" string
256  std::string speciesSPName(int k) const;
257 
258  //! Return a const reference to the vector of species names
259  const std::vector<std::string>& speciesNames() const;
260 
261  /// Returns the number of species in the phase
262  size_t nSpecies() const {
263  return m_kk;
264  }
265 
266  //! Check that the specified species index is in range.
267  //! Throws an exception if k is greater than nSpecies()-1
268  void checkSpeciesIndex(size_t k) const;
269 
270  //! Check that an array size is at least nSpecies().
271  //! Throws an exception if kk is less than nSpecies(). Used before calls
272  //! which take an array pointer.
273  void checkSpeciesArraySize(size_t kk) const;
274 
275  //!@} end group Element and Species Information
276 
277  //! Save the current internal state of the phase.
278  //! Write to vector 'state' the current internal state.
279  //! @param state output vector. Will be resized to nSpecies() + 2.
280  void saveState(vector_fp& state) const;
281 
282  //! Write to array 'state' the current internal state.
283  //! @param lenstate length of the state array. Must be >= nSpecies()+2
284  //! @param state output vector. Must be of length nSpecies() + 2 or
285  //! greater.
286  void saveState(size_t lenstate, doublereal* state) const;
287 
288  //! Restore a state saved on a previous call to saveState.
289  //! @param state State vector containing the previously saved state.
290  void restoreState(const vector_fp& state);
291 
292  //! Restore the state of the phase from a previously saved state vector.
293  //! @param lenstate Length of the state vector
294  //! @param state Vector of state conditions.
295  void restoreState(size_t lenstate, const doublereal* state);
296 
297  /*! @name Set thermodynamic state
298  * Set the internal thermodynamic state by setting the internally stored
299  * temperature, density and species composition. Note that the composition
300  * is always set first.
301  *
302  * Temperature and density are held constant if not explicitly set.
303  */
304  //!@{
305 
306  //! Set the species mole fractions by name.
307  //! Species not listed by name in \c xMap are set to zero.
308  //! @param xMap map from species names to mole fraction values.
309  void setMoleFractionsByName(const compositionMap& xMap);
310 
311  //! Set the mole fractions of a group of species by name. Species which
312  //! are not listed by name in the composition map are set to zero.
313  //! @param x string x in the form of a composition map
314  void setMoleFractionsByName(const std::string& x);
315 
316  //! Set the species mass fractions by name.
317  //! Species not listed by name in \c yMap are set to zero.
318  //! @param yMap map from species names to mass fraction values.
319  void setMassFractionsByName(const compositionMap& yMap);
320 
321  //! Set the species mass fractions by name.
322  //! Species not listed by name in \c x are set to zero.
323  //! @param x String containing a composition map
324  void setMassFractionsByName(const std::string& x);
325 
326  //! Set the internally stored temperature (K), density, and mole fractions.
327  //! @param t Temperature in kelvin
328  //! @param dens Density (kg/m^3)
329  //! @param x vector of species mole fractions, length m_kk
330  void setState_TRX(doublereal t, doublereal dens, const doublereal* x);
331 
332  //! Set the internally stored temperature (K), density, and mole fractions.
333  //! @param t Temperature in kelvin
334  //! @param dens Density (kg/m^3)
335  //! @param x Composition Map containing the mole fractions.
336  //! Species not included in the map are assumed to have
337  //! a zero mole fraction.
338  void setState_TRX(doublereal t, doublereal dens, const compositionMap& x);
339 
340  //! Set the internally stored temperature (K), density, and mass fractions.
341  //! @param t Temperature in kelvin
342  //! @param dens Density (kg/m^3)
343  //! @param y vector of species mass fractions, length m_kk
344  void setState_TRY(doublereal t, doublereal dens, const doublereal* y);
345 
346  //! Set the internally stored temperature (K), density, and mass fractions.
347  //! @param t Temperature in kelvin
348  //! @param dens Density (kg/m^3)
349  //! @param y Composition Map containing the mass fractions.
350  //! Species not included in the map are assumed to have
351  //! a zero mass fraction.
352  void setState_TRY(doublereal t, doublereal dens, const compositionMap& y);
353 
354  //! Set the internally stored temperature (K), molar density (kmol/m^3), and
355  //! mole fractions.
356  //! @param t Temperature in kelvin
357  //! @param n molar density (kmol/m^3)
358  //! @param x vector of species mole fractions, length m_kk
359  void setState_TNX(doublereal t, doublereal n, const doublereal* x);
360 
361  //! Set the internally stored temperature (K) and density (kg/m^3)
362  //! @param t Temperature in kelvin
363  //! @param rho Density (kg/m^3)
364  void setState_TR(doublereal t, doublereal rho);
365 
366  //! Set the internally stored temperature (K) and mole fractions.
367  //! @param t Temperature in kelvin
368  //! @param x vector of species mole fractions, length m_kk
369  void setState_TX(doublereal t, doublereal* x);
370 
371  //! Set the internally stored temperature (K) and mass fractions.
372  //! @param t Temperature in kelvin
373  //! @param y vector of species mass fractions, length m_kk
374  void setState_TY(doublereal t, doublereal* y);
375 
376  //! Set the density (kg/m^3) and mole fractions.
377  //! @param rho Density (kg/m^3)
378  //! @param x vector of species mole fractions, length m_kk
379  void setState_RX(doublereal rho, doublereal* x);
380 
381  //! Set the density (kg/m^3) and mass fractions.
382  //! @param rho Density (kg/m^3)
383  //! @param y vector of species mass fractions, length m_kk
384  void setState_RY(doublereal rho, doublereal* y);
385 
386  //!@} end group set thermo state
387 
388  //! Molecular weight of species \c k.
389  //! @param k index of species \c k
390  //! @returns the molecular weight of species \c k.
391  doublereal molecularWeight(size_t k) const;
392 
393  //! Copy the vector of molecular weights into vector weights.
394  //! @param weights Output vector of molecular weights (kg/kmol)
395  void getMolecularWeights(vector_fp& weights) const;
396 
397  //! Copy the vector of molecular weights into array weights.
398  //! @param weights Output array of molecular weights (kg/kmol)
399  void getMolecularWeights(doublereal* weights) const;
400 
401  //! Return a const reference to the internal vector of molecular weights.
402  //! units = kg / kmol
403  const vector_fp& molecularWeights() const;
404 
405  //! This routine returns the size of species k
406  //!
407  //! The meaning and dimensions are model-dependent. For surface phases, the
408  //! size is the number of sites occupied by one molecule of the species
409  //! [nondimensional]. For models which utilize the species partial molar
410  //! volumes, this is the molar volume of the species in its reference state.
411  //! For other models, this value may have no meaning.
412  //! @param k index of the species
413  //! @return The size of the species
414  doublereal size(size_t k) const {
415  return m_speciesSize[k];
416  }
417 
418  /// @name Composition
419  //@{
420 
421  //! Get the mole fractions by name.
422  //! @param threshold Exclude species with mole fractions less than or
423  //! equal to this threshold.
424  //! @return Map of species names to mole fractions
425  compositionMap getMoleFractionsByName(double threshold=0.0) const;
426 
427  //! Return the mole fraction of a single species
428  //! @param k species index
429  //! @return Mole fraction of the species
430  doublereal moleFraction(size_t k) const;
431 
432  //! Return the mole fraction of a single species
433  //! @param name String name of the species
434  //! @return Mole fraction of the species
435  doublereal moleFraction(const std::string& name) const;
436 
437  //! Get the mass fractions by name.
438  //! @param threshold Exclude species with mass fractions less than or
439  //! equal to this threshold.
440  //! @return Map of species names to mass fractions
441  compositionMap getMassFractionsByName(double threshold=0.0) const;
442 
443  //! Return the mass fraction of a single species
444  //! @param k species index
445  //! @return Mass fraction of the species
446  doublereal massFraction(size_t k) const;
447 
448  //! Return the mass fraction of a single species
449  //! @param name String name of the species
450  //! @return Mass Fraction of the species
451  doublereal massFraction(const std::string& name) const;
452 
453  //! Get the species mole fraction vector.
454  //! @param x On return, x contains the mole fractions. Must have a
455  //! length greater than or equal to the number of species.
456  void getMoleFractions(doublereal* const x) const;
457 
458  //! Set the mole fractions to the specified values.
459  //! There is no restriction on the sum of the mole fraction vector.
460  //! Internally, the Phase object will normalize this vector before storing
461  //! its contents.
462  //! @param x Array of unnormalized mole fraction values (input). Must
463  //! have a length greater than or equal to the number of species, m_kk.
464  virtual void setMoleFractions(const doublereal* const x);
465 
466  //! Set the mole fractions to the specified values without normalizing.
467  //! This is useful when the normalization condition is being handled by
468  //! some other means, for example by a constraint equation as part of a
469  //! larger set of equations.
470  //! @param x Input vector of mole fractions. Length is m_kk.
471  virtual void setMoleFractions_NoNorm(const doublereal* const x);
472 
473  //! Get the species mass fractions.
474  //! @param[out] y Array of mass fractions, length nSpecies()
475  void getMassFractions(doublereal* const y) const;
476 
477  //! Return a const pointer to the mass fraction array
478  const doublereal* massFractions() const {
479  return &m_y[0];
480  }
481 
482  //! Set the mass fractions to the specified values and normalize them.
483  //! @param[in] y Array of unnormalized mass fraction values. Length
484  //! must be greater than or equal to the number of
485  //! species. The Phase object will normalize this vector
486  //! before storing its contents.
487  virtual void setMassFractions(const doublereal* const y);
488 
489  //! Set the mass fractions to the specified values without normalizing.
490  //! This is useful when the normalization condition is being handled by
491  //! some other means, for example by a constraint equation as part of a
492  //! larger set of equations.
493  //! @param y Input vector of mass fractions. Length is m_kk.
494  virtual void setMassFractions_NoNorm(const doublereal* const y);
495 
496  //! Get the species concentrations (kmol/m^3).
497  /*!
498  * @param[out] c The vector of species concentrations. Units are
499  * kmol/m^3. The length of the vector must be greater than
500  * or equal to the number of species within the phase.
501  */
502  void getConcentrations(doublereal* const c) const;
503 
504  //! Concentration of species k.
505  //! If k is outside the valid range, an exception will be thrown.
506  /*!
507  * @param[in] k Index of the species within the phase.
508  *
509  * @returns the concentration of species k (kmol m-3).
510  */
511  doublereal concentration(const size_t k) const;
512 
513  //! Set the concentrations to the specified values within the phase.
514  //! We set the concentrations here and therefore we set the overall density
515  //! of the phase. We hold the temperature constant during this operation.
516  //! Therefore, we have possibly changed the pressure of the phase by
517  //! calling this routine.
518  //! @param[in] conc Array of concentrations in dimensional units. For
519  //! bulk phases c[k] is the concentration of the kth
520  //! species in kmol/m3. For surface phases, c[k] is the
521  //! concentration in kmol/m2. The length of the vector
522  //! is the number of species in the phase.
523  virtual void setConcentrations(const doublereal* const conc);
524 
525  //! Set the concentrations without ignoring negative concentrations
526  virtual void setConcentrationsNoNorm(const double* const conc);
527 
528  //! Elemental mass fraction of element m
529  /*!
530  * The elemental mass fraction \f$Z_{\mathrm{mass},m}\f$ of element \f$m\f$
531  * is defined as
532  * \f[
533  * Z_{\mathrm{mass},m} = \sum_k \frac{a_{m,k} M_m}{M_k} Y_k
534  * \f]
535  * with \f$a_{m,k}\f$ being the number of atoms of element \f$m\f$ in
536  * species \f$k\f$, \f$M_m\f$ the atomic weight of element \f$m\f$,
537  * \f$M_k\f$ the molecular weight of species \f$k\f$, and \f$Y_k\f$
538  * the mass fraction of species \f$k\f$.
539  *
540  * @param[in] m Index of the element within the phase. If m is outside
541  * the valid range, an exception will be thrown.
542  *
543  * @return the elemental mass fraction of element m.
544  */
545  doublereal elementalMassFraction(const size_t m) const;
546 
547  //! Elemental mole fraction of element m
548  /*!
549  * The elemental mole fraction \f$Z_{\mathrm{mole},m}\f$ of element \f$m\f$
550  * is the number of atoms of element *m* divided by the total number of
551  * atoms. It is defined as:
552  *
553  * \f[
554  * Z_{\mathrm{mole},m} = \frac{\sum_k a_{m,k} X_k}
555  * {\sum_k \sum_j a_{j,k} X_k}
556  * \f]
557  * with \f$a_{m,k}\f$ being the number of atoms of element \f$m\f$ in
558  * species \f$k\f$, \f$\sum_j\f$ being a sum over all elements, and
559  * \f$X_k\f$ being the mole fraction of species \f$k\f$.
560  *
561  * @param[in] m Index of the element within the phase. If m is outside the
562  * valid range, an exception will be thrown.
563  * @return the elemental mole fraction of element m.
564  */
565  doublereal elementalMoleFraction(const size_t m) const;
566 
567  //! Returns a const pointer to the start of the moleFraction/MW array.
568  //! This array is the array of mole fractions, each divided by the mean
569  //! molecular weight.
570  const doublereal* moleFractdivMMW() const;
571 
572  //@}
573 
574  //! Dimensionless electrical charge of a single molecule of species k
575  //! The charge is normalized by the the magnitude of the electron charge
576  //! @param k species index
577  doublereal charge(size_t k) const {
578  return m_speciesCharge[k];
579  }
580 
581  //! Charge density [C/m^3].
582  doublereal chargeDensity() const;
583 
584  //! Returns the number of spatial dimensions (1, 2, or 3)
585  size_t nDim() const {
586  return m_ndim;
587  }
588 
589  //! Set the number of spatial dimensions (1, 2, or 3). The number of
590  //! spatial dimensions is used for vector involving directions.
591  //! @param ndim Input number of dimensions.
592  void setNDim(size_t ndim) {
593  m_ndim = ndim;
594  }
595 
596  //! @name Thermodynamic Properties
597  //!@{
598 
599  //! Temperature (K).
600  //! @return The temperature of the phase
601  doublereal temperature() const {
602  return m_temp;
603  }
604 
605  //! Density (kg/m^3).
606  //! @return The density of the phase
607  virtual doublereal density() const {
608  return m_dens;
609  }
610 
611  //! Molar density (kmol/m^3).
612  //! @return The molar density of the phase
613  doublereal molarDensity() const;
614 
615  //! Molar volume (m^3/kmol).
616  //! @return The molar volume of the phase
617  doublereal molarVolume() const;
618 
619  //! Set the internally stored density (kg/m^3) of the phase.
620  //! Note the density of a phase is an independent variable.
621  //! @param[in] density_ density (kg/m^3).
622  virtual void setDensity(const doublereal density_) {
623  if (density_ > 0.0) {
624  m_dens = density_;
625  } else {
626  throw CanteraError("Phase::setDensity()",
627  "density must be positive");
628  }
629  }
630 
631  //! Set the internally stored molar density (kmol/m^3) of the phase.
632  //! @param[in] molarDensity Input molar density (kmol/m^3).
633  virtual void setMolarDensity(const doublereal molarDensity);
634 
635  //! Set the internally stored temperature of the phase (K).
636  //! @param temp Temperature in Kelvin
637  virtual void setTemperature(const doublereal temp) {
638  if (temp > 0) {
639  m_temp = temp;
640  } else {
641  throw CanteraError("Phase::setTemperature",
642  "temperature must be positive");
643  }
644  }
645  //@}
646 
647  //! @name Mean Properties
648  //!@{
649 
650  //! Evaluate the mole-fraction-weighted mean of an array Q.
651  //! \f[ \sum_k X_k Q_k. \f]
652  //! Q should contain pure-species molar property values.
653  //! @param[in] Q Array of length m_kk that is to be averaged.
654  //! @return mole-fraction-weighted mean of Q
655  doublereal mean_X(const doublereal* const Q) const;
656 
657  //! @copydoc Phase::mean_X(const doublereal* const Q) const
658  doublereal mean_X(const vector_fp& Q) const;
659 
660  //! The mean molecular weight. Units: (kg/kmol)
661  doublereal meanMolecularWeight() const {
662  return m_mmw;
663  }
664 
665  //! Evaluate \f$ \sum_k X_k \log X_k \f$.
666  //! @return The indicated sum. Dimensionless.
667  doublereal sum_xlogx() const;
668 
669  //@}
670 
671  //! @name Adding Elements and Species
672  //! These methods are used to add new elements or species. These are not
673  //! usually called by user programs.
674  //!
675  //! Since species are checked to insure that they are only composed of
676  //! declared elements, it is necessary to first add all elements before
677  //! adding any species.
678  //!@{
679 
680  //! Add an element.
681  //! @param symbol Atomic symbol std::string.
682  //! @param weight Atomic mass in amu.
683  //! @param atomicNumber Atomic number of the element (unitless)
684  //! @param entropy298 Entropy of the element at 298 K and 1 bar in its
685  //! most stable form. The default is the value ENTROPY298_UNKNOWN,
686  //! which is interpreted as an unknown, and if used will cause
687  //! %Cantera to throw an error.
688  //! @param elem_type Specifies the type of the element constraint
689  //! equation. This defaults to CT_ELEM_TYPE_ABSPOS, i.e., an element.
690  //! @return index of the element added
691  size_t addElement(const std::string& symbol, doublereal weight=-12345.0,
692  int atomicNumber=0, doublereal entropy298=ENTROPY298_UNKNOWN,
693  int elem_type=CT_ELEM_TYPE_ABSPOS);
694 
695  //! Add a Species to this Phase. Returns `true` if the species was
696  //! successfully added, or `false` if the species was ignored.
697  //!
698  //! Derived classes which need to size arrays according to the number of
699  //! species should overload this method. The derived class implementation
700  //! should call the base class method, and, if this returns `true`
701  //! (indicating that the species has been added), adjust their array sizes
702  //! accordingly.
703  //!
704  //! @see ignoreUndefinedElements addUndefinedElements throwUndefinedElements
705  virtual bool addSpecies(shared_ptr<Species> spec);
706 
707  //! Modify the thermodynamic data associated with a species.
708  /*!
709  * The species name, elemental composition, and type of thermo
710  * parameterization must be unchanged. If there are Kinetics objects that
711  * depend on this phase, Kinetics::invalidateCache() should be called on
712  * those objects after calling this function.
713  */
714  virtual void modifySpecies(size_t k, shared_ptr<Species> spec);
715 
716  //! Return the Species object for the named species.
717  shared_ptr<Species> species(const std::string& name) const;
718 
719  //! Return the Species object for species whose index is *k*.
720  shared_ptr<Species> species(size_t k) const;
721 
722  //! Set behavior when adding a species containing undefined elements to just
723  //! skip the species.
725 
726  //! Set behavior when adding a species containing undefined elements to add
727  //! those elements to the phase.
728  void addUndefinedElements();
729 
730  //! Set the behavior when adding a species containing undefined elements to
731  //! throw an exception. This is the default behavior.
732  void throwUndefinedElements();
733 
734  struct UndefElement { enum behavior {
735  error, ignore, add
736  }; };
737 
738  //!@} end group adding species and elements
739 
740  //! Returns a bool indicating whether the object is ready for use
741  /*!
742  * @returns true if the object is ready for calculation, false otherwise.
743  */
744  virtual bool ready() const;
745 
746  //! Return the State Mole Fraction Number
747  int stateMFNumber() const {
748  return m_stateNum;
749  }
750 
751  //! Invalidate any cached values which are normally updated only when a
752  //! change in state is detected
753  virtual void invalidateCache();
754 
755 protected:
756  //! Cached for saved calculations within each ThermoPhase.
757  /*!
758  * For more information on how to use this, see examples within the source
759  * code and documentation for this within ValueCache class itself.
760  */
762 
763  //! Set the molecular weight of a single species to a given value.
764  //!
765  //! Used by phases where the equation of state is defined for a specific
766  //! value of the molecular weight which may not exactly correspond to the
767  //! value computed from the chemical formula.
768  //! @param k id of the species
769  //! @param mw Molecular Weight (kg kmol-1)
770  void setMolecularWeight(const int k, const double mw) {
771  m_molwts[k] = mw;
772  m_rmolwts[k] = 1.0/mw;
773  }
774 
775  //! Apply changes to the state which are needed after the composition
776  //! changes. This function is called after any call to setMassFractions(),
777  //! setMoleFractions(), or similar. For phases which need to execute a
778  //! callback after any change to the composition, it should be done by
779  //! overriding this function rather than overriding all of the composition-
780  //! setting functions. Derived class implementations of compositionChanged()
781  //! should call the parent class method as well.
782  virtual void compositionChanged();
783 
784  size_t m_kk; //!< Number of species in the phase.
785 
786  //! Dimensionality of the phase. Volumetric phases have dimensionality 3
787  //! and surface phases have dimensionality 2.
788  size_t m_ndim;
789 
790  //! Atomic composition of the species. The number of atoms of element i
791  //! in species k is equal to m_speciesComp[k * m_mm + i]
792  //! The length of this vector is equal to m_kk * m_mm
794 
795  //!Vector of species sizes. length m_kk. Used in some equations of state
796  //! which employ the constant partial molar volume approximation, and for
797  //! surface phases.
799 
800  vector_fp m_speciesCharge; //!< Vector of species charges. length m_kk.
801 
802  std::map<std::string, shared_ptr<Species> > m_species;
803 
804  //! Flag determining behavior when adding species with an undefined element
805  UndefElement::behavior m_undefinedElementBehavior;
806 
807 private:
808  XML_Node* m_xml; //!< XML node containing the XML info for this phase
809 
810  //! ID of the phase. This is the value of the ID attribute of the XML
811  //! phase node. The field will stay that way even if the name is changed.
812  std::string m_id;
813 
814  //! Name of the phase.
815  //! Initially, this is the value of the ID attribute of the XML phase node.
816  //! It may be changed to another value during the course of a calculation.
817  std::string m_name;
818 
819  doublereal m_temp; //!< Temperature (K). This is an independent variable
820 
821  //! Density (kg m-3). This is an independent variable except in the
822  //! incompressible degenerate case. Thus, the pressure is determined from
823  //! this variable rather than other way round.
824  doublereal m_dens;
825 
826  doublereal m_mmw; //!< mean molecular weight of the mixture (kg kmol-1)
827 
828  //! m_ym[k] = mole fraction of species k divided by the mean molecular
829  //! weight of mixture.
830  mutable vector_fp m_ym;
831 
832  //! Mass fractions of the species
833  /*!
834  * Note, this vector
835  * Length is m_kk
836  */
837  mutable vector_fp m_y;
838 
839  vector_fp m_molwts; //!< species molecular weights (kg kmol-1)
840 
841  vector_fp m_rmolwts; //!< inverse of species molecular weights (kmol kg-1)
842 
843  //! State Change variable. Whenever the mole fraction vector changes,
844  //! this int is incremented.
846 
847  //! Vector of the species names
848  std::vector<std::string> m_speciesNames;
849 
850  //! Map of species names to indices
851  std::map<std::string, size_t> m_speciesIndices;
852 
853  size_t m_mm; //!< Number of elements.
854  vector_fp m_atomicWeights; //!< element atomic weights (kg kmol-1)
855  vector_int m_atomicNumbers; //!< element atomic numbers
856  std::vector<std::string> m_elementNames; //!< element names
857  vector_int m_elem_type; //!< Vector of element types
858 
859  //! Entropy at 298.15 K and 1 bar of stable state pure elements (J kmol-1)
861 };
862 
863 }
864 
865 #endif
Phase()
Default constructor.
Definition: Phase.cpp:20
std::map< std::string, doublereal > compositionMap
Map connecting a string name with a double.
Definition: ct_defs.h:149
const vector_fp & molecularWeights() const
Return a const reference to the internal vector of molecular weights.
Definition: Phase.cpp:513
virtual void setMoleFractions_NoNorm(const doublereal *const x)
Set the mole fractions to the specified values without normalizing.
Definition: Phase.cpp:359
vector_fp m_y
Mass fractions of the species.
Definition: Phase.h:837
void setState_RY(doublereal rho, doublereal *y)
Set the density (kg/m^3) and mass fractions.
Definition: Phase.cpp:490
doublereal molarVolume() const
Molar volume (m^3/kmol).
Definition: Phase.cpp:676
std::vector< std::string > m_elementNames
element names
Definition: Phase.h:856
size_t nElements() const
Number of elements.
Definition: Phase.cpp:161
void restoreState(const vector_fp &state)
Restore a state saved on a previous call to saveState.
Definition: Phase.cpp:309
void setState_TRY(doublereal t, doublereal dens, const doublereal *y)
Set the internally stored temperature (K), density, and mass fractions.
Definition: Phase.cpp:452
void getMassFractions(doublereal *const y) const
Get the species mass fractions.
Definition: Phase.cpp:584
vector_fp m_atomicWeights
element atomic weights (kg kmol-1)
Definition: Phase.h:854
doublereal temperature() const
Temperature (K).
Definition: Phase.h:601
const doublereal * massFractions() const
Return a const pointer to the mass fraction array.
Definition: Phase.h:478
size_t speciesIndex(const std::string &name) const
Returns the index of a species named &#39;name&#39; within the Phase object.
Definition: Phase.cpp:251
void saveState(vector_fp &state) const
Save the current internal state of the phase.
Definition: Phase.cpp:297
doublereal moleFraction(size_t k) const
Return the mole fraction of a single species.
Definition: Phase.cpp:547
vector_fp m_speciesCharge
Vector of species charges. length m_kk.
Definition: Phase.h:800
size_t elementIndex(const std::string &name) const
Return the index of element named &#39;name&#39;.
Definition: Phase.cpp:186
doublereal sum_xlogx() const
Evaluate .
Definition: Phase.cpp:700
vector_fp m_entropy298
Entropy at 298.15 K and 1 bar of stable state pure elements (J kmol-1)
Definition: Phase.h:860
compositionMap getMoleFractionsByName(double threshold=0.0) const
Get the mole fractions by name.
Definition: Phase.cpp:518
Class XML_Node is a tree-based representation of the contents of an XML file.
Definition: xml.h:97
Contains the getElementWeight function and the definitions of element constraint types.
std::map< std::string, size_t > m_speciesIndices
Map of species names to indices.
Definition: Phase.h:851
void checkSpeciesIndex(size_t k) const
Check that the specified species index is in range.
Definition: Phase.cpp:278
vector_fp m_speciesSize
Vector of species sizes.
Definition: Phase.h:798
size_t nSpecies() const
Returns the number of species in the phase.
Definition: Phase.h:262
virtual doublereal density() const
Density (kg/m^3).
Definition: Phase.h:607
doublereal m_mmw
mean molecular weight of the mixture (kg kmol-1)
Definition: Phase.h:826
Class Phase is the base class for phases of matter, managing the species and elements in a phase...
Definition: Phase.h:96
#define CT_ELEM_TYPE_ABSPOS
Normal element constraint consisting of positive coefficients for the formula matrix.
Definition: Elements.h:37
virtual bool ready() const
Returns a bool indicating whether the object is ready for use.
Definition: Phase.cpp:890
UndefElement::behavior m_undefinedElementBehavior
Flag determining behavior when adding species with an undefined element.
Definition: Phase.h:805
std::string m_name
Name of the phase.
Definition: Phase.h:817
const std::vector< std::string > & speciesNames() const
Return a const reference to the vector of species names.
Definition: Phase.cpp:273
doublereal mean_X(const doublereal *const Q) const
Evaluate the mole-fraction-weighted mean of an array Q.
Definition: Phase.cpp:690
int atomicNumber(size_t m) const
Atomic number of element m.
Definition: Phase.cpp:220
void setName(const std::string &nm)
Sets the string name for the phase.
Definition: Phase.cpp:156
std::vector< int > vector_int
Vector of ints.
Definition: ct_defs.h:159
doublereal elementalMoleFraction(const size_t m) const
Elemental mole fraction of element m.
Definition: Phase.cpp:648
void throwUndefinedElements()
Set the behavior when adding a species containing undefined elements to throw an exception.
Definition: Phase.cpp:886
int m_stateNum
State Change variable.
Definition: Phase.h:845
size_t m_ndim
Dimensionality of the phase.
Definition: Phase.h:788
doublereal atomicWeight(size_t m) const
Atomic weight of element m.
Definition: Phase.cpp:201
void ignoreUndefinedElements()
Set behavior when adding a species containing undefined elements to just skip the species...
Definition: Phase.cpp:878
std::vector< std::string > m_speciesNames
Vector of the species names.
Definition: Phase.h:848
void setState_TRX(doublereal t, doublereal dens, const doublereal *x)
Set the internally stored temperature (K), density, and mole fractions.
Definition: Phase.cpp:431
void getAtoms(size_t k, double *atomArray) const
Get a vector containing the atomic composition of species k.
Definition: Phase.cpp:244
void checkElementArraySize(size_t mm) const
Check that an array size is at least nElements().
Definition: Phase.cpp:173
doublereal molarDensity() const
Molar density (kmol/m^3).
Definition: Phase.cpp:666
shared_ptr< Species > species(const std::string &name) const
Return the Species object for the named species.
Definition: Phase.cpp:868
std::string speciesName(size_t k) const
Name of the species with index k.
Definition: Phase.cpp:267
virtual void compositionChanged()
Apply changes to the state which are needed after the composition changes.
Definition: Phase.cpp:899
doublereal entropyElement298(size_t m) const
Entropy of the element in its standard state at 298 K and 1 bar.
Definition: Phase.cpp:206
std::string speciesSPName(int k) const
Returns the expanded species name of a species, including the phase name This is guaranteed to be uni...
Definition: Phase.cpp:292
const doublereal * moleFractdivMMW() const
Returns a const pointer to the start of the moleFraction/MW array.
Definition: Phase.cpp:563
void setMoleFractionsByName(const compositionMap &xMap)
Set the species mole fractions by name.
Definition: Phase.cpp:368
ValueCache m_cache
Cached for saved calculations within each ThermoPhase.
Definition: Phase.h:761
void setMassFractionsByName(const compositionMap &yMap)
Set the species mass fractions by name.
Definition: Phase.cpp:412
virtual void modifySpecies(size_t k, shared_ptr< Species > spec)
Modify the thermodynamic data associated with a species.
Definition: Phase.cpp:851
Base class for exceptions thrown by Cantera classes.
Definition: ctexceptions.h:65
void setState_TNX(doublereal t, doublereal n, const doublereal *x)
Set the internally stored temperature (K), molar density (kmol/m^3), and mole fractions.
Definition: Phase.cpp:438
int stateMFNumber() const
Return the State Mole Fraction Number.
Definition: Phase.h:747
virtual void invalidateCache()
Invalidate any cached values which are normally updated only when a change in state is detected...
Definition: Phase.cpp:895
virtual void setConcentrations(const doublereal *const conc)
Set the concentrations to the specified values within the phase.
Definition: Phase.cpp:600
virtual void setMoleFractions(const doublereal *const x)
Set the mole fractions to the specified values.
Definition: Phase.cpp:327
doublereal concentration(const size_t k) const
Concentration of species k.
Definition: Phase.cpp:589
virtual bool addSpecies(shared_ptr< Species > spec)
Add a Species to this Phase.
Definition: Phase.cpp:761
void getMoleFractions(doublereal *const x) const
Get the species mole fraction vector.
Definition: Phase.cpp:542
virtual void setMolarDensity(const doublereal molarDensity)
Set the internally stored molar density (kmol/m^3) of the phase.
Definition: Phase.cpp:671
XML_Node & xml() const
Returns a const reference to the XML_Node that describes the phase.
Definition: Phase.cpp:117
virtual void setConcentrationsNoNorm(const double *const conc)
Set the concentrations without ignoring negative concentrations.
Definition: Phase.cpp:620
compositionMap getMassFractionsByName(double threshold=0.0) const
Get the mass fractions by name.
Definition: Phase.cpp:530
#define ENTROPY298_UNKNOWN
Number indicating we don&#39;t know the entropy of the element in its most stable state at 298...
Definition: Elements.h:87
virtual void setMassFractions_NoNorm(const doublereal *const y)
Set the mass fractions to the specified values without normalizing.
Definition: Phase.cpp:401
void setNDim(size_t ndim)
Set the number of spatial dimensions (1, 2, or 3).
Definition: Phase.h:592
void setState_TX(doublereal t, doublereal *x)
Set the internally stored temperature (K) and mole fractions.
Definition: Phase.cpp:472
vector_int m_elem_type
Vector of element types.
Definition: Phase.h:857
doublereal elementalMassFraction(const size_t m) const
Elemental mass fraction of element m.
Definition: Phase.cpp:637
void setMolecularWeight(const int k, const double mw)
Set the molecular weight of a single species to a given value.
Definition: Phase.h:770
doublereal chargeDensity() const
Charge density [C/m^3].
Definition: Phase.cpp:681
std::string id() const
Return the string id for the phase.
Definition: Phase.cpp:141
size_t addElement(const std::string &symbol, doublereal weight=-12345.0, int atomicNumber=0, doublereal entropy298=ENTROPY298_UNKNOWN, int elem_type=CT_ELEM_TYPE_ABSPOS)
Add an element.
Definition: Phase.cpp:705
void checkSpeciesArraySize(size_t kk) const
Check that an array size is at least nSpecies().
Definition: Phase.cpp:285
vector_fp m_molwts
species molecular weights (kg kmol-1)
Definition: Phase.h:839
void setState_TY(doublereal t, doublereal *y)
Set the internally stored temperature (K) and mass fractions.
Definition: Phase.cpp:478
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
doublereal meanMolecularWeight() const
The mean molecular weight. Units: (kg/kmol)
Definition: Phase.h:661
void addUndefinedElements()
Set behavior when adding a species containing undefined elements to add those elements to the phase...
Definition: Phase.cpp:882
size_t m_mm
Number of elements.
Definition: Phase.h:853
virtual void setTemperature(const doublereal temp)
Set the internally stored temperature of the phase (K).
Definition: Phase.h:637
std::string m_id
ID of the phase.
Definition: Phase.h:812
doublereal m_dens
Density (kg m-3).
Definition: Phase.h:824
vector_fp m_ym
m_ym[k] = mole fraction of species k divided by the mean molecular weight of mixture.
Definition: Phase.h:830
int elementType(size_t m) const
Return the element constraint type Possible types include:
Definition: Phase.cpp:225
void getMolecularWeights(vector_fp &weights) const
Copy the vector of molecular weights into vector weights.
Definition: Phase.cpp:502
XML_Node * m_xml
XML node containing the XML info for this phase.
Definition: Phase.h:808
std::string name() const
Return the name of the phase.
Definition: Phase.cpp:151
vector_int m_atomicNumbers
element atomic numbers
Definition: Phase.h:855
size_t nDim() const
Returns the number of spatial dimensions (1, 2, or 3)
Definition: Phase.h:585
doublereal size(size_t k) const
This routine returns the size of species k.
Definition: Phase.h:414
doublereal m_temp
Temperature (K). This is an independent variable.
Definition: Phase.h:819
virtual void setMassFractions(const doublereal *const y)
Set the mass fractions to the specified values and normalize them.
Definition: Phase.cpp:387
vector_fp m_rmolwts
inverse of species molecular weights (kmol kg-1)
Definition: Phase.h:841
vector_fp m_speciesComp
Atomic composition of the species.
Definition: Phase.h:793
void setXMLdata(XML_Node &xmlPhase)
Stores the XML tree information for the current phase.
Definition: Phase.cpp:122
void setState_RX(doublereal rho, doublereal *x)
Set the density (kg/m^3) and mole fractions.
Definition: Phase.cpp:484
size_t m_kk
Number of species in the phase.
Definition: Phase.h:784
int changeElementType(int m, int elem_type)
Change the element type of the mth constraint Reassigns an element type.
Definition: Phase.cpp:230
void setID(const std::string &id)
Set the string id for the phase.
Definition: Phase.cpp:146
Declaration for class Cantera::Species.
doublereal molecularWeight(size_t k) const
Molecular weight of species k.
Definition: Phase.cpp:496
Namespace for the Cantera kernel.
Definition: application.cpp:29
void getConcentrations(doublereal *const c) const
Get the species concentrations (kmol/m^3).
Definition: Phase.cpp:595
doublereal nAtoms(size_t k, size_t m) const
Number of atoms of element m in species k.
Definition: Phase.cpp:237
void setState_TR(doublereal t, doublereal rho)
Set the internally stored temperature (K) and density (kg/m^3)
Definition: Phase.cpp:466
const std::vector< std::string > & elementNames() const
Return a read-only reference to the vector of element names.
Definition: Phase.cpp:196
const vector_fp & atomicWeights() const
Return a read-only reference to the vector of atomic weights.
Definition: Phase.cpp:215
doublereal charge(size_t k) const
Dimensionless electrical charge of a single molecule of species k The charge is normalized by the the...
Definition: Phase.h:577
Definitions for the classes that are thrown when Cantera experiences an error condition (also contain...
std::string elementName(size_t m) const
Name of the element with index m.
Definition: Phase.cpp:180
void checkElementIndex(size_t m) const
Check that the specified element index is in range.
Definition: Phase.cpp:166
virtual void setDensity(const doublereal density_)
Set the internally stored density (kg/m^3) of the phase.
Definition: Phase.h:622
doublereal massFraction(size_t k) const
Return the mass fraction of a single species.
Definition: Phase.cpp:568