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