32 #include "libmesh/quadrature.h" 
   36   template<
typename Chemistry>
 
   38                                                                            SharedPtr<CatalycityBase>& gamma,
 
   39                                                                            const std::vector<VariableIndex>& species_vars,
 
   42                                                                            unsigned int reactant_species_idx,
 
   43                                                                            unsigned int product_species_idx )
 
   44   : 
CatalyticWallBase<Chemistry>(chem,gamma,species_vars,T_var,p0,reactant_species_idx),
 
   45     _reactant_species_idx(reactant_species_idx),
 
   46     _reactant_var_idx(species_vars[reactant_species_idx]),
 
   47     _product_species_idx(product_species_idx),
 
   48     _product_var_idx(species_vars[product_species_idx])
 
   51   template<
typename Chemistry>
 
   54                                                                            const unsigned int reactant_species_idx,
 
   55                                                                            const unsigned int product_species_idx )
 
   57     _reactant_species_idx(reactant_species_idx),
 
   58     _product_species_idx(product_species_idx)
 
   63   template<
typename Chemistry>
 
   66     libmesh_do_once(libmesh_deprecated());
 
   68     const std::string r_var_name = std::string(
"w_"+this->_chemistry.species_name( this->_reactant_species_idx ) );
 
   70     const std::string p_var_name = std::string(
"w_"+this->_chemistry.species_name( this->_product_species_idx ) );
 
   72     libmesh_assert( system.has_variable( r_var_name ) );
 
   73     libmesh_assert( system.has_variable( p_var_name ) );
 
   75     this->_reactant_var_idx = system.variable_number( r_var_name );
 
   77     this->_product_var_idx = system.variable_number( p_var_name );
 
   80   template<
typename Chemistry>
 
   84                                                             bool is_axisymmetric )
 
   86     libMesh::FEGenericBase<libMesh::Real>* side_fe = NULL;
 
   87     context.get_side_fe( _reactant_var_idx, side_fe );
 
   90     const unsigned int n_var_dofs = context.get_dof_indices(_reactant_var_idx).size();
 
   92     libmesh_assert_equal_to( n_var_dofs, context.get_dof_indices(_product_var_idx).size() );
 
   95     const std::vector<libMesh::Real> &JxW_side = side_fe->get_JxW();
 
   98     const std::vector<std::vector<libMesh::Real> >& var_phi_side = side_fe->get_phi();
 
  101     const std::vector<libMesh::Point>& var_qpoint = side_fe->get_xyz();
 
  104     libMesh::DenseSubVector<libMesh::Number> &F_r_var = context.get_elem_residual(_reactant_var_idx);
 
  107     libMesh::DenseSubVector<libMesh::Number> &F_p_var = context.get_elem_residual(_product_var_idx);
 
  109     unsigned int n_qpoints = context.get_side_qrule().n_points();
 
  111     for (
unsigned int qp=0; qp != n_qpoints; qp++)
 
  113         libMesh::Real jac = JxW_side[qp];
 
  117             const libMesh::Number r = var_qpoint[qp](0);
 
  121         std::vector<libMesh::Real> mass_fractions(this->_chem_ptr->n_species());
 
  122         for( 
unsigned int s = 0; s < this->_chem_ptr->n_species(); s++ )
 
  123           mass_fractions[s] = context.side_value(this->_species_vars[s], qp);
 
  125         libMesh::Real Y_r = mass_fractions[this->_reactant_species_idx];
 
  126         libMesh::Real T =  context.side_value(this->_T_var, qp);
 
  127         libMesh::Real R_mix = this->_chem_ptr->R_mix(mass_fractions);
 
  128         libMesh::Real rho = this->rho( T, this->_p0, R_mix );
 
  130         const libMesh::Real r_value = this->compute_reactant_mass_flux(rho, Y_r, T);
 
  132         const libMesh::Real p_value = -r_value;
 
  134         for (
unsigned int i=0; i != n_var_dofs; i++)
 
  136             F_r_var(i) += sign*r_value*var_phi_side[i][qp]*jac;
 
  138             F_p_var(i) += sign*p_value*var_phi_side[i][qp]*jac;
 
  140             if( compute_jacobian )
 
  141               libmesh_not_implemented();
 
  149   template<
typename Chemistry>
 
  152                                                                const bool request_jacobian )
 
  154     libmesh_do_once(libmesh_deprecated());
 
  156     libMesh::FEGenericBase<libMesh::Real>* side_fe = NULL;
 
  157     context.get_side_fe( _reactant_var_idx, side_fe );
 
  160     const unsigned int n_var_dofs = context.get_dof_indices(_reactant_var_idx).size();
 
  162     libmesh_assert_equal_to( n_var_dofs, context.get_dof_indices(_product_var_idx).size() );
 
  165     const std::vector<libMesh::Real> &JxW_side = side_fe->get_JxW();
 
  168     const std::vector<std::vector<libMesh::Real> >& var_phi_side = side_fe->get_phi();
 
  171     const std::vector<libMesh::Point>& var_qpoint = side_fe->get_xyz();
 
  174     libMesh::DenseSubVector<libMesh::Number> &F_r_var = context.get_elem_residual(_reactant_var_idx);
 
  177     libMesh::DenseSubVector<libMesh::Number> &F_p_var = context.get_elem_residual(_product_var_idx);
 
  179     unsigned int n_qpoints = context.get_side_qrule().n_points();
 
  181     for (
unsigned int qp=0; qp != n_qpoints; qp++)
 
  183         libMesh::Real jac = JxW_side[qp];
 
  187             const libMesh::Number r = var_qpoint[qp](0);
 
  197         const libMesh::Real r_value = this->compute_reactant_mass_flux(rho, Y_r, T);
 
  199         const libMesh::Real p_value = -r_value;
 
  201         for (
unsigned int i=0; i != n_var_dofs; i++)
 
  203             F_r_var(i) += r_value*var_phi_side[i][qp]*jac;
 
  205             F_p_var(i) += p_value*var_phi_side[i][qp]*jac;
 
  207             if( request_jacobian )
 
  209                 libmesh_not_implemented();
 
static bool is_axisymmetric()
unsigned int VariableIndex
More descriptive name of the type used for variable indices. 
virtual void init(const libMesh::FEMSystem &system)
Deprecated. 
virtual void apply_fluxes(AssemblyContext &context, const CachedValues &cache, const bool request_jacobian)
Deprecated. 
const std::vector< std::vector< libMesh::Number > > & get_cached_vector_values(unsigned int quantity) const 
virtual bool eval_flux(bool compute_jacobian, AssemblyContext &context, libMesh::Real sign, bool is_axisymmetric)
GasRecombinationCatalyticWall(SharedPtr< Chemistry > &chem, SharedPtr< CatalycityBase > &gamma, const std::vector< VariableIndex > &species_vars, VariableIndex T_var, libMesh::Real p0, unsigned int reactant_species_idx, unsigned int product_species_idx)
const std::vector< libMesh::Number > & get_cached_values(unsigned int quantity) const