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

Wrapper class for evaluating thermo properties using Cantera. More...

#include <cantera_thermo.h>

Collaboration diagram for GRINS::CanteraThermodynamics:
Collaboration graph
[legend]

Public Member Functions

 CanteraThermodynamics (CanteraMixture &mixture)
 
 ~CanteraThermodynamics ()
 
libMesh::Real cp (const CachedValues &cache, unsigned int qp) const
 
libMesh::Real cv (const CachedValues &cache, unsigned int qp) const
 
libMesh::Real h (const CachedValues &cache, unsigned int qp, unsigned int species) const
 
void h (const CachedValues &cache, unsigned int qp, std::vector< libMesh::Real > &h) const
 
libMesh::Real h (const libMesh::Real &T, unsigned int species) const
 

Protected Attributes

CanteraMixture_cantera_mixture
 
Cantera::IdealGasMix & _cantera_gas
 

Private Member Functions

 CanteraThermodynamics ()
 

Detailed Description

Wrapper class for evaluating thermo properties using Cantera.

This class is expected to be constructed after threads have been forked and will only live during the lifetime of the thread. Note that this documentation will always be built regardless if Cantera is included in the GRINS build or not. Check configure output to confirm that Cantera was included in the build if you wish to use it.

Definition at line 61 of file cantera_thermo.h.

Constructor & Destructor Documentation

GRINS::CanteraThermodynamics::CanteraThermodynamics ( CanteraMixture mixture)

Definition at line 43 of file cantera_thermo.C.

44  : _cantera_mixture(mixture),
45  _cantera_gas(mixture.get_chemistry())
46  {
47  return;
48  }
Cantera::IdealGasMix & _cantera_gas
CanteraMixture & _cantera_mixture
GRINS::CanteraThermodynamics::~CanteraThermodynamics ( )

Definition at line 50 of file cantera_thermo.C.

51  {
52  return;
53  }
GRINS::CanteraThermodynamics::CanteraThermodynamics ( )
private

Member Function Documentation

libMesh::Real GRINS::CanteraThermodynamics::cp ( const CachedValues cache,
unsigned int  qp 
) 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 55 of file cantera_thermo.C.

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

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

56  {
57  const libMesh::Real T = cache.get_cached_values(Cache::TEMPERATURE)[qp];
58  const libMesh::Real P = cache.get_cached_values(Cache::THERMO_PRESSURE)[qp];
59  const std::vector<libMesh::Real>& Y = cache.get_cached_vector_values(Cache::MASS_FRACTIONS)[qp];
60 
61  libmesh_assert_equal_to( Y.size(), _cantera_gas.nSpecies() );
62 
63  libMesh::Real cp = 0.0;
64 
65  {
68  libMesh::Threads::spin_mutex::scoped_lock lock(cantera_mutex);
69 
70  try
71  {
72  _cantera_gas.setState_TPY( T, P, &Y[0] );
73 
74  cp = _cantera_gas.cp_mass();
75  }
76  catch(Cantera::CanteraError)
77  {
78  Cantera::showErrors(std::cerr);
79  libmesh_error();
80  }
81 
82  }
83 
84  return cp;
85  }
libMesh::Real cp(const CachedValues &cache, unsigned int qp) const
Cantera::IdealGasMix & _cantera_gas
libMesh::Real GRINS::CanteraThermodynamics::cv ( const CachedValues cache,
unsigned int  qp 
) 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 87 of file cantera_thermo.C.

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

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

88  {
89  const libMesh::Real T = cache.get_cached_values(Cache::TEMPERATURE)[qp];
90  const libMesh::Real P = cache.get_cached_values(Cache::THERMO_PRESSURE)[qp];
91  const std::vector<libMesh::Real>& Y = cache.get_cached_vector_values(Cache::MASS_FRACTIONS)[qp];
92 
93  libmesh_assert_equal_to( Y.size(), _cantera_gas.nSpecies() );
94 
95  libMesh::Real cv = 0.0;
96 
97  {
100  libMesh::Threads::spin_mutex::scoped_lock lock(cantera_mutex);
101 
102  try
103  {
104  _cantera_gas.setState_TPY( T, P, &Y[0] );
105 
106  cv = _cantera_gas.cv_mass();
107  }
108  catch(Cantera::CanteraError)
109  {
110  Cantera::showErrors(std::cerr);
111  libmesh_error();
112  }
113 
114  }
115 
116  return cv;
117  }
Cantera::IdealGasMix & _cantera_gas
libMesh::Real cv(const CachedValues &cache, unsigned int qp) const
libMesh::Real GRINS::CanteraThermodynamics::h ( const CachedValues cache,
unsigned int  qp,
unsigned int  species 
) 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 119 of file cantera_thermo.C.

References _cantera_gas, _cantera_mixture, GRINS::CachedValues::get_cached_values(), GRINS::CachedValues::get_cached_vector_values(), GRINS::Cache::MASS_FRACTIONS, GRINS::CanteraMixture::R(), GRINS::Cache::TEMPERATURE, and GRINS::Cache::THERMO_PRESSURE.

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

121  {
122  const libMesh::Real T = cache.get_cached_values(Cache::TEMPERATURE)[qp];
123  const libMesh::Real P = cache.get_cached_values(Cache::THERMO_PRESSURE)[qp];
124  const std::vector<libMesh::Real>& Y = cache.get_cached_vector_values(Cache::MASS_FRACTIONS)[qp];
125 
126  libmesh_assert_equal_to( Y.size(), _cantera_gas.nSpecies() );
127 
128  std::vector<libMesh::Real> h_RT( Y.size(), 0.0 );
129 
130  {
133  libMesh::Threads::spin_mutex::scoped_lock lock(cantera_mutex);
134 
135  try
136  {
137  _cantera_gas.setState_TPY( T, P, &Y[0] );
138 
139  _cantera_gas.getEnthalpy_RT( &h_RT[0] );
140  }
141  catch(Cantera::CanteraError)
142  {
143  Cantera::showErrors(std::cerr);
144  libmesh_error();
145  }
146 
147  }
148 
149  return h_RT[species]*_cantera_mixture.R(species)*T;
150  }
libMesh::Real R(unsigned int species) const
Cantera::IdealGasMix & _cantera_gas
CanteraMixture & _cantera_mixture
void GRINS::CanteraThermodynamics::h ( const CachedValues cache,
unsigned int  qp,
std::vector< libMesh::Real > &  h 
) 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 152 of file cantera_thermo.C.

References _cantera_gas, _cantera_mixture, GRINS::CachedValues::get_cached_values(), GRINS::CachedValues::get_cached_vector_values(), GRINS::Cache::MASS_FRACTIONS, GRINS::CanteraMixture::R(), GRINS::Cache::TEMPERATURE, and GRINS::Cache::THERMO_PRESSURE.

154  {
155  const libMesh::Real T = cache.get_cached_values(Cache::TEMPERATURE)[qp];
156  const libMesh::Real P = cache.get_cached_values(Cache::THERMO_PRESSURE)[qp];
157  const std::vector<libMesh::Real>& Y = cache.get_cached_vector_values(Cache::MASS_FRACTIONS)[qp];
158 
159  libmesh_assert_equal_to( Y.size(), h.size() );
160  libmesh_assert_equal_to( Y.size(), _cantera_gas.nSpecies() );
161 
162  {
165  libMesh::Threads::spin_mutex::scoped_lock lock(cantera_mutex);
166 
167  try
168  {
169  _cantera_gas.setState_TPY( T, P, &Y[0] );
170 
171  _cantera_gas.getEnthalpy_RT( &h[0] );
172  }
173  catch(Cantera::CanteraError)
174  {
175  Cantera::showErrors(std::cerr);
176  libmesh_error();
177  }
178 
179  for( unsigned int s = 0; s < h.size(); s++ )
180  {
181  h[s] *= _cantera_mixture.R(s)*T;
182  }
183 
184  }
185 
186  return;
187  }
libMesh::Real R(unsigned int species) const
Cantera::IdealGasMix & _cantera_gas
libMesh::Real h(const CachedValues &cache, unsigned int qp, unsigned int species) const
CanteraMixture & _cantera_mixture
libMesh::Real GRINS::CanteraThermodynamics::h ( const libMesh::Real &  T,
unsigned int  species 
) const

Definition at line 189 of file cantera_thermo.C.

References _cantera_gas, _cantera_mixture, and GRINS::CanteraMixture::R().

190  {
191  std::vector<libMesh::Real> h_RT( _cantera_gas.nSpecies(), 0.0 );
192 
193  try
194  {
195  _cantera_gas.setTemperature( T );
196 
197  _cantera_gas.getEnthalpy_RT( &h_RT[0] );
198  }
199  catch(Cantera::CanteraError)
200  {
201  Cantera::showErrors(std::cerr);
202  libmesh_error();
203  }
204 
205  return h_RT[species]*_cantera_mixture.R(species)*T;
206  }
libMesh::Real R(unsigned int species) const
Cantera::IdealGasMix & _cantera_gas
CanteraMixture & _cantera_mixture

Member Data Documentation

Cantera::IdealGasMix& GRINS::CanteraThermodynamics::_cantera_gas
protected

Definition at line 82 of file cantera_thermo.h.

Referenced by cp(), cv(), and h().

CanteraMixture& GRINS::CanteraThermodynamics::_cantera_mixture
protected

Definition at line 80 of file cantera_thermo.h.

Referenced by h().


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