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