32 #include "libmesh/quadrature.h"
36 template<
typename Chemistry>
38 SharedPtr<CatalycityBase>& gamma,
39 const std::vector<VariableIndex>& species_vars,
42 unsigned int reactant_gas_species_idx,
43 unsigned int reactant_solid_species_idx,
44 unsigned int product_species_idx)
45 :
CatalyticWallBase<Chemistry>(chem,gamma,species_vars,T_var,p0,reactant_gas_species_idx),
46 _reactant_gas_species_idx(reactant_gas_species_idx),
47 _reactant_gas_var_idx(species_vars[reactant_gas_species_idx]),
48 _reactant_solid_species_idx(reactant_solid_species_idx),
49 _product_species_idx(product_species_idx),
50 _product_var_idx(species_vars[product_species_idx])
53 template<
typename Chemistry>
56 const unsigned int reactant_gas_species_idx,
57 const unsigned int reactant_solid_species_idx,
58 const unsigned int product_species_idx )
60 _reactant_gas_species_idx(reactant_gas_species_idx),
61 _reactant_solid_species_idx(reactant_solid_species_idx),
62 _product_species_idx(product_species_idx)
67 template<
typename Chemistry>
70 libmesh_do_once(libmesh_deprecated());
72 const std::string r_var_name = std::string(
"w_"+this->_chemistry.species_name( this->_reactant_gas_species_idx ) );
74 const std::string p_var_name = std::string(
"w_"+this->_chemistry.species_name( this->_product_species_idx ) );
76 libmesh_assert( system.has_variable( r_var_name ) );
77 libmesh_assert( system.has_variable( p_var_name ) );
79 this->_reactant_gas_var_idx = system.variable_number( r_var_name );
81 this->_product_var_idx = system.variable_number( p_var_name );
86 template<
typename Chemistry>
90 bool is_axisymmetric )
92 libMesh::FEGenericBase<libMesh::Real>* side_fe = NULL;
93 context.get_side_fe( _reactant_gas_var_idx, side_fe );
96 const unsigned int n_var_dofs = context.get_dof_indices(_reactant_gas_var_idx).size();
98 libmesh_assert_equal_to( n_var_dofs, context.get_dof_indices(_product_var_idx).size() );
101 const std::vector<libMesh::Real> &JxW_side = side_fe->get_JxW();
104 const std::vector<std::vector<libMesh::Real> >& var_phi_side = side_fe->get_phi();
107 const std::vector<libMesh::Point>& var_qpoint = side_fe->get_xyz();
110 libMesh::DenseSubVector<libMesh::Number> &F_r_var = context.get_elem_residual(_reactant_gas_var_idx);
113 libMesh::DenseSubVector<libMesh::Number> &F_p_var = context.get_elem_residual(_product_var_idx);
115 unsigned int n_qpoints = context.get_side_qrule().n_points();
117 for (
unsigned int qp=0; qp != n_qpoints; qp++)
119 libMesh::Real jac = JxW_side[qp];
123 const libMesh::Number r = var_qpoint[qp](0);
127 std::vector<libMesh::Real> mass_fractions(this->_chem_ptr->n_species());
128 for(
unsigned int s = 0; s < this->_chem_ptr->n_species(); s++ )
129 mass_fractions[s] = context.side_value(this->_species_vars[s], qp);
131 libMesh::Real T = context.side_value(this->_T_var, qp);
132 libMesh::Real R_mix = this->_chem_ptr->R_mix(mass_fractions);
133 libMesh::Real rho = this->rho( T, this->_p0, R_mix );
135 libMesh::Real Y_r = mass_fractions[this->_reactant_gas_species_idx];
137 const libMesh::Real r_value = this->compute_reactant_gas_mass_flux(rho, Y_r, T);
139 const libMesh::Real p_value = this->compute_product_mass_flux(rho, Y_r, T);
141 for (
unsigned int i=0; i != n_var_dofs; i++)
143 F_r_var(i) += sign*r_value*var_phi_side[i][qp]*jac;
145 F_p_var(i) += sign*p_value*var_phi_side[i][qp]*jac;
147 if( compute_jacobian )
149 libmesh_not_implemented();
158 template<
typename Chemistry>
161 const bool request_jacobian )
163 libmesh_do_once(libmesh_deprecated());
165 libMesh::FEGenericBase<libMesh::Real>* side_fe = NULL;
166 context.get_side_fe( _reactant_gas_var_idx, side_fe );
169 const unsigned int n_var_dofs = context.get_dof_indices(_reactant_gas_var_idx).size();
171 libmesh_assert_equal_to( n_var_dofs, context.get_dof_indices(_product_var_idx).size() );
174 const std::vector<libMesh::Real> &JxW_side = side_fe->get_JxW();
177 const std::vector<std::vector<libMesh::Real> >& var_phi_side = side_fe->get_phi();
180 const std::vector<libMesh::Point>& var_qpoint = side_fe->get_xyz();
183 libMesh::DenseSubVector<libMesh::Number> &F_r_var = context.get_elem_residual(_reactant_gas_var_idx);
186 libMesh::DenseSubVector<libMesh::Number> &F_p_var = context.get_elem_residual(_product_var_idx);
188 unsigned int n_qpoints = context.get_side_qrule().n_points();
190 for (
unsigned int qp=0; qp != n_qpoints; qp++)
192 libMesh::Real jac = JxW_side[qp];
196 const libMesh::Number r = var_qpoint[qp](0);
206 const libMesh::Real r_value = this->compute_reactant_gas_mass_flux(rho, Y_r, T);
208 const libMesh::Real p_value = this->compute_product_mass_flux(rho, Y_r, T);
210 for (
unsigned int i=0; i != n_var_dofs; i++)
212 F_r_var(i) += r_value*var_phi_side[i][qp]*jac;
214 F_p_var(i) += p_value*var_phi_side[i][qp]*jac;
216 if( request_jacobian )
218 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.
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)
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)
virtual void apply_fluxes(AssemblyContext &context, const CachedValues &cache, const bool request_jacobian)
Deprecated.
const std::vector< libMesh::Number > & get_cached_values(unsigned int quantity) const