GRINS-0.6.0
gas_recombination_catalytic_wall.C
Go to the documentation of this file.
1 //-----------------------------------------------------------------------bl-
2 //--------------------------------------------------------------------------
3 //
4 // GRINS - General Reacting Incompressible Navier-Stokes
5 //
6 // Copyright (C) 2014-2015 Paul T. Bauman, Roy H. Stogner
7 // Copyright (C) 2010-2013 The PECOS Development Team
8 //
9 // This library is free software; you can redistribute it and/or
10 // modify it under the terms of the Version 2.1 GNU Lesser General
11 // Public License as published by the Free Software Foundation.
12 //
13 // This library is distributed in the hope that it will be useful,
14 // but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // Lesser General Public License for more details.
17 //
18 // You should have received a copy of the GNU Lesser General Public
19 // License along with this library; if not, write to the Free Software
20 // Foundation, Inc. 51 Franklin Street, Fifth Floor,
21 // Boston, MA 02110-1301 USA
22 //
23 //-----------------------------------------------------------------------el-
24 
25 // This class
27 
28 // libMesh
29 #include "libmesh/quadrature.h"
30 
31 namespace GRINS
32 {
33  template<typename Chemistry>
35  CatalycityBase& gamma,
36  const unsigned int reactant_species_idx,
37  const unsigned int product_species_idx )
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  }
44 
45  template<typename Chemistry>
47  {
48  return;
49  }
50 
51  template<typename Chemistry>
52  void GasRecombinationCatalyticWall<Chemistry>::init( const libMesh::FEMSystem& system )
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  }
67 
68  template<typename Chemistry>
70  const CachedValues& cache,
71  const bool request_jacobian )
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  }
133 
134 } // end namespace GRINS
virtual void init(const libMesh::FEMSystem &system)
virtual void apply_fluxes(AssemblyContext &context, const CachedValues &cache, const bool request_jacobian)
GRINS namespace.
const std::vector< std::vector< libMesh::Number > > & get_cached_vector_values(unsigned int quantity) const
const std::vector< libMesh::Number > & get_cached_values(unsigned int quantity) const
Definition: cached_values.C:99
GasRecombinationCatalyticWall(const Chemistry &chem_mixture, CatalycityBase &gamma, const unsigned int reactant_species_idx, const unsigned int product_species_idx)

Generated on Mon Jun 22 2015 21:32:20 for GRINS-0.6.0 by  doxygen 1.8.9.1