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