Cantera  2.2.1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes | List of all members

A species thermodynamic property manager for the NASA polynomial parameterization with two temperature ranges. More...

#include <NasaThermo.h>

Inheritance diagram for NasaThermo:
[legend]
Collaboration diagram for NasaThermo:
[legend]

Public Member Functions

 NasaThermo (const NasaThermo &right)
 
NasaThermooperator= (const NasaThermo &right)
 
virtual SpeciesThermoduplMyselfAsSpeciesThermo () const
 Duplication routine for objects derived from SpeciesThermo. More...
 
virtual void install (const std::string &name, size_t index, int type, const doublereal *c, doublereal min_temp, doublereal max_temp, doublereal ref_pressure)
 install a new species thermodynamic property parameterization for one species. More...
 
virtual void install_STIT (size_t index, shared_ptr< SpeciesThermoInterpType > stit_ptr)
 Install a new species thermodynamic property parameterization for one species. More...
 
virtual void update_one (size_t k, doublereal t, doublereal *cp_R, doublereal *h_RT, doublereal *s_R) const
 Like update(), but only updates the single species k. More...
 
virtual void update (doublereal t, doublereal *cp_R, doublereal *h_RT, doublereal *s_R) const
 Compute the reference-state properties for all species. More...
 
virtual doublereal minTemp (size_t k=npos) const
 Minimum temperature. More...
 
virtual doublereal maxTemp (size_t k=npos) const
 Maximum temperature. More...
 
virtual doublereal refPressure (size_t k=npos) const
 The reference-state pressure for species k. More...
 
virtual int reportType (size_t index) const
 This utility function reports the type of parameterization used for the species with index number index. More...
 
virtual void reportParams (size_t index, int &type, doublereal *const c, doublereal &minTemp, doublereal &maxTemp, doublereal &refPressure) const
 
virtual doublereal reportOneHf298 (const size_t k) const
 Report the 298 K Heat of Formation of the standard state of one species (J kmol-1) More...
 
virtual void modifyOneHf298 (const size_t k, const doublereal Hf298New)
 Modify the value of the 298 K Heat of Formation of the standard state of one species in the phase (J kmol-1) More...
 
- Public Member Functions inherited from SpeciesThermo
 SpeciesThermo ()
 Constructor. More...
 
virtual ~SpeciesThermo ()
 Destructor. More...
 
bool ready (size_t nSpecies)
 Check if data for all species (0 through nSpecies-1) has been installed. More...
 

Public Attributes

const int ID
 Initialized to the type of parameterization. More...
 

Protected Member Functions

doublereal cp_R (double t, const doublereal *c)
 Compute the nondimensional heat capacity using the given NASA polynomial. More...
 
doublereal enthalpy_RT (double t, const doublereal *c)
 Compute the nondimensional enthalpy using the given NASA polynomial. More...
 
doublereal entropy_R (double t, const doublereal *c)
 Compute the nondimensional entropy using the given NASA polynomial. More...
 
double checkContinuity (const std::string &name, double tmid, doublereal *clow, doublereal *chigh)
 Adjust polynomials to be continuous at the midpoint temperature. More...
 
void fixDiscontinuities (doublereal Tlow, doublereal Tmid, doublereal Thigh, doublereal *clow, doublereal *chigh)
 Adjust polynomials to be continuous at the midpoint temperature. More...
 
- Protected Member Functions inherited from SpeciesThermo
void markInstalled (size_t k)
 Mark species k as having its thermodynamic data installed. More...
 

Protected Attributes

std::vector< std::vector
< NasaPoly1 > > 
m_high
 Vector of vector of NasaPoly1's for the high temp region. More...
 
std::vector< std::vector
< NasaPoly1 > > 
m_low
 Vector of vector of NasaPoly1's for the low temp region. More...
 
std::map< int, int > m_index
 Map between the midpoint temperature, as an int, to the group number. More...
 
vector_fp m_tmid
 Vector of log temperature limits. More...
 
doublereal m_tlow_max
 Maximum value of the low temperature limit. More...
 
doublereal m_thigh_min
 Minimum value of the high temperature limit. More...
 
vector_fp m_tlow
 Vector of low temperature limits (species index) More...
 
vector_fp m_thigh
 Vector of low temperature limits (species index) More...
 
doublereal m_p0
 Reference pressure (Pa) More...
 
int m_ngroups
 number of groups More...
 
vector_fp m_t
 Vector of temperature polynomials. More...
 
std::map< size_t, size_t > m_group_map
 
std::map< size_t, size_t > m_posInGroup_map
 
std::map< size_t, std::string > m_name
 Species name as a function of the species index. More...
 

Detailed Description

A species thermodynamic property manager for the NASA polynomial parameterization with two temperature ranges.

This class is designed to efficiently evaluate the properties of a large number of species with the NASA parameterization.

The original NASA polynomial parameterization expressed the heat capacity as a fourth-order polynomial in temperature, with separate coefficients for each of two temperature ranges. (The newer NASA format adds coefficients for 1/T and 1/T^2, and allows multiple temperature ranges.) This class is designed for use with the original parameterization, which is used, for example, by the Chemkin software package.

In many cases, the midpoint temperature is the same for many species. To take advantage of this, class NasaThermo groups species with a common midpoint temperature, so that checking which range the desired temperature is in need be done only once for each group.

Note
There is a special CTML element for entering the coefficients of this parameterization.
See Also
importCTML
Deprecated:
To be removed after Cantera 2.2. Use GeneralSpeciesThermo instead.

Definition at line 46 of file NasaThermo.h.

Member Function Documentation

virtual SpeciesThermo* duplMyselfAsSpeciesThermo ( ) const
inlinevirtual

Duplication routine for objects derived from SpeciesThermo.

This function can be used to duplicate objects derived from SpeciesThermo even if the application only has a pointer to SpeciesThermo to work with.

Implements SpeciesThermo.

Definition at line 55 of file NasaThermo.h.

void install ( const std::string &  name,
size_t  index,
int  type,
const doublereal *  c,
doublereal  min_temp,
doublereal  max_temp,
doublereal  ref_pressure 
)
virtual

install a new species thermodynamic property parameterization for one species.

Parameters
nameName of the species
indexThe 'update' method will update the property values for this species at position i index in the property arrays.
typeint flag specifying the type of parameterization to be installed.
cvector of coefficients for the parameterization.
  • c[0] midpoint temperature
  • c[1] - c[7] coefficients for low T range
  • c[8] - c[14] coefficients for high T range
min_tempminimum temperature for which this parameterization is valid.
max_tempmaximum temperature for which this parameterization is valid.
ref_pressurestandard-state pressure for this parameterization.
See Also
speciesThermoTypes.h

Implements SpeciesThermo.

Definition at line 61 of file NasaThermo.cpp.

References NasaThermo::checkContinuity(), Cantera::fp2str(), Cantera::int2str(), NasaThermo::m_group_map, NasaThermo::m_high, NasaThermo::m_index, NasaThermo::m_low, NasaThermo::m_name, NasaThermo::m_ngroups, NasaThermo::m_p0, NasaThermo::m_posInGroup_map, NasaThermo::m_thigh, NasaThermo::m_thigh_min, NasaThermo::m_tlow, NasaThermo::m_tlow_max, NasaThermo::m_tmid, SpeciesThermo::markInstalled(), NASA, and Cantera::writelog().

virtual void install_STIT ( size_t  index,
shared_ptr< SpeciesThermoInterpType stit 
)
inlinevirtual

Install a new species thermodynamic property parameterization for one species.

Parameters
indexIndex of the species being installed
stitPointer to the SpeciesThermoInterpType object This will set up the thermo for one species

Implements SpeciesThermo.

Definition at line 85 of file NasaThermo.h.

void update_one ( size_t  k,
doublereal  t,
doublereal *  cp_R,
doublereal *  h_RT,
doublereal *  s_R 
) const
virtual

Like update(), but only updates the single species k.

Parameters
kspecies index
tTemperature (Kelvin)
cp_RVector of Dimensionless heat capacities. (length m_kk).
h_RTVector of Dimensionless enthalpies. (length m_kk).
s_RVector of Dimensionless entropies. (length m_kk).

Reimplemented from SpeciesThermo.

Definition at line 124 of file NasaThermo.cpp.

References Cantera::getValue(), NasaThermo::m_group_map, NasaThermo::m_high, NasaThermo::m_low, NasaThermo::m_posInGroup_map, NasaThermo::m_t, SpeciesThermoInterpType::maxTemp(), and NasaPoly1::updateProperties().

void update ( doublereal  T,
doublereal *  cp_R,
doublereal *  h_RT,
doublereal *  s_R 
) const
virtual

Compute the reference-state properties for all species.

Given temperature T in K, this method updates the values of the non- dimensional heat capacity at constant pressure, enthalpy, and entropy, at the reference pressure, Pref of each of the standard states.

Parameters
TTemperature (Kelvin)
cp_RVector of Dimensionless heat capacities. (length m_kk).
h_RTVector of Dimensionless enthalpies. (length m_kk).
s_RVector of Dimensionless entropies. (length m_kk).

Implements SpeciesThermo.

Definition at line 149 of file NasaThermo.cpp.

References NasaThermo::m_high, NasaThermo::m_low, NasaThermo::m_ngroups, NasaThermo::m_t, and NasaThermo::m_tmid.

virtual doublereal minTemp ( size_t  k = npos) const
inlinevirtual

Minimum temperature.

If no argument is supplied, this method returns the minimum temperature for which all parameterizations are valid. If an integer index k is supplied, then the value returned is the minimum temperature for species k in the phase.

Parameters
kSpecies index

Implements SpeciesThermo.

Definition at line 103 of file NasaThermo.h.

References NasaThermo::m_tlow, NasaThermo::m_tlow_max, and Cantera::npos.

virtual doublereal maxTemp ( size_t  k = npos) const
inlinevirtual

Maximum temperature.

If no argument is supplied, this method returns the maximum temperature for which all parameterizations are valid. If an integer index k is supplied, then the value returned is the maximum temperature for parameterization k.

Parameters
kSpecies Index

Implements SpeciesThermo.

Definition at line 111 of file NasaThermo.h.

References NasaThermo::m_thigh, NasaThermo::m_thigh_min, and Cantera::npos.

virtual doublereal refPressure ( size_t  k = npos) const
inlinevirtual

The reference-state pressure for species k.

Returns the reference state pressure in Pascals for species k. If k is left out of the argument list, it returns the reference state pressure for the first species. Note that some SpeciesThermo implementations, such as those for ideal gases, require that all species in the same phase have the same reference state pressures.

Parameters
kSpecies Index

Implements SpeciesThermo.

Definition at line 119 of file NasaThermo.h.

References NasaThermo::m_p0.

virtual int reportType ( size_t  index) const
inlinevirtual

This utility function reports the type of parameterization used for the species with index number index.

Parameters
indexSpecies index

Implements SpeciesThermo.

Definition at line 123 of file NasaThermo.h.

References NASA.

Referenced by NasaThermo::reportParams().

void reportParams ( size_t  index,
int &  type,
doublereal *const  c,
doublereal &  minTemp,
doublereal &  maxTemp,
doublereal &  refPressure 
) const
virtual

This utility function reports back the type of parameterization and all of the parameters for the species, index.

Parameters
indexSpecies index
typeInteger type of the standard type
cVector of coefficients used to set the parameters for the standard state. For the NASA object, there are 15 coefficients.
minTempoutput - Minimum temperature
maxTempoutput - Maximum temperature
refPressureoutput - reference pressure (Pa).

Implements SpeciesThermo.

Definition at line 176 of file NasaThermo.cpp.

References Cantera::getValue(), NasaThermo::m_group_map, NasaThermo::m_high, NasaThermo::m_low, NasaThermo::m_posInGroup_map, SpeciesThermoInterpType::maxTemp(), NASA, NASA1, NasaPoly1::reportParameters(), and NasaThermo::reportType().

doublereal reportOneHf298 ( const size_t  k) const
virtual

Report the 298 K Heat of Formation of the standard state of one species (J kmol-1)

The 298K Heat of Formation is defined as the enthalpy change to create the standard state of the species from its constituent elements in their standard states at 298 K and 1 bar.

Parameters
kspecies index
Returns
Returns the current value of the Heat of Formation at 298K and 1 bar

Implements SpeciesThermo.

Definition at line 218 of file NasaThermo.cpp.

References Cantera::getValue(), NasaThermo::m_group_map, NasaThermo::m_high, NasaThermo::m_low, NasaThermo::m_posInGroup_map, SpeciesThermoInterpType::maxTemp(), and NasaPoly1::reportHf298().

Referenced by NasaThermo::modifyOneHf298().

void modifyOneHf298 ( const size_t  k,
const doublereal  Hf298New 
)
virtual

Modify the value of the 298 K Heat of Formation of the standard state of one species in the phase (J kmol-1)

The 298K heat of formation is defined as the enthalpy change to create the standard state of the species from its constituent elements in their standard states at 298 K and 1 bar.

Parameters
kIndex of the species
Hf298NewSpecify the new value of the Heat of Formation at 298K and 1 bar. units = J/kmol.

Implements SpeciesThermo.

Definition at line 234 of file NasaThermo.cpp.

References Cantera::getValue(), NasaThermo::m_group_map, NasaThermo::m_high, NasaThermo::m_low, NasaThermo::m_posInGroup_map, SpeciesThermoInterpType::maxTemp(), NasaPoly1::modifyOneHf298(), NasaPoly1::reportHf298(), and NasaThermo::reportOneHf298().

doublereal cp_R ( double  t,
const doublereal *  c 
)
protected

Compute the nondimensional heat capacity using the given NASA polynomial.

Parameters
ttemperature
ccoefficient array

Definition at line 259 of file NasaThermo.cpp.

References Cantera::poly4().

Referenced by NasaThermo::checkContinuity(), and NasaThermo::fixDiscontinuities().

doublereal enthalpy_RT ( double  t,
const doublereal *  c 
)
protected

Compute the nondimensional enthalpy using the given NASA polynomial.

Parameters
ttemperature
ccoefficient array

Definition at line 264 of file NasaThermo.cpp.

Referenced by NasaThermo::checkContinuity(), and NasaThermo::fixDiscontinuities().

doublereal entropy_R ( double  t,
const doublereal *  c 
)
protected

Compute the nondimensional entropy using the given NASA polynomial.

Parameters
ttemperature
ccoefficient array

Definition at line 270 of file NasaThermo.cpp.

Referenced by NasaThermo::checkContinuity(), and NasaThermo::fixDiscontinuities().

doublereal checkContinuity ( const std::string &  name,
double  tmid,
doublereal *  clow,
doublereal *  chigh 
)
protected

Adjust polynomials to be continuous at the midpoint temperature.

Check to see if the provided coefficients are nearly continuous. Adjust the values to get more precise continuity to avoid convergence issues with algorithms that expect these quantities to be continuous.

Parameters
namestring name of species
tmidMid temperature, between the two temperature regions
clowcoefficients for lower temperature region
chighcoefficients for higher temperature region

Definition at line 276 of file NasaThermo.cpp.

References NasaThermo::cp_R(), NasaThermo::enthalpy_RT(), NasaThermo::entropy_R(), Cantera::fp2str(), and Cantera::writelog().

Referenced by NasaThermo::install().

void fixDiscontinuities ( doublereal  Tlow,
doublereal  Tmid,
doublereal  Thigh,
doublereal *  clow,
doublereal *  chigh 
)
protected

Adjust polynomials to be continuous at the midpoint temperature.

We seek a set of coefficients for the low- and high-temperature polynomials which are continuous in Cp, H, and S at the midpoint while minimizing the difference between the values in Cp, H, and S over the entire valid temperature range. To do this, we formulate a linear least-squares problem to be solved for 11 of the 14 coefficients, with the remaining 3 coefficients eliminated in the process of satisfying the continuity constraints.

Parameters
TlowMinimum temperature at which the low-T polynomial is valid
TmidMid temperature, between the two temperature regions
ThighMaximum temperature at which the high-T polynomial is valid
clowcoefficients for lower temperature region
chighcoefficients for higher temperature region

Definition at line 327 of file NasaThermo.cpp.

References AssertTrace, NasaThermo::cp_R(), NasaThermo::enthalpy_RT(), and NasaThermo::entropy_R().

Member Data Documentation

const int ID

Initialized to the type of parameterization.

Note, this value is used in some template functions

Definition at line 154 of file NasaThermo.h.

std::vector<std::vector<NasaPoly1> > m_high
protected

Vector of vector of NasaPoly1's for the high temp region.

This is the high temp region representation. The first Length is equal to the number of groups. The second vector is equal to the number of species in that particular group.

Definition at line 164 of file NasaThermo.h.

Referenced by NasaThermo::install(), NasaThermo::modifyOneHf298(), NasaThermo::reportOneHf298(), NasaThermo::reportParams(), NasaThermo::update(), and NasaThermo::update_one().

std::vector<std::vector<NasaPoly1> > m_low
protected

Vector of vector of NasaPoly1's for the low temp region.

This is the low temp region representation. The first Length is equal to the number of groups. The second vector is equal to the number of species in that particular group.

Definition at line 173 of file NasaThermo.h.

Referenced by NasaThermo::install(), NasaThermo::modifyOneHf298(), NasaThermo::reportOneHf298(), NasaThermo::reportParams(), NasaThermo::update(), and NasaThermo::update_one().

std::map<int, int> m_index
protected

Map between the midpoint temperature, as an int, to the group number.

Length is equal to the number of groups. Only used in the setup.

Definition at line 179 of file NasaThermo.h.

Referenced by NasaThermo::install().

vector_fp m_tmid
protected

Vector of log temperature limits.

Length is equal to the number of groups.

Definition at line 185 of file NasaThermo.h.

Referenced by NasaThermo::install(), and NasaThermo::update().

doublereal m_tlow_max
protected

Maximum value of the low temperature limit.

Definition at line 188 of file NasaThermo.h.

Referenced by NasaThermo::install(), and NasaThermo::minTemp().

doublereal m_thigh_min
protected

Minimum value of the high temperature limit.

Definition at line 191 of file NasaThermo.h.

Referenced by NasaThermo::install(), and NasaThermo::maxTemp().

vector_fp m_tlow
protected

Vector of low temperature limits (species index)

Length is equal to number of species

Definition at line 197 of file NasaThermo.h.

Referenced by NasaThermo::install(), and NasaThermo::minTemp().

vector_fp m_thigh
protected

Vector of low temperature limits (species index)

Length is equal to number of species

Definition at line 203 of file NasaThermo.h.

Referenced by NasaThermo::install(), and NasaThermo::maxTemp().

doublereal m_p0
protected

Reference pressure (Pa)

all species must have the same reference pressure.

Definition at line 209 of file NasaThermo.h.

Referenced by NasaThermo::install(), and NasaThermo::refPressure().

int m_ngroups
protected

number of groups

Definition at line 212 of file NasaThermo.h.

Referenced by NasaThermo::install(), and NasaThermo::update().

vector_fp m_t
mutableprotected

Vector of temperature polynomials.

Definition at line 215 of file NasaThermo.h.

Referenced by NasaThermo::update(), and NasaThermo::update_one().

std::map<size_t, size_t> m_group_map
protected

This map takes as its index, the species index in the phase. It returns the group index, where the temperature polynomials for that species are stored. group indices start at 1, so a decrement is always performed to access vectors.

Definition at line 223 of file NasaThermo.h.

Referenced by NasaThermo::install(), NasaThermo::modifyOneHf298(), NasaThermo::reportOneHf298(), NasaThermo::reportParams(), and NasaThermo::update_one().

std::map<size_t, size_t> m_posInGroup_map
protected

This map takes as its index, the species index in the phase. It returns the position index within the group, where the temperature polynomials for that species are stored.

Definition at line 230 of file NasaThermo.h.

Referenced by NasaThermo::install(), NasaThermo::modifyOneHf298(), NasaThermo::reportOneHf298(), NasaThermo::reportParams(), and NasaThermo::update_one().

std::map<size_t, std::string> m_name
protected

Species name as a function of the species index.

Definition at line 233 of file NasaThermo.h.

Referenced by NasaThermo::install().


The documentation for this class was generated from the following files: