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