25 #ifndef GRINS_GAS_SOLID_CATALYTIC_WALL_H 
   26 #define GRINS_GAS_SOLID_CATALYTIC_WALL_H 
   35   template<
typename Chemistry>
 
   41                           SharedPtr<CatalycityBase>& gamma,
 
   42                           const std::vector<VariableIndex>& species_vars,
 
   45                           unsigned int reactant_gas_species_idx,
 
   46                           unsigned int reactant_solid_species_idx,
 
   47                           unsigned int product_species_idx);
 
   52                            const unsigned int reactant_gas_species_idx,
 
   53                            const unsigned int reactant_solid_species_idx,
 
   54                            const unsigned int product_species_idx );
 
   59     virtual void init( 
const libMesh::FEMSystem& system );
 
   61     virtual bool eval_flux( 
bool compute_jacobian,
 
   64                             bool is_axisymmetric );
 
   69                                const bool request_jacobian );
 
   72                                                   const libMesh::Real Y_r,
 
   73                                                   const libMesh::Real T );
 
   77                                                            const libMesh::Real Y_r,
 
   78                                                            const libMesh::Real T );
 
   81                                              const libMesh::Real Y_r,
 
   82                                              const libMesh::Real T );
 
   85                                                               const libMesh::Real Y_r,
 
   86                                                               const libMesh::Real T );
 
   89                                                                const std::vector<libMesh::Real> Y,
 
   90                                                                const libMesh::Real T );
 
  104   template<
typename Chemistry>
 
  107                                                                                   const libMesh::Real Y_r,
 
  108                                                                                   const libMesh::Real T )
 
  110     const libMesh::Real rho_r = rho*Y_r;
 
  112     const libMesh::Real omega_dot = this->omega_dot( rho_r, T );
 
  117   template<
typename Chemistry>
 
  120                                                                                            const libMesh::Real Y_r,
 
  121                                                                                            const libMesh::Real T )
 
  123     const libMesh::Real rho_r = rho*Y_r;
 
  125     const libMesh::Real M_r = this->_chemistry.M(_reactant_gas_species_idx);
 
  127     const libMesh::Real M_solid = this->_chemistry.M(_reactant_solid_species_idx);
 
  129     const libMesh::Real omega_dot = this->omega_dot( rho_r, T )*M_solid/M_r;
 
  134   template<
typename Chemistry>
 
  137                                                                              const libMesh::Real Y_r,
 
  138                                                                              const libMesh::Real T )
 
  140     const libMesh::Real rho_r = rho*Y_r;
 
  142     const libMesh::Real M_r = this->_chemistry.M(_reactant_gas_species_idx);
 
  144     const libMesh::Real M_p = this->_chemistry.M(_product_species_idx);
 
  146     const libMesh::Real omega_dot = this->omega_dot( rho_r, T )*M_p/M_r;
 
  151   template<
typename Chemistry>
 
  154                                                                                               const libMesh::Real Y_r,
 
  155                                                                                               const libMesh::Real T )
 
  157     const libMesh::Real rho_r = rho*Y_r;
 
  159     const libMesh::Real M_r = this->_chemistry.M(_reactant_gas_species_idx);
 
  161     const libMesh::Real M_solid = this->_chemistry.M(_reactant_solid_species_idx);
 
  163     const libMesh::Real domega_dot_dT = this->domega_dot_dT( rho_r, T )*M_solid/M_r;
 
  165     return -domega_dot_dT;
 
  168   template<
typename Chemistry>
 
  171                                                                                                const std::vector<libMesh::Real> Y,
 
  172                                                                                                const libMesh::Real T )
 
  174     const libMesh::Real Y_r = Y[_reactant_gas_species_idx];
 
  176     const libMesh::Real rho_r = rho*Y_r;
 
  178     const libMesh::Real M_r = this->_chemistry.M(_reactant_gas_species_idx);
 
  180     const libMesh::Real M_solid = this->_chemistry.M(_reactant_solid_species_idx);
 
  182     const libMesh::Real R = this->_chemistry.R_mix(Y);
 
  184     const libMesh::Real domega_dot_dYs = this->domega_dot_dws( rho_r, Y_r, T, R )*M_solid/M_r;
 
  186     return -domega_dot_dYs;
 
  191 #endif // GRINS_GAS_SOLID_CATALYTIC_WALL_H 
unsigned int VariableIndex
More descriptive name of the type used for variable indices. 
 
const unsigned int _reactant_gas_species_idx
 
virtual void init(const libMesh::FEMSystem &system)
Deprecated. 
 
libMesh::Real compute_reactant_gas_mass_flux(const libMesh::Real rho, const libMesh::Real Y_r, const libMesh::Real T)
 
virtual ~GasSolidCatalyticWall()
 
const unsigned int _reactant_solid_species_idx
 
libMesh::Real compute_reactant_solid_mass_consumption(const libMesh::Real rho, const libMesh::Real Y_r, const libMesh::Real T)
 
const unsigned int _product_species_idx
 
libMesh::Real compute_reactant_solid_mass_consumption_dT(const libMesh::Real rho, const libMesh::Real Y_r, const libMesh::Real T)
 
VariableIndex _product_var_idx
 
virtual bool eval_flux(bool compute_jacobian, AssemblyContext &context, libMesh::Real sign, bool is_axisymmetric)
 
GasSolidCatalyticWall(SharedPtr< Chemistry > &chem, SharedPtr< CatalycityBase > &gamma, const std::vector< VariableIndex > &species_vars, VariableIndex T_var, libMesh::Real p0, unsigned int reactant_gas_species_idx, unsigned int reactant_solid_species_idx, unsigned int product_species_idx)
 
VariableIndex _reactant_gas_var_idx
 
virtual void apply_fluxes(AssemblyContext &context, const CachedValues &cache, const bool request_jacobian)
Deprecated. 
 
libMesh::Real compute_reactant_solid_mass_consumption_dYs(const libMesh::Real rho, const std::vector< libMesh::Real > Y, const libMesh::Real T)
 
libMesh::Real compute_product_mass_flux(const libMesh::Real rho, const libMesh::Real Y_r, const libMesh::Real T)
 
libMesh::Real rho(libMesh::Real T, libMesh::Real p0, libMesh::Real R_mix) const