GRINS-0.6.0
Public Member Functions | Protected Attributes | Private Member Functions | List of all members
GRINS::CanteraKinetics Class Reference

#include <cantera_kinetics.h>

Public Member Functions

 CanteraKinetics (CanteraMixture &mixture)
 
 ~CanteraKinetics ()
 
void omega_dot (const CachedValues &cache, unsigned int qp, std::vector< libMesh::Real > &omega_dot) const
 
void omega_dot_TPY (const libMesh::Real T, const libMesh::Real P, const std::vector< libMesh::Real > &mass_fractions, std::vector< libMesh::Real > &omega_dot) const
 
void omega_dot_TRY (const libMesh::Real &T, const libMesh::Real rho, const std::vector< libMesh::Real > &mass_fractions, std::vector< libMesh::Real > &omega_dot) const
 

Protected Attributes

Cantera::IdealGasMix & _cantera_gas
 

Private Member Functions

 CanteraKinetics ()
 

Detailed Description

Definition at line 54 of file cantera_kinetics.h.

Constructor & Destructor Documentation

GRINS::CanteraKinetics::CanteraKinetics ( CanteraMixture mixture)

Definition at line 48 of file cantera_kinetics.C.

49  : _cantera_gas( mixture.get_chemistry() )
50  {
51  return;
52  }
Cantera::IdealGasMix & _cantera_gas
GRINS::CanteraKinetics::~CanteraKinetics ( )

Definition at line 54 of file cantera_kinetics.C.

55  {
56  return;
57  }
GRINS::CanteraKinetics::CanteraKinetics ( )
private

Member Function Documentation

void GRINS::CanteraKinetics::omega_dot ( const CachedValues cache,
unsigned int  qp,
std::vector< libMesh::Real > &  omega_dot 
) const

Definition at line 59 of file cantera_kinetics.C.

References GRINS::CachedValues::get_cached_values(), GRINS::CachedValues::get_cached_vector_values(), GRINS::Cache::MASS_FRACTIONS, omega_dot_TPY(), GRINS::Cache::TEMPERATURE, and GRINS::Cache::THERMO_PRESSURE.

Referenced by main(), and GRINS::CanteraEvaluator::omega_dot().

62  {
63  const libMesh::Real T = cache.get_cached_values(Cache::TEMPERATURE)[qp];
64  const libMesh::Real P = cache.get_cached_values(Cache::THERMO_PRESSURE)[qp];
65  const std::vector<libMesh::Real>& Y = cache.get_cached_vector_values(Cache::MASS_FRACTIONS)[qp];
66 
67  this->omega_dot_TPY( T, P, Y, omega_dot );
68 
69  return;
70  }
void omega_dot(const CachedValues &cache, unsigned int qp, std::vector< libMesh::Real > &omega_dot) const
void omega_dot_TPY(const libMesh::Real T, const libMesh::Real P, const std::vector< libMesh::Real > &mass_fractions, std::vector< libMesh::Real > &omega_dot) const
void GRINS::CanteraKinetics::omega_dot_TPY ( const libMesh::Real  T,
const libMesh::Real  P,
const std::vector< libMesh::Real > &  mass_fractions,
std::vector< libMesh::Real > &  omega_dot 
) const
Todo:
Need to make sure this will work in a threaded environment. Not sure if we will get thread lock here or not.

Definition at line 72 of file cantera_kinetics.C.

References _cantera_gas.

Referenced by main(), and omega_dot().

75  {
76  libmesh_assert_equal_to( mass_fractions.size(), omega_dot.size() );
77  libmesh_assert_equal_to( mass_fractions.size(), _cantera_gas.nSpecies() );
78  libmesh_assert_greater(T,0.0);
79  libmesh_assert_greater(P,0.0);
80 
81  {
84  libMesh::Threads::spin_mutex::scoped_lock lock(cantera_mutex);
85 
86  try
87  {
88  _cantera_gas.setState_TPY(T, P, &mass_fractions[0]);
89  _cantera_gas.getNetProductionRates(&omega_dot[0]);
90  }
91  catch(Cantera::CanteraError)
92  {
93  Cantera::showErrors(std::cerr);
94  libmesh_error();
95  }
96  }
97 
98 #ifdef DEBUG
99  for( unsigned int s = 0; s < omega_dot.size(); s++ )
100  {
101  if( libMesh::libmesh_isnan(omega_dot[s]) )
102  {
103  std::cout << "T = " << T << std::endl
104  << "P = " << P << std::endl;
105  for( unsigned int s = 0; s < omega_dot.size(); s++ )
106  {
107  std::cout << "Y[" << s << "] = " << mass_fractions[s] << std::endl;
108  }
109  for( unsigned int s = 0; s < omega_dot.size(); s++ )
110  {
111  std::cout << "omega_dot[" << s << "] = " << omega_dot[s] << std::endl;
112  }
113 
114  libmesh_error();
115  }
116  }
117 #endif
118 
119  for( unsigned int s = 0; s < omega_dot.size(); s++ )
120  {
121  // convert [kmol/m^3-s] to [kg/m^3-s]
122  omega_dot[s] *= this->_cantera_gas.molecularWeight(s);
123  }
124 
125  return;
126  }
void omega_dot(const CachedValues &cache, unsigned int qp, std::vector< libMesh::Real > &omega_dot) const
Cantera::IdealGasMix & _cantera_gas
void GRINS::CanteraKinetics::omega_dot_TRY ( const libMesh::Real &  T,
const libMesh::Real  rho,
const std::vector< libMesh::Real > &  mass_fractions,
std::vector< libMesh::Real > &  omega_dot 
) const
Todo:
Need to make sure this will work in a threaded environment. Not sure if we will get thread lock here or not.

Definition at line 128 of file cantera_kinetics.C.

References _cantera_gas.

Referenced by GRINS::CanteraEvaluator::omega_dot().

131  {
132  libmesh_assert_equal_to( mass_fractions.size(), omega_dot.size() );
133  libmesh_assert_equal_to( mass_fractions.size(), _cantera_gas.nSpecies() );
134  libmesh_assert_greater(T,0.0);
135  libmesh_assert_greater(rho,0.0);
136 
137  {
140  libMesh::Threads::spin_mutex::scoped_lock lock(cantera_mutex);
141 
142  try
143  {
144  _cantera_gas.setState_TRY(T, rho, &mass_fractions[0]);
145  _cantera_gas.getNetProductionRates(&omega_dot[0]);
146  }
147  catch(Cantera::CanteraError)
148  {
149  Cantera::showErrors(std::cerr);
150  libmesh_error();
151  }
152  }
153 
154  for( unsigned int s = 0; s < omega_dot.size(); s++ )
155  {
156  // convert [kmol/m^3-s] to [kg/m^3-s]
157  omega_dot[s] *= this->_cantera_gas.molecularWeight(s);
158  }
159 
160  return;
161  }
void omega_dot(const CachedValues &cache, unsigned int qp, std::vector< libMesh::Real > &omega_dot) const
Cantera::IdealGasMix & _cantera_gas

Member Data Documentation

Cantera::IdealGasMix& GRINS::CanteraKinetics::_cantera_gas
protected

Definition at line 74 of file cantera_kinetics.h.

Referenced by omega_dot_TPY(), and omega_dot_TRY().


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

Generated on Mon Jun 22 2015 21:32:22 for GRINS-0.6.0 by  doxygen 1.8.9.1