antioch-0.4.0
List of all members | Public Member Functions | Protected Attributes
Antioch::KineticsEvaluator< CoeffType, StateType > Class Template Reference

Class to handle computing mass source terms for a given ReactionSet. More...

#include <kinetics_evaluator.h>

Public Member Functions

 KineticsEvaluator (const ReactionSet< CoeffType > &reaction_set, const StateType &example)
 Constructor. Requires a reaction set to be evaluated later,. More...
 
 ~KineticsEvaluator ()
 
const ReactionSet< CoeffType > & reaction_set () const
 
template<typename VectorStateType , typename KC >
void compute_mass_sources (const KC &conditions, const VectorStateType &molar_densities, const VectorStateType &h_RT_minus_s_R, VectorStateType &mass_sources)
 Compute species production/destruction rates per unit volume. More...
 
template<typename VectorStateType , typename KC >
void compute_mass_sources_and_derivs (const KC &conditions, const VectorStateType &molar_densities, const VectorStateType &h_RT_minus_s_R, const VectorStateType &dh_RT_minus_s_R_dT, VectorStateType &mass_sources, VectorStateType &dmass_dT, std::vector< VectorStateType > &dmass_drho_s)
 Compute species production/destruction rate derivatives. More...
 
template<typename VectorStateType , typename KC >
void compute_mole_sources (const KC &conditions, const VectorStateType &molar_densities, const VectorStateType &h_RT_minus_s_R, VectorStateType &mole_sources)
 Compute species molar production/destruction rates per unit volume. More...
 
template<typename VectorStateType , typename KC >
void compute_mole_sources_and_derivs (const KC &conditions, const VectorStateType &molar_densities, const VectorStateType &h_RT_minus_s_R, const VectorStateType &dh_RT_minus_s_R_dT, VectorStateType &mole_sources, VectorStateType &dmole_dT, std::vector< VectorStateType > &dmole_dX_s)
 Compute species production/destruction rate derivatives. More...
 
unsigned int n_species () const
 
unsigned int n_reactions () const
 

Protected Attributes

const ReactionSet< CoeffType > & _reaction_set
 
const ChemicalMixture
< CoeffType > & 
_chem_mixture
 
std::vector< StateType > _net_reaction_rates
 
std::vector< StateType > _dnet_rate_dT
 
std::vector< std::vector
< StateType > > 
_dnet_rate_dX_s
 

Detailed Description

template<typename CoeffType = double, typename StateType = CoeffType>
class Antioch::KineticsEvaluator< CoeffType, StateType >

Class to handle computing mass source terms for a given ReactionSet.

This class preallocates work arrays and so must be created within a spawned thread, if running in a threaded environment. It takes a reference to an already created ReactionSet, so there's little construction penalty.

Definition at line 53 of file kinetics_evaluator.h.

Constructor & Destructor Documentation

template<typename CoeffType , typename StateType >
Antioch::KineticsEvaluator< CoeffType, StateType >::KineticsEvaluator ( const ReactionSet< CoeffType > &  reaction_set,
const StateType &  example 
)
inline

Constructor. Requires a reaction set to be evaluated later,.

Definition at line 152 of file kinetics_evaluator.h.

References Antioch::ReactionSet< CoeffType >::n_species().

155  _chem_mixture( reaction_set.chemical_mixture() ),
156  _net_reaction_rates( reaction_set.n_reactions(), example ),
157  _dnet_rate_dT( reaction_set.n_reactions(), example ),
158  _dnet_rate_dX_s( reaction_set.n_reactions() )
159  {
160 
161  for( unsigned int r = 0; r < this->n_reactions(); r++ )
162  {
163  _dnet_rate_dX_s[r].resize( reaction_set.n_species(), example );
164  }
165 
166  return;
167  }
unsigned int n_reactions() const
const ChemicalMixture< CoeffType > & _chem_mixture
std::vector< StateType > _dnet_rate_dT
std::vector< StateType > _net_reaction_rates
const ReactionSet< CoeffType > & reaction_set() const
const ReactionSet< CoeffType > & _reaction_set
std::vector< std::vector< StateType > > _dnet_rate_dX_s
template<typename CoeffType , typename StateType >
Antioch::KineticsEvaluator< CoeffType, StateType >::~KineticsEvaluator ( )
inline

Definition at line 172 of file kinetics_evaluator.h.

173  {
174  return;
175  }

Member Function Documentation

template<typename CoeffType , typename StateType >
template<typename VectorStateType , typename KC >
void Antioch::KineticsEvaluator< CoeffType, StateType >::compute_mass_sources ( const KC &  conditions,
const VectorStateType &  molar_densities,
const VectorStateType &  h_RT_minus_s_R,
VectorStateType &  mass_sources 
)
inline

Compute species production/destruction rates per unit volume.

$ \left(kg/sec/m^3\right)$

Definition at line 234 of file kinetics_evaluator.h.

Referenced by tester(), and vectester().

238  {
239  // Quantities asserted in compute_mole_sources call
240  this->compute_mole_sources( conditions, molar_densities, h_RT_minus_s_R, mass_sources );
241 
242  // finally scale by molar mass
243  for (unsigned int s=0; s < this->n_species(); s++)
244  {
245  mass_sources[s] *= _chem_mixture.M(s);
246  }
247 
248  return;
249  }
const ChemicalMixture< CoeffType > & _chem_mixture
void compute_mole_sources(const KC &conditions, const VectorStateType &molar_densities, const VectorStateType &h_RT_minus_s_R, VectorStateType &mole_sources)
Compute species molar production/destruction rates per unit volume.
unsigned int n_species() const
template<typename CoeffType , typename StateType >
template<typename VectorStateType , typename KC >
void Antioch::KineticsEvaluator< CoeffType, StateType >::compute_mass_sources_and_derivs ( const KC &  conditions,
const VectorStateType &  molar_densities,
const VectorStateType &  h_RT_minus_s_R,
const VectorStateType &  dh_RT_minus_s_R_dT,
VectorStateType &  mass_sources,
VectorStateType &  dmass_dT,
std::vector< VectorStateType > &  dmass_drho_s 
)
inline

Compute species production/destruction rate derivatives.

In mass units, e.g. $ \frac{\partial \dot{\omega}}{dT} [\left(kg/sec/m^3/K\right)]$

Definition at line 358 of file kinetics_evaluator.h.

Referenced by tester(), and vectester().

365  {
366  // Asserts are in compute_mole_sources
367  this->compute_mole_sources_and_derivs( conditions, molar_densities, h_RT_minus_s_R, dh_RT_minus_s_R_dT,
368  mass_sources, dmass_dT, dmass_drho_s );
369 
370  // Convert from mole units to mass units
371  for (unsigned int s=0; s < this->n_species(); s++)
372  {
373  mass_sources[s] *= _chem_mixture.M(s);
374  dmass_dT[s] *= _chem_mixture.M(s);
375 
376  for (unsigned int t=0; t < this->n_species(); t++)
377  {
378  dmass_drho_s[s][t] *= _chem_mixture.M(s)/_chem_mixture.M(t);
379  }
380 
381  }
382 
383  return;
384  }
const ChemicalMixture< CoeffType > & _chem_mixture
void compute_mole_sources_and_derivs(const KC &conditions, const VectorStateType &molar_densities, const VectorStateType &h_RT_minus_s_R, const VectorStateType &dh_RT_minus_s_R_dT, VectorStateType &mole_sources, VectorStateType &dmole_dT, std::vector< VectorStateType > &dmole_dX_s)
Compute species production/destruction rate derivatives.
unsigned int n_species() const
template<typename CoeffType , typename StateType >
template<typename VectorStateType , typename KC >
void Antioch::KineticsEvaluator< CoeffType, StateType >::compute_mole_sources ( const KC &  conditions,
const VectorStateType &  molar_densities,
const VectorStateType &  h_RT_minus_s_R,
VectorStateType &  mole_sources 
)
inline

Compute species molar production/destruction rates per unit volume.

$ \left(mole/sec/m^3\right)$

Todo:
Make these assertions vector-compatible
Todo:
Do we need to really initialize this?

Definition at line 180 of file kinetics_evaluator.h.

References antioch_assert_equal_to, Antioch::Reaction< CoeffType, VectorCoeffType >::n_products(), Antioch::Reaction< CoeffType, VectorCoeffType >::n_reactants(), Antioch::Reaction< CoeffType, VectorCoeffType >::product_id(), Antioch::Reaction< CoeffType, VectorCoeffType >::product_stoichiometric_coefficient(), Antioch::Reaction< CoeffType, VectorCoeffType >::reactant_id(), Antioch::Reaction< CoeffType, VectorCoeffType >::reactant_stoichiometric_coefficient(), and Antioch::set_zero().

184  {
186  // antioch_assert_greater(T, 0.0);
187  antioch_assert_equal_to( molar_densities.size(), this->n_species() );
188  antioch_assert_equal_to( h_RT_minus_s_R.size(), this->n_species() );
189  antioch_assert_equal_to( mole_sources.size(), this->n_species() );
190 
193 
194  Antioch::set_zero(mole_sources);
195 
196  typename constructor_or_reference<const KineticsConditions<StateType,VectorStateType>, const KC>::type //either (KineticsConditions<> &) or (KineticsConditions<>)
197  kinetics_conditions(conditions);
198  // compute the requisite reaction rates
199  this->_reaction_set.compute_reaction_rates( kinetics_conditions, molar_densities,
200  h_RT_minus_s_R, _net_reaction_rates );
201 
202  // compute the actual mole sources in kmol/sec/m^3
203  for (unsigned int rxn = 0; rxn < this->n_reactions(); rxn++)
204  {
205  const Reaction<CoeffType>& reaction = this->_reaction_set.reaction(rxn);
206  const StateType rate = _net_reaction_rates[rxn];
207 
208  // reactant contributions
209  for (unsigned int r = 0; r < reaction.n_reactants(); r++)
210  {
211  const unsigned int r_id = reaction.reactant_id(r);
212  const unsigned int r_stoich = reaction.reactant_stoichiometric_coefficient(r);
213 
214  mole_sources[r_id] -= (static_cast<CoeffType>(r_stoich)*rate);
215  }
216 
217  // product contributions
218  for (unsigned int p=0; p < reaction.n_products(); p++)
219  {
220  const unsigned int p_id = reaction.product_id(p);
221  const unsigned int p_stoich = reaction.product_stoichiometric_coefficient(p);
222 
223  mole_sources[p_id] += (static_cast<CoeffType>(p_stoich)*rate);
224  }
225  }
226 
227  return;
228  }
unsigned int n_reactions() const
#define antioch_assert_equal_to(expr1, expr2)
unsigned int n_species() const
std::vector< StateType > _net_reaction_rates
const ReactionSet< CoeffType > & _reaction_set
void set_zero(_Matrix< _Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols > &a)
Definition: eigen_utils.h:217
template<typename CoeffType , typename StateType >
template<typename VectorStateType , typename KC >
void Antioch::KineticsEvaluator< CoeffType, StateType >::compute_mole_sources_and_derivs ( const KC &  conditions,
const VectorStateType &  molar_densities,
const VectorStateType &  h_RT_minus_s_R,
const VectorStateType &  dh_RT_minus_s_R_dT,
VectorStateType &  mole_sources,
VectorStateType &  dmole_dT,
std::vector< VectorStateType > &  dmole_dX_s 
)
inline

Compute species production/destruction rate derivatives.

In mass units, e.g. $ \frac{\partial \dot{\omega}}{dT} [\left(mole/sec/m^3/K\right)]$

Todo:
Make these assertions vector-compatible
Todo:
Do we need to really initialize these?
Todo:
Do we need to really initialize this?
Todo:
Are these going to get optimized out? Should we remove them?

Definition at line 254 of file kinetics_evaluator.h.

References antioch_assert_equal_to, Antioch::Reaction< CoeffType, VectorCoeffType >::n_products(), Antioch::Reaction< CoeffType, VectorCoeffType >::n_reactants(), Antioch::Reaction< CoeffType, VectorCoeffType >::product_id(), Antioch::Reaction< CoeffType, VectorCoeffType >::product_stoichiometric_coefficient(), Antioch::Reaction< CoeffType, VectorCoeffType >::reactant_id(), Antioch::Reaction< CoeffType, VectorCoeffType >::reactant_stoichiometric_coefficient(), and Antioch::set_zero().

261  {
263  // antioch_assert_greater(T, 0.0);
264  antioch_assert_equal_to( molar_densities.size(), this->n_species() );
265  antioch_assert_equal_to( h_RT_minus_s_R.size(), this->n_species() );
266  antioch_assert_equal_to( dh_RT_minus_s_R_dT.size(), this->n_species() );
267  antioch_assert_equal_to( mole_sources.size(), this->n_species() );
268  antioch_assert_equal_to( dmole_dT.size(), this->n_species() );
269  antioch_assert_equal_to( dmole_dX_s.size(), this->n_species() );
270 #ifdef NDEBUG
271 #else
272  for (unsigned int s=0; s < this->n_species(); s++)
273  {
274  antioch_assert_equal_to( dmole_dX_s[s].size(), this->n_species() );
275  }
276 #endif
277 
281 
282  Antioch::set_zero(mole_sources);
283  Antioch::set_zero(dmole_dT);
284  for (unsigned int s=0; s < this->n_species(); s++)
285  {
286  Antioch::set_zero(dmole_dX_s[s]);
287  }
288 
289  for (unsigned int rxn=0; rxn < this->n_reactions(); rxn++)
290  {
293  }
294 
295  typename constructor_or_reference<const KineticsConditions<StateType,VectorStateType>, const KC>::type //either (KineticsConditions<> &) or (KineticsConditions<>)
296  kinetics_conditions(conditions);
297  // compute the requisite reaction rates
298  this->_reaction_set.compute_reaction_rates_and_derivs( kinetics_conditions, molar_densities,
299  h_RT_minus_s_R, dh_RT_minus_s_R_dT,
301  _dnet_rate_dT,
302  _dnet_rate_dX_s );
303  // compute the actual mole sources in kmol/sec/m^3
304  for (unsigned int rxn = 0; rxn < this->n_reactions(); rxn++)
305  {
306 
307  const Reaction<CoeffType>& reaction = this->_reaction_set.reaction(rxn);
308 
310  const StateType rate = _net_reaction_rates[rxn];
311  const StateType drate_dT = _dnet_rate_dT[rxn];
312  const VectorStateType drate_dX_s = _dnet_rate_dX_s[rxn];
313 
314  // reactant contributions
315  for (unsigned int r = 0; r < reaction.n_reactants(); r++)
316  {
317  const unsigned int r_id = reaction.reactant_id(r);
318  const unsigned int r_stoich = reaction.reactant_stoichiometric_coefficient(r);
319 
320  mole_sources[r_id] -= (static_cast<CoeffType>(r_stoich)*rate);
321 
322  // d/dT rate contributions
323  dmole_dT[r_id] -= (static_cast<CoeffType>(r_stoich)*drate_dT);
324 
325  // d(.m)/dX_s rate contributions
326  for (unsigned int s=0; s < this->n_species(); s++)
327  {
328  dmole_dX_s[r_id][s] -= (static_cast<CoeffType>(r_stoich)*drate_dX_s[s]);
329  }
330  }
331 
332  // product contributions
333  for (unsigned int p=0; p < reaction.n_products(); p++)
334  {
335  const unsigned int p_id = reaction.product_id(p);
336  const unsigned int p_stoich = reaction.product_stoichiometric_coefficient(p);
337 
338  mole_sources[p_id] += (static_cast<CoeffType>(p_stoich)*rate);
339 
340  // d/dT rate contributions
341  dmole_dT[p_id] += (static_cast<CoeffType>(p_stoich)*drate_dT);
342 
343  // d/dX_s rate contributions
344  for (unsigned int s=0; s < this->n_species(); s++)
345  {
346  dmole_dX_s[p_id][s] += (static_cast<CoeffType>(p_stoich)*drate_dX_s[s]);
347  }
348  }
349  }
350 
351  return;
352  }
unsigned int n_reactions() const
#define antioch_assert_equal_to(expr1, expr2)
unsigned int n_species() const
std::vector< StateType > _dnet_rate_dT
std::vector< StateType > _net_reaction_rates
const ReactionSet< CoeffType > & _reaction_set
void set_zero(_Matrix< _Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols > &a)
Definition: eigen_utils.h:217
std::vector< std::vector< StateType > > _dnet_rate_dX_s
template<typename CoeffType , typename StateType >
unsigned int Antioch::KineticsEvaluator< CoeffType, StateType >::n_reactions ( ) const
inline

Definition at line 143 of file kinetics_evaluator.h.

144  {
145  return _reaction_set.n_reactions();
146  }
const ReactionSet< CoeffType > & _reaction_set
template<typename CoeffType , typename StateType >
unsigned int Antioch::KineticsEvaluator< CoeffType, StateType >::n_species ( ) const
inline

Definition at line 136 of file kinetics_evaluator.h.

137  {
138  return _chem_mixture.n_species();
139  }
const ChemicalMixture< CoeffType > & _chem_mixture
template<typename CoeffType , typename StateType >
const ReactionSet< CoeffType > & Antioch::KineticsEvaluator< CoeffType, StateType >::reaction_set ( ) const
inline

Definition at line 129 of file kinetics_evaluator.h.

130  {
131  return _reaction_set;
132  }
const ReactionSet< CoeffType > & _reaction_set

Member Data Documentation

template<typename CoeffType = double, typename StateType = CoeffType>
const ChemicalMixture<CoeffType>& Antioch::KineticsEvaluator< CoeffType, StateType >::_chem_mixture
protected

Definition at line 117 of file kinetics_evaluator.h.

template<typename CoeffType = double, typename StateType = CoeffType>
std::vector<StateType> Antioch::KineticsEvaluator< CoeffType, StateType >::_dnet_rate_dT
protected

Definition at line 121 of file kinetics_evaluator.h.

template<typename CoeffType = double, typename StateType = CoeffType>
std::vector<std::vector<StateType> > Antioch::KineticsEvaluator< CoeffType, StateType >::_dnet_rate_dX_s
protected

Definition at line 123 of file kinetics_evaluator.h.

template<typename CoeffType = double, typename StateType = CoeffType>
std::vector<StateType> Antioch::KineticsEvaluator< CoeffType, StateType >::_net_reaction_rates
protected

Definition at line 119 of file kinetics_evaluator.h.

template<typename CoeffType = double, typename StateType = CoeffType>
const ReactionSet<CoeffType>& Antioch::KineticsEvaluator< CoeffType, StateType >::_reaction_set
protected

Definition at line 115 of file kinetics_evaluator.h.


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

Generated on Thu Jul 7 2016 11:09:49 for antioch-0.4.0 by  doxygen 1.8.8