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