GRINS-0.6.0
Public Member Functions | Protected Attributes | List of all members
GRINS::GasRecombinationCatalyticWall< Chemistry > Class Template Reference

#include <gas_recombination_catalytic_wall.h>

Inheritance diagram for GRINS::GasRecombinationCatalyticWall< Chemistry >:
Inheritance graph
[legend]
Collaboration diagram for GRINS::GasRecombinationCatalyticWall< Chemistry >:
Collaboration graph
[legend]

Public Member Functions

 GasRecombinationCatalyticWall (const Chemistry &chem_mixture, CatalycityBase &gamma, const unsigned int reactant_species_idx, const unsigned int product_species_idx)
 
virtual ~GasRecombinationCatalyticWall ()
 
virtual void init (const libMesh::FEMSystem &system)
 
virtual void apply_fluxes (AssemblyContext &context, const CachedValues &cache, const bool request_jacobian)
 
libMesh::Real compute_reactant_mass_flux (const libMesh::Real rho, const libMesh::Real Y_r, const libMesh::Real T)
 
libMesh::Real compute_product_mass_flux (const libMesh::Real rho, const libMesh::Real Y_r, const libMesh::Real T)
 
void set_axisymmetric (bool is_axisymmetric)
 
libMesh::Real omega_dot (const libMesh::Real rho_s, const libMesh::Real T) const
 $ \rho_s \gamma \sqrt{ \frac{R_s T}{2\pi} } $ More...
 
libMesh::Real domega_dot_dws (const libMesh::Real rho_s, const libMesh::Real w_s, const libMesh::Real T, const libMesh::Real R) const
 
libMesh::Real domega_dot_dT (const libMesh::Real rho_s, const libMesh::Real T) const
 
void set_catalycity_params (const std::vector< libMesh::Real > &params)
 

Protected Attributes

const unsigned int _reactant_species_idx
 
VariableIndex _reactant_var_idx
 
const unsigned int _product_species_idx
 
VariableIndex _product_var_idx
 
const Chemistry & _chemistry
 
boost::scoped_ptr< CatalycityBase_gamma_s
 
const libMesh::Real _C
 $ \sqrt{ \frac{R_s}{2\pi} } $ More...
 
bool _is_axisymmetric
 

Detailed Description

template<typename Chemistry>
class GRINS::GasRecombinationCatalyticWall< Chemistry >

Definition at line 36 of file gas_recombination_catalytic_wall.h.

Constructor & Destructor Documentation

template<typename Chemistry >
GRINS::GasRecombinationCatalyticWall< Chemistry >::GasRecombinationCatalyticWall ( const Chemistry &  chem_mixture,
CatalycityBase gamma,
const unsigned int  reactant_species_idx,
const unsigned int  product_species_idx 
)

Definition at line 34 of file gas_recombination_catalytic_wall.C.

38  : CatalyticWallBase<Chemistry>(chem_mixture,gamma,reactant_species_idx),
39  _reactant_species_idx(reactant_species_idx),
40  _product_species_idx(product_species_idx)
41  {
42  return;
43  }
template<typename Chemistry >
GRINS::GasRecombinationCatalyticWall< Chemistry >::~GasRecombinationCatalyticWall ( )
virtual

Definition at line 46 of file gas_recombination_catalytic_wall.C.

47  {
48  return;
49  }

Member Function Documentation

template<typename Chemistry >
void GRINS::GasRecombinationCatalyticWall< Chemistry >::apply_fluxes ( AssemblyContext context,
const CachedValues cache,
const bool  request_jacobian 
)
virtual

Implements GRINS::CatalyticWallBase< Chemistry >.

Definition at line 69 of file gas_recombination_catalytic_wall.C.

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

72  {
73  libMesh::FEGenericBase<libMesh::Real>* side_fe = NULL;
74  context.get_side_fe( _reactant_var_idx, side_fe );
75 
76  // The number of local degrees of freedom in each variable.
77  const unsigned int n_var_dofs = context.get_dof_indices(_reactant_var_idx).size();
78 
79  libmesh_assert_equal_to( n_var_dofs, context.get_dof_indices(_product_var_idx).size() );
80 
81  // Element Jacobian * quadrature weight for side integration.
82  const std::vector<libMesh::Real> &JxW_side = side_fe->get_JxW();
83 
84  // The var shape functions at side quadrature points.
85  const std::vector<std::vector<libMesh::Real> >& var_phi_side = side_fe->get_phi();
86 
87  // Physical location of the quadrature points
88  const std::vector<libMesh::Point>& var_qpoint = side_fe->get_xyz();
89 
90  // reactant residual
91  libMesh::DenseSubVector<libMesh::Number> &F_r_var = context.get_elem_residual(_reactant_var_idx);
92 
93  // product residual
94  libMesh::DenseSubVector<libMesh::Number> &F_p_var = context.get_elem_residual(_product_var_idx);
95 
96  unsigned int n_qpoints = context.get_side_qrule().n_points();
97 
98  for (unsigned int qp=0; qp != n_qpoints; qp++)
99  {
100  libMesh::Real jac = JxW_side[qp];
101 
102  if( this->_is_axisymmetric )
103  {
104  const libMesh::Number r = var_qpoint[qp](0);
105  jac *= r;
106  }
107 
108  const libMesh::Real rho = cache.get_cached_values(Cache::MIXTURE_DENSITY)[qp];
109 
110  const libMesh::Real Y_r = cache.get_cached_vector_values(Cache::MASS_FRACTIONS)[qp][this->_reactant_species_idx];
111 
112  const libMesh::Real T = cache.get_cached_values(Cache::TEMPERATURE)[qp];
113 
114  const libMesh::Real r_value = this->compute_reactant_mass_flux(rho, Y_r, T);
115 
116  const libMesh::Real p_value = -r_value;
117 
118  for (unsigned int i=0; i != n_var_dofs; i++)
119  {
120  F_r_var(i) += r_value*var_phi_side[i][qp]*jac;
121 
122  F_p_var(i) += p_value*var_phi_side[i][qp]*jac;
123 
124  if( request_jacobian )
125  {
126  libmesh_not_implemented();
127  }
128  }
129  }
130 
131  return;
132  }
libMesh::Real compute_reactant_mass_flux(const libMesh::Real rho, const libMesh::Real Y_r, const libMesh::Real T)
template<typename Chemistry >
libMesh::Real GRINS::GasRecombinationCatalyticWall< Chemistry >::compute_product_mass_flux ( const libMesh::Real  rho,
const libMesh::Real  Y_r,
const libMesh::Real  T 
)
inline

Definition at line 86 of file gas_recombination_catalytic_wall.h.

89  {
90  return -(this->compute_reactant_mass_flux(rho,Y_r,T));
91  }
libMesh::Real compute_reactant_mass_flux(const libMesh::Real rho, const libMesh::Real Y_r, const libMesh::Real T)
template<typename Chemistry >
libMesh::Real GRINS::GasRecombinationCatalyticWall< Chemistry >::compute_reactant_mass_flux ( const libMesh::Real  rho,
const libMesh::Real  Y_r,
const libMesh::Real  T 
)
inline

Definition at line 73 of file gas_recombination_catalytic_wall.h.

76  {
77  const libMesh::Real rho_r = rho*Y_r;
78 
79  const libMesh::Real omega_dot = this->omega_dot( rho_r, T );
80 
81  return -omega_dot;
82  }
libMesh::Real omega_dot(const libMesh::Real rho_s, const libMesh::Real T) const
template<typename Chemistry >
libMesh::Real GRINS::CatalyticWallBase< Chemistry >::domega_dot_dT ( const libMesh::Real  rho_s,
const libMesh::Real  T 
) const
inlineinherited

Definition at line 108 of file catalytic_wall_base.h.

Referenced by test().

109  {
110  libMesh::Real sqrtT = std::sqrt(T);
111 
112  return rho_s*_C*( 0.5/sqrtT*(*_gamma_s)(T) + sqrtT*(*_gamma_s).dT(T) );
113  }
const libMesh::Real _C
template<typename Chemistry >
libMesh::Real GRINS::CatalyticWallBase< Chemistry >::domega_dot_dws ( const libMesh::Real  rho_s,
const libMesh::Real  w_s,
const libMesh::Real  T,
const libMesh::Real  R 
) const
inlineinherited

Definition at line 100 of file catalytic_wall_base.h.

Referenced by test().

102  {
103  return (1.0/w_s - rho_s/R)*(this->omega_dot( rho_s, T ));
104  }
libMesh::Real omega_dot(const libMesh::Real rho_s, const libMesh::Real T) const
template<typename Chemistry >
void GRINS::GasRecombinationCatalyticWall< Chemistry >::init ( const libMesh::FEMSystem &  system)
virtual

Reimplemented from GRINS::CatalyticWallBase< Chemistry >.

Definition at line 52 of file gas_recombination_catalytic_wall.C.

53  {
54  const std::string r_var_name = std::string("w_"+this->_chemistry.species_name( this->_reactant_species_idx ) );
55 
56  const std::string p_var_name = std::string("w_"+this->_chemistry.species_name( this->_product_species_idx ) );
57 
58  libmesh_assert( system.has_variable( r_var_name ) );
59  libmesh_assert( system.has_variable( p_var_name ) );
60 
61  this->_reactant_var_idx = system.variable_number( r_var_name );
62 
63  this->_product_var_idx = system.variable_number( p_var_name );
64 
65  return;
66  }
const Chemistry & _chemistry
template<typename Chemistry >
libMesh::Real GRINS::CatalyticWallBase< Chemistry >::omega_dot ( const libMesh::Real  rho_s,
const libMesh::Real  T 
) const
inlineinherited

$ \rho_s \gamma \sqrt{ \frac{R_s T}{2\pi} } $

Definition at line 93 of file catalytic_wall_base.h.

Referenced by test().

94  {
95  return rho_s*(*_gamma_s)(T)*_C*std::sqrt(T);
96  }
const libMesh::Real _C
template<typename Chemistry >
void GRINS::CatalyticWallBase< Chemistry >::set_axisymmetric ( bool  is_axisymmetric)
inherited

Definition at line 64 of file catalytic_wall_base.C.

65  {
66  this->_is_axisymmetric = is_axisymmetric;
67 
68  return;
69  }
template<typename Chemistry >
void GRINS::CatalyticWallBase< Chemistry >::set_catalycity_params ( const std::vector< libMesh::Real > &  params)
inherited

Definition at line 73 of file catalytic_wall_base.C.

74  {
75  _gamma_s->set_params( params );
76  return;
77  }
boost::scoped_ptr< CatalycityBase > _gamma_s

Member Data Documentation

template<typename Chemistry>
const libMesh::Real GRINS::CatalyticWallBase< Chemistry >::_C
protectedinherited

$ \sqrt{ \frac{R_s}{2\pi} } $

Definition at line 84 of file catalytic_wall_base.h.

template<typename Chemistry>
const Chemistry& GRINS::CatalyticWallBase< Chemistry >::_chemistry
protectedinherited

Definition at line 79 of file catalytic_wall_base.h.

template<typename Chemistry>
boost::scoped_ptr<CatalycityBase> GRINS::CatalyticWallBase< Chemistry >::_gamma_s
protectedinherited

Definition at line 81 of file catalytic_wall_base.h.

template<typename Chemistry>
bool GRINS::CatalyticWallBase< Chemistry >::_is_axisymmetric
protectedinherited

Definition at line 86 of file catalytic_wall_base.h.

template<typename Chemistry>
const unsigned int GRINS::GasRecombinationCatalyticWall< Chemistry >::_product_species_idx
protected

Definition at line 66 of file gas_recombination_catalytic_wall.h.

template<typename Chemistry>
VariableIndex GRINS::GasRecombinationCatalyticWall< Chemistry >::_product_var_idx
protected

Definition at line 67 of file gas_recombination_catalytic_wall.h.

template<typename Chemistry>
const unsigned int GRINS::GasRecombinationCatalyticWall< Chemistry >::_reactant_species_idx
protected

Definition at line 63 of file gas_recombination_catalytic_wall.h.

template<typename Chemistry>
VariableIndex GRINS::GasRecombinationCatalyticWall< Chemistry >::_reactant_var_idx
protected

Definition at line 64 of file gas_recombination_catalytic_wall.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