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

Base class for falloff processes coupled with efficiencies. More...

#include <falloff_threebody_reaction.h>

Inheritance diagram for Antioch::FalloffThreeBodyReaction< CoeffType, FalloffType >:
Antioch::Reaction< CoeffType >

Public Member Functions

 FalloffThreeBodyReaction (const unsigned int n_species, const std::string &equation, const bool &reversible=true, const ReactionType::ReactionType &falloffType=ReactionType::LINDEMANN_FALLOFF_THREE_BODY, const KineticsModel::KineticsModel kin=KineticsModel::KOOIJ)
 Construct a single reaction mechanism. More...
 
 ~FalloffThreeBodyReaction ()
 
template<typename StateType , typename VectorStateType >
StateType compute_forward_rate_coefficient (const VectorStateType &molar_densities, const KineticsConditions< StateType, VectorStateType > &conditions) const
 
template<typename StateType , typename VectorStateType >
void compute_forward_rate_coefficient_and_derivatives (const VectorStateType &molar_densities, const KineticsConditions< StateType, VectorStateType > &conditions, StateType &kfwd, StateType &dkfwd_dT, VectorStateType &dkfkwd_dX) const
 
const FalloffType & F () const
 Return const reference to the falloff object. More...
 
FalloffType & F ()
 Return writeable reference to the falloff object. More...
 
unsigned int n_species () const
 
const std::string & equation () const
 
const std::string & id () const
 
void set_id (const std::string &id)
 set the reaction id. More...
 
ReactionType::ReactionType type () const
 Type of reaction. More...
 
void set_type (const ReactionType::ReactionType type)
 Set the type of reaction. More...
 
void set_reversibility (const bool reversible)
 Set the reversibility of reaction. More...
 
KineticsModel::KineticsModel kinetics_model () const
 Model of kinetics. More...
 
void set_kinetics_model (const KineticsModel::KineticsModel kin)
 Set the model of kinetics. More...
 
bool initialized () const
 
void set_parameter_of_rate (KineticsModel::Parameters parameter, CoeffType new_value, unsigned int n_reaction=0, const std::string &unit="SI")
 reset a parameter from the rate constant More...
 
void set_parameter_of_rate (KineticsModel::Parameters parameter, CoeffType new_value, unsigned int n_reaction, int l, const std::string &unit="SI")
 reset a parameter from the rate constant, vector parameters More...
 
CoeffType get_parameter_of_rate (KineticsModel::Parameters parameter, unsigned int n_reaction=0, const std::string &unit="SI") const
 get a parameter from the rate constant More...
 
CoeffType get_parameter_of_rate (KineticsModel::Parameters parameter, unsigned int n_reaction, const std::string &unit, int l) const
 get a parameter from the rate constant, vectorized version More...
 
void set_parameter_of_chemical_process (ReactionType::Parameters parameter, CoeffType new_value, unsigned int species=std::numeric_limits< unsigned int >::max())
 reset a parameter from the chemical process More...
 
CoeffType get_parameter_of_chemical_process (ReactionType::Parameters parameter, unsigned int species=std::numeric_limits< unsigned int >::max()) const
 get a parameter from the chemical process More...
 
bool reversible () const
 
unsigned int n_reactants () const
 
unsigned int n_products () const
 
const std::string & reactant_name (const unsigned int r) const
 
const std::string & product_name (const unsigned int p) const
 
unsigned int reactant_id (const unsigned int r) const
 
unsigned int product_id (const unsigned int p) const
 
unsigned int reactant_stoichiometric_coefficient (const unsigned int r) const
 
unsigned int product_stoichiometric_coefficient (const unsigned int p) const
 
CoeffType reactant_partial_order (const unsigned int r) const
 
CoeffType product_partial_order (const unsigned int p) const
 
void add_reactant (const std::string &name, const unsigned int r_id, const unsigned int stoichiometric_coeff, const CoeffType partial_order=std::numeric_limits< CoeffType >::infinity())
 
void add_product (const std::string &name, const unsigned int p_id, const unsigned int stoichiometric_coeff, const CoeffType partial_order=std::numeric_limits< CoeffType >::infinity())
 
void clear_reactant ()
 
void clear_product ()
 
void set_efficiency (const std::string &, const unsigned int s, const CoeffType efficiency)
 
CoeffType get_efficiency (const unsigned int s) const
 
CoeffType efficiency (const unsigned int s) const
 
void initialize (unsigned int index=0)
 Computes derived quantities. More...
 
int gamma () const
 
StateType equilibrium_constant (const StateType &P0_RT, const VectorStateType &h_RT_minus_s_R) const
 
void equilibrium_constant_and_derivative (const StateType &T, const StateType &P0_RT, const VectorStateType &h_RT_minus_s_R, const VectorStateType &ddT_h_RT_minus_s_R, StateType &keq, StateType &dkeq_dT) const
 
StateType compute_forward_rate_coefficient (const VectorStateType &molar_densities, const StateType &temp) const
 
void compute_forward_rate_coefficient_and_derivatives (const VectorStateType &molar_densities, const StateType &temp, StateType &kfwd, StateType &dkfwd_dT, VectorStateType &dkfwd_dX) const
 
StateType compute_rate_of_progress (const VectorStateType &molar_densities, const KineticsConditions< StateType, VectorStateType > &conditions, const StateType &P0_RT, const VectorStateType &h_RT_minus_s_R) const
 
StateType compute_rate_of_progress (const VectorStateType &molar_densities, const StateType &temp, const StateType &P0_RT, const VectorStateType &h_RT_minus_s_R) const
 
void compute_rate_of_progress_and_derivatives (const VectorStateType &molar_densities, const ChemicalMixture< CoeffType > &, const KineticsConditions< StateType, VectorStateType > &conditions, const StateType &P0_RT, const VectorStateType &h_RT_minus_s_R, const VectorStateType &dh_RT_minus_s_R_dT, StateType &net_reaction_rate, StateType &dnet_rate_dT, VectorStateType &dnet_rate_dX_s) const
 
void compute_rate_of_progress_and_derivatives (const VectorStateType &molar_densities, const ChemicalMixture< CoeffType > &chem_mixture, const StateType &temp, const StateType &P0_RT, const VectorStateType &h_RT_minus_s_R, const VectorStateType &dh_RT_minus_s_R_dT, StateType &net_reaction_rate, StateType &dnet_rate_dT, VectorStateType &dnet_rate_dX_s) const
 
const KineticsType< CoeffType,
VectorCoeffType > & 
forward_rate (unsigned int ir=0) const
 Return const reference to the forward rate object. More...
 
KineticsType< CoeffType,
VectorCoeffType > & 
forward_rate (unsigned int ir=0)
 Return writeable reference to the forward rate object. More...
 
FalloffType & falloff ()
 Return writeable reference to the falloff object, test type. More...
 
void add_forward_rate (KineticsType< CoeffType, VectorCoeffType > *rate)
 Add a forward rate object. More...
 
void swap_forward_rates (unsigned int irate, unsigned int jrate)
 Swap two forward rates object. More...
 
unsigned int n_rate_constants () const
 Return the number of rate constant objects. More...
 
void print (std::ostream &os=std::cout) const
 Formatted print, by default to std::cout. More...
 

Protected Attributes

FalloffType _F
 
unsigned int _n_species
 
std::string _id
 
std::string _equation
 
std::vector< std::string > _reactant_names
 
std::vector< std::string > _product_names
 
std::vector< unsigned int > _reactant_ids
 
std::vector< unsigned int > _product_ids
 
std::vector< unsigned int > _reactant_stoichiometry
 
std::vector< unsigned int > _product_stoichiometry
 
std::vector< unsigned int > _species_reactant_stoichiometry
 
std::vector< unsigned int > _species_product_stoichiometry
 
std::vector< CoeffType > _species_reactant_partial_order
 
std::vector< CoeffType > _species_product_partial_order
 
std::vector< int > _species_delta_stoichiometry
 
int _gamma
 
bool _initialized
 
bool _reversible
 
ReactionType::ReactionType _type
 
KineticsModel::KineticsModel _kintype
 
std::vector< KineticsType
< CoeffType, VectorCoeffType > * > 
_forward_rate
 The forward reaction rate modified Arrhenius form. More...
 
std::vector< CoeffType > _efficiencies
 efficiencies for three body reactions More...
 

Detailed Description

template<typename CoeffType = double, typename FalloffType = LindemannFalloff<CoeffType>>
class Antioch::FalloffThreeBodyReaction< CoeffType, FalloffType >

Base class for falloff processes coupled with efficiencies.

This class encapsulates a model between falloff and three-body reaction.

Sylvain: I strongly disapprove using this model, see section ``A twist in the physics: three-body falloff" in section 2.2.2 Kinetics computing of the <a href="https://github.com/libantioch/model_doc">model documentation.

It performs the common operations. A falloff rate constant is defined by the equation

\[ k(T,[M]) = \frac{[M] k_0(T)}{1+[M]\frac{k_0(T)}{k_\infty(T)}}\times F \]

with

\[ \begin{split} k_0(T) & = \lim_{[M] \rightarrow 0} k(T,[M])\\ k_\infty(T) & = \lim_{[M] \rightarrow \infty} k(T,[M]) \end{split} \]

$k_0$ and $k_\infty$ being respectively the low and high pressure rate constant limits, considered elementary (as pressure in these conditions is constant). Thus

\[ \begin{split} k_0(T) & = \alpha_0(T) \\ k_\infty(T) & = \alpha_\infty(T) \end{split} \]

with $\alpha_{i,\,i=0,\infty}(T)$ any kinetics model (see base class KineticsType), and $[M]$ the mixture concentration (or pressure, it's equivalent, $[M] = \frac{P}{\mathrm{R}T}$ in the ideal gas model). All reactions are assumed to be reversible, the kinetics models are assumed to be the same. The three-body part is for the computations of $[M]$

\[ [M] = \sum_{s=1}^{n} \epsilon_sc_n \]

We have:

\[ \begin{split} \frac{\partial k(T,[M])}{\partial T} & = k(T,[M]) F \left( \frac{\partial k_0(T)}{\partial T}\frac{1}{k_0(T)} - \frac{\partial k_0(T)}{\partial T}\frac{1}{k_0(T) + \frac{k_\infty(T)}{[M]}} + \frac{\partial k_\infty(T)}{\partial T} \frac{k_0(T)}{k_\infty \left(k_0(T) + \frac{k_\infty(T)}{[M]}\right)} \right) + k(T,[M]) \frac{\partial F}{\partial T} \\[10pt] \frac{\partial k(T,[M])}{\partial c_i} & = F \epsilon_i\frac{k(T,[M])}{[M] + [M]^2\frac{k_0(T)}{k_\infty(T)}} + k(T,[M]) \frac{\partial F}{\partial c_i} \end{split} \]

By default, the falloff is LindemannFalloff and the kinetics model KooijRate.

Definition at line 102 of file falloff_threebody_reaction.h.

Constructor & Destructor Documentation

template<typename CoeffType , typename FalloffType >
Antioch::FalloffThreeBodyReaction< CoeffType, FalloffType >::FalloffThreeBodyReaction ( const unsigned int  n_species,
const std::string &  equation,
const bool &  reversible = true,
const ReactionType::ReactionType falloffType = ReactionType::LINDEMANN_FALLOFF_THREE_BODY,
const KineticsModel::KineticsModel  kin = KineticsModel::KOOIJ 
)
inline

Construct a single reaction mechanism.

Definition at line 144 of file falloff_threebody_reaction.h.

149  :Reaction<CoeffType>(n_species,equation,reversible,falloffType,kin),
150  _F(n_species)
151 
152  {
155  }
const std::string & equation() const
unsigned int n_species() const
std::vector< CoeffType > _efficiencies
efficiencies for three body reactions
Definition: reaction.h:379
template<typename CoeffType , typename FalloffType >
Antioch::FalloffThreeBodyReaction< CoeffType, FalloffType >::~FalloffThreeBodyReaction ( )
inline

Definition at line 160 of file falloff_threebody_reaction.h.

161  {
162  return;
163  }

Member Function Documentation

void Antioch::Reaction< CoeffType, VectorCoeffType >::add_forward_rate ( KineticsType< CoeffType, VectorCoeffType > *  rate)
inherited

Add a forward rate object.

Referenced by tester().

void Antioch::Reaction< CoeffType, VectorCoeffType >::add_product ( const std::string &  name,
const unsigned int  p_id,
const unsigned int  stoichiometric_coeff,
const CoeffType  partial_order = std::numeric_limits<CoeffType>::infinity() 
)
inherited
void Antioch::Reaction< CoeffType, VectorCoeffType >::add_reactant ( const std::string &  name,
const unsigned int  r_id,
const unsigned int  stoichiometric_coeff,
const CoeffType  partial_order = std::numeric_limits<CoeffType>::infinity() 
)
inherited
void Antioch::Reaction< CoeffType, VectorCoeffType >::clear_product ( )
inherited
void Antioch::Reaction< CoeffType, VectorCoeffType >::clear_reactant ( )
inherited
template<typename CoeffType , typename FalloffType >
template<typename StateType , typename VectorStateType >
StateType Antioch::FalloffThreeBodyReaction< CoeffType, FalloffType >::compute_forward_rate_coefficient ( const VectorStateType &  molar_densities,
const KineticsConditions< StateType, VectorStateType > &  conditions 
) const
inline

Definition at line 182 of file falloff_threebody_reaction.h.

References Antioch::KineticsConditions< StateType, VectorStateType >::T(), and Antioch::zero_clone().

Referenced by tester().

184  {
185 //falloff is k(T,[M]) = k0*[M]/(1 + [M]*k0/kinf) * F = k0 * ([M]^-1 + k0 * kinf^-1)^-1 * F
186  StateType M = Antioch::zero_clone(conditions.T());
187  for(unsigned int s = 0; s < molar_densities.size(); s++)
188  {
189  M += this->efficiency(s) * molar_densities[s];
190  }
191 
192  const StateType k0 = (*this->_forward_rate[0])(conditions);
193  const StateType kinf = (*this->_forward_rate[1])(conditions);
194 
195  return k0 / (ant_pow(M,-1) + k0 / kinf) * _F(conditions.T(),M,k0,kinf);
196  }
std::vector< KineticsType< CoeffType, VectorCoeffType > * > _forward_rate
The forward reaction rate modified Arrhenius form.
Definition: reaction.h:376
_Matrix< _Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols > zero_clone(const _Matrix< _Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols > &ex)
Definition: eigen_utils.h:145
CoeffType efficiency(const unsigned int s) const
StateType Antioch::Reaction< CoeffType, VectorCoeffType >::compute_forward_rate_coefficient ( const VectorStateType &  molar_densities,
const StateType &  temp 
) const
inherited
template<typename CoeffType , typename FalloffType >
template<typename StateType , typename VectorStateType >
void Antioch::FalloffThreeBodyReaction< CoeffType, FalloffType >::compute_forward_rate_coefficient_and_derivatives ( const VectorStateType &  molar_densities,
const KineticsConditions< StateType, VectorStateType > &  conditions,
StateType &  kfwd,
StateType &  dkfwd_dT,
VectorStateType &  dkfkwd_dX 
) const
inline

Definition at line 201 of file falloff_threebody_reaction.h.

References Antioch::KineticsConditions< StateType, VectorStateType >::T(), and Antioch::zero_clone().

Referenced by tester().

206  {
207  //variables, k0,kinf and derivatives
208  StateType k0 = Antioch::zero_clone(conditions.T());
209  StateType dk0_dT = Antioch::zero_clone(conditions.T());
210  StateType kinf = Antioch::zero_clone(conditions.T());
211  StateType dkinf_dT = Antioch::zero_clone(conditions.T());
212 
213  this->_forward_rate[0]->compute_rate_and_derivative(conditions,k0,dk0_dT);
214  this->_forward_rate[1]->compute_rate_and_derivative(conditions,kinf,dkinf_dT);
215 
216  StateType M = Antioch::zero_clone(conditions.T());
217  for(unsigned int s = 0; s < molar_densities.size(); s++)
218  {
219  M += this->efficiency(s) * molar_densities[s];
220  }
221 
222  //F
223  StateType f = Antioch::zero_clone(conditions.T());
224  StateType df_dT = Antioch::zero_clone(conditions.T());
225  VectorStateType df_dX = Antioch::zero_clone(molar_densities);
226  _F.F_and_derivatives(conditions.T(),M,k0,dk0_dT,kinf,dkinf_dT,f,df_dT,df_dX);
227 
228 // k(T,[M]) = k0*[M]/(1 + [M]*k0/kinf) * F = k0 * ([M]^-1 + k0 * kinf^-1)^-1 * F
229  kfwd = k0 / (ant_pow(M,-1) + k0/kinf); //temp variable here for calculations dk_d{T,X}
230 
231 //dk_dT = F * dkfwd_dT + kfwd * dF_dT
232 // = F * kfwd * (dk0_dT/k0 - dk0_dT/(kinf/[M] + k0) + k0 * dkinf_dT/(kinf * (kinf/[M] + k0) ) )
233 // + dF_dT * kfwd
234  StateType temp = (kinf/M + k0);
235  dkfwd_dT = f * kfwd * (dk0_dT/k0 - dk0_dT/temp + dkinf_dT * k0/(kinf * temp))
236  + df_dT * kfwd;
237 
238  dkfwd_dX.resize(this->n_species(), kfwd);
239 
240  StateType tmp = f * kfwd / (M + ant_pow(M,2) * k0/kinf);
241 //dkfwd_dX = F * dkfwd_dX + kfwd * dF_dX
242 // = F * epsilon_i * kfwd / ([M] + [M]^2 k0/kinf) + kfwd * dF_dX
243  for(unsigned int ic = 0; ic < this->n_species(); ic++)
244  {
245  dkfwd_dX[ic] = this->efficiency(ic) * ( tmp + df_dX[ic] * kfwd );
246  }
247 
248  kfwd *= f; //finalize
249 
250  return;
251  }
unsigned int n_species() const
std::vector< KineticsType< CoeffType, VectorCoeffType > * > _forward_rate
The forward reaction rate modified Arrhenius form.
Definition: reaction.h:376
_Matrix< _Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols > zero_clone(const _Matrix< _Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols > &ex)
Definition: eigen_utils.h:145
CoeffType efficiency(const unsigned int s) const
void Antioch::Reaction< CoeffType, VectorCoeffType >::compute_forward_rate_coefficient_and_derivatives ( const VectorStateType &  molar_densities,
const StateType &  temp,
StateType &  kfwd,
StateType &  dkfwd_dT,
VectorStateType &  dkfwd_dX 
) const
inherited
StateType Antioch::Reaction< CoeffType, VectorCoeffType >::compute_rate_of_progress ( const VectorStateType &  molar_densities,
const KineticsConditions< StateType, VectorStateType > &  conditions,
const StateType &  P0_RT,
const VectorStateType &  h_RT_minus_s_R 
) const
inherited
StateType Antioch::Reaction< CoeffType, VectorCoeffType >::compute_rate_of_progress ( const VectorStateType &  molar_densities,
const StateType &  temp,
const StateType &  P0_RT,
const VectorStateType &  h_RT_minus_s_R 
) const
inherited
void Antioch::Reaction< CoeffType, VectorCoeffType >::compute_rate_of_progress_and_derivatives ( const VectorStateType &  molar_densities,
const ChemicalMixture< CoeffType > &  ,
const KineticsConditions< StateType, VectorStateType > &  conditions,
const StateType &  P0_RT,
const VectorStateType &  h_RT_minus_s_R,
const VectorStateType &  dh_RT_minus_s_R_dT,
StateType &  net_reaction_rate,
StateType &  dnet_rate_dT,
VectorStateType &  dnet_rate_dX_s 
) const
inherited
void Antioch::Reaction< CoeffType, VectorCoeffType >::compute_rate_of_progress_and_derivatives ( const VectorStateType &  molar_densities,
const ChemicalMixture< CoeffType > &  chem_mixture,
const StateType &  temp,
const StateType &  P0_RT,
const VectorStateType &  h_RT_minus_s_R,
const VectorStateType &  dh_RT_minus_s_R_dT,
StateType &  net_reaction_rate,
StateType &  dnet_rate_dT,
VectorStateType &  dnet_rate_dX_s 
) const
inherited
CoeffType Antioch::Reaction< CoeffType, VectorCoeffType >::efficiency ( const unsigned int  s) const
inherited
const std::string& Antioch::Reaction< CoeffType, VectorCoeffType >::equation ( ) const
inherited
Returns
the equation for this reaction.

Referenced by Antioch::read_reaction_set_data().

StateType Antioch::Reaction< CoeffType, VectorCoeffType >::equilibrium_constant ( const StateType &  P0_RT,
const VectorStateType &  h_RT_minus_s_R 
) const
inherited
void Antioch::Reaction< CoeffType, VectorCoeffType >::equilibrium_constant_and_derivative ( const StateType &  T,
const StateType &  P0_RT,
const VectorStateType &  h_RT_minus_s_R,
const VectorStateType &  ddT_h_RT_minus_s_R,
StateType &  keq,
StateType &  dkeq_dT 
) const
inherited
template<typename CoeffType , typename FalloffType >
const FalloffType & Antioch::FalloffThreeBodyReaction< CoeffType, FalloffType >::F ( ) const
inline

Return const reference to the falloff object.

Definition at line 174 of file falloff_threebody_reaction.h.

Referenced by tester().

175  {
176  return _F;
177  }
template<typename CoeffType , typename FalloffType >
FalloffType & Antioch::FalloffThreeBodyReaction< CoeffType, FalloffType >::F ( )
inline

Return writeable reference to the falloff object.

Definition at line 167 of file falloff_threebody_reaction.h.

168  {
169  return _F;
170  }
FalloffType& Antioch::Reaction< CoeffType, VectorCoeffType >::falloff ( )
inherited

Return writeable reference to the falloff object, test type.

const KineticsType<CoeffType,VectorCoeffType>& Antioch::Reaction< CoeffType, VectorCoeffType >::forward_rate ( unsigned int  ir = 0) const
inherited

Return const reference to the forward rate object.

KineticsType<CoeffType,VectorCoeffType>& Antioch::Reaction< CoeffType, VectorCoeffType >::forward_rate ( unsigned int  ir = 0)
inherited

Return writeable reference to the forward rate object.

int Antioch::Reaction< CoeffType, VectorCoeffType >::gamma ( ) const
inherited
CoeffType Antioch::Reaction< CoeffType, VectorCoeffType >::get_efficiency ( const unsigned int  s) const
inherited
CoeffType Antioch::Reaction< CoeffType, VectorCoeffType >::get_parameter_of_chemical_process ( ReactionType::Parameters  parameter,
unsigned int  species = std::numeric_limits<unsigned int>::max() 
) const
inherited

get a parameter from the chemical process

CoeffType Antioch::Reaction< CoeffType, VectorCoeffType >::get_parameter_of_rate ( KineticsModel::Parameters  parameter,
unsigned int  n_reaction = 0,
const std::string &  unit = "SI" 
) const
inherited

get a parameter from the rate constant

CoeffType Antioch::Reaction< CoeffType, VectorCoeffType >::get_parameter_of_rate ( KineticsModel::Parameters  parameter,
unsigned int  n_reaction,
const std::string &  unit,
int  l 
) const
inherited

get a parameter from the rate constant, vectorized version

const std::string& Antioch::Reaction< CoeffType, VectorCoeffType >::id ( ) const
inherited
Returns
the reaction id.
void Antioch::Reaction< CoeffType, VectorCoeffType >::initialize ( unsigned int  index = 0)
inherited

Computes derived quantities.

bool Antioch::Reaction< CoeffType, VectorCoeffType >::initialized ( ) const
inherited
KineticsModel::KineticsModel Antioch::Reaction< CoeffType, VectorCoeffType >::kinetics_model ( ) const
inherited

Model of kinetics.

unsigned int Antioch::Reaction< CoeffType, VectorCoeffType >::n_products ( ) const
inherited
Returns
the number of products.
unsigned int Antioch::Reaction< CoeffType, VectorCoeffType >::n_rate_constants ( ) const
inherited

Return the number of rate constant objects.

unsigned int Antioch::Reaction< CoeffType, VectorCoeffType >::n_reactants ( ) const
inherited
Returns
the number of reactants.
unsigned int Antioch::Reaction< CoeffType, VectorCoeffType >::n_species ( ) const
inherited
void Antioch::Reaction< CoeffType, VectorCoeffType >::print ( std::ostream &  os = std::cout) const
inherited

Formatted print, by default to std::cout.

unsigned int Antioch::Reaction< CoeffType, VectorCoeffType >::product_id ( const unsigned int  p) const
inherited
const std::string& Antioch::Reaction< CoeffType, VectorCoeffType >::product_name ( const unsigned int  p) const
inherited
Returns
the name of the p th product.
CoeffType Antioch::Reaction< CoeffType, VectorCoeffType >::product_partial_order ( const unsigned int  p) const
inherited
unsigned int Antioch::Reaction< CoeffType, VectorCoeffType >::product_stoichiometric_coefficient ( const unsigned int  p) const
inherited
unsigned int Antioch::Reaction< CoeffType, VectorCoeffType >::reactant_id ( const unsigned int  r) const
inherited
const std::string& Antioch::Reaction< CoeffType, VectorCoeffType >::reactant_name ( const unsigned int  r) const
inherited
Returns
the name of the r th reactant.
CoeffType Antioch::Reaction< CoeffType, VectorCoeffType >::reactant_partial_order ( const unsigned int  r) const
inherited
unsigned int Antioch::Reaction< CoeffType, VectorCoeffType >::reactant_stoichiometric_coefficient ( const unsigned int  r) const
inherited
bool Antioch::Reaction< CoeffType, VectorCoeffType >::reversible ( ) const
inherited
Returns
the reversibility state of reaction.
void Antioch::Reaction< CoeffType, VectorCoeffType >::set_efficiency ( const std::string &  ,
const unsigned int  s,
const CoeffType  efficiency 
)
inherited

Referenced by tester().

void Antioch::Reaction< CoeffType, VectorCoeffType >::set_id ( const std::string &  id)
inherited

set the reaction id.

void Antioch::Reaction< CoeffType, VectorCoeffType >::set_kinetics_model ( const KineticsModel::KineticsModel  kin)
inherited

Set the model of kinetics.

void Antioch::Reaction< CoeffType, VectorCoeffType >::set_parameter_of_chemical_process ( ReactionType::Parameters  parameter,
CoeffType  new_value,
unsigned int  species = std::numeric_limits<unsigned int>::max() 
)
inherited

reset a parameter from the chemical process

void Antioch::Reaction< CoeffType, VectorCoeffType >::set_parameter_of_rate ( KineticsModel::Parameters  parameter,
CoeffType  new_value,
unsigned int  n_reaction = 0,
const std::string &  unit = "SI" 
)
inherited

reset a parameter from the rate constant

void Antioch::Reaction< CoeffType, VectorCoeffType >::set_parameter_of_rate ( KineticsModel::Parameters  parameter,
CoeffType  new_value,
unsigned int  n_reaction,
int  l,
const std::string &  unit = "SI" 
)
inherited

reset a parameter from the rate constant, vector parameters

void Antioch::Reaction< CoeffType, VectorCoeffType >::set_reversibility ( const bool  reversible)
inherited

Set the reversibility of reaction.

void Antioch::Reaction< CoeffType, VectorCoeffType >::set_type ( const ReactionType::ReactionType  type)
inherited

Set the type of reaction.

reversible reactions are considered.

void Antioch::Reaction< CoeffType, VectorCoeffType >::swap_forward_rates ( unsigned int  irate,
unsigned int  jrate 
)
inherited

Swap two forward rates object.

ReactionType::ReactionType Antioch::Reaction< CoeffType, VectorCoeffType >::type ( ) const
inherited

Type of reaction.

reversible reactions are considered.

Member Data Documentation

std::vector<CoeffType> Antioch::Reaction< CoeffType, VectorCoeffType >::_efficiencies
protectedinherited

efficiencies for three body reactions

Definition at line 379 of file reaction.h.

std::string Antioch::Reaction< CoeffType, VectorCoeffType >::_equation
protectedinherited

Definition at line 357 of file reaction.h.

template<typename CoeffType = double, typename FalloffType = LindemannFalloff<CoeffType>>
FalloffType Antioch::FalloffThreeBodyReaction< CoeffType, FalloffType >::_F
protected

Definition at line 137 of file falloff_threebody_reaction.h.

std::vector<KineticsType<CoeffType,VectorCoeffType>* > Antioch::Reaction< CoeffType, VectorCoeffType >::_forward_rate
protectedinherited

The forward reaction rate modified Arrhenius form.

Definition at line 376 of file reaction.h.

int Antioch::Reaction< CoeffType, VectorCoeffType >::_gamma
protectedinherited

Definition at line 369 of file reaction.h.

std::string Antioch::Reaction< CoeffType, VectorCoeffType >::_id
protectedinherited

Definition at line 356 of file reaction.h.

bool Antioch::Reaction< CoeffType, VectorCoeffType >::_initialized
protectedinherited

Definition at line 370 of file reaction.h.

KineticsModel::KineticsModel Antioch::Reaction< CoeffType, VectorCoeffType >::_kintype
protectedinherited

Definition at line 373 of file reaction.h.

unsigned int Antioch::Reaction< CoeffType, VectorCoeffType >::_n_species
protectedinherited

Definition at line 355 of file reaction.h.

std::vector<unsigned int> Antioch::Reaction< CoeffType, VectorCoeffType >::_product_ids
protectedinherited

Definition at line 361 of file reaction.h.

std::vector<std::string> Antioch::Reaction< CoeffType, VectorCoeffType >::_product_names
protectedinherited

Definition at line 359 of file reaction.h.

std::vector<unsigned int> Antioch::Reaction< CoeffType, VectorCoeffType >::_product_stoichiometry
protectedinherited

Definition at line 363 of file reaction.h.

std::vector<unsigned int> Antioch::Reaction< CoeffType, VectorCoeffType >::_reactant_ids
protectedinherited

Definition at line 360 of file reaction.h.

std::vector<std::string> Antioch::Reaction< CoeffType, VectorCoeffType >::_reactant_names
protectedinherited

Definition at line 358 of file reaction.h.

std::vector<unsigned int> Antioch::Reaction< CoeffType, VectorCoeffType >::_reactant_stoichiometry
protectedinherited

Definition at line 362 of file reaction.h.

bool Antioch::Reaction< CoeffType, VectorCoeffType >::_reversible
protectedinherited

Definition at line 371 of file reaction.h.

std::vector<int> Antioch::Reaction< CoeffType, VectorCoeffType >::_species_delta_stoichiometry
protectedinherited

Definition at line 368 of file reaction.h.

std::vector<CoeffType> Antioch::Reaction< CoeffType, VectorCoeffType >::_species_product_partial_order
protectedinherited

Definition at line 367 of file reaction.h.

std::vector<unsigned int> Antioch::Reaction< CoeffType, VectorCoeffType >::_species_product_stoichiometry
protectedinherited

Definition at line 365 of file reaction.h.

std::vector<CoeffType> Antioch::Reaction< CoeffType, VectorCoeffType >::_species_reactant_partial_order
protectedinherited

Definition at line 366 of file reaction.h.

std::vector<unsigned int> Antioch::Reaction< CoeffType, VectorCoeffType >::_species_reactant_stoichiometry
protectedinherited

Definition at line 364 of file reaction.h.

ReactionType::ReactionType Antioch::Reaction< CoeffType, VectorCoeffType >::_type
protectedinherited

Definition at line 372 of file reaction.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