GRINS-0.8.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-2017 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 // GRINS
29 #include "grins/physics.h"
30 
31 // libMesh
32 #include "libmesh/quadrature.h"
33 
34 namespace GRINS
35 {
36  template<typename Chemistry>
38  SharedPtr<CatalycityBase>& gamma,
39  const std::vector<VariableIndex>& species_vars,
40  VariableIndex T_var,
41  libMesh::Real p0,
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])
49  {}
50 
51  template<typename Chemistry>
53  CatalycityBase& gamma,
54  const unsigned int reactant_species_idx,
55  const unsigned int product_species_idx )
56  : CatalyticWallBase<Chemistry>(chem_mixture,gamma,reactant_species_idx),
57  _reactant_species_idx(reactant_species_idx),
58  _product_species_idx(product_species_idx)
59  {
60  libmesh_deprecated();
61  }
62 
63  template<typename Chemistry>
64  void GasRecombinationCatalyticWall<Chemistry>::init( const libMesh::FEMSystem& system )
65  {
66  libmesh_do_once(libmesh_deprecated());
67 
68  const std::string r_var_name = std::string("w_"+this->_chemistry.species_name( this->_reactant_species_idx ) );
69 
70  const std::string p_var_name = std::string("w_"+this->_chemistry.species_name( this->_product_species_idx ) );
71 
72  libmesh_assert( system.has_variable( r_var_name ) );
73  libmesh_assert( system.has_variable( p_var_name ) );
74 
75  this->_reactant_var_idx = system.variable_number( r_var_name );
76 
77  this->_product_var_idx = system.variable_number( p_var_name );
78  }
79 
80  template<typename Chemistry>
82  AssemblyContext& context,
83  libMesh::Real sign,
84  bool is_axisymmetric )
85  {
86  libMesh::FEGenericBase<libMesh::Real>* side_fe = NULL;
87  context.get_side_fe( _reactant_var_idx, side_fe );
88 
89  // The number of local degrees of freedom in each variable.
90  const unsigned int n_var_dofs = context.get_dof_indices(_reactant_var_idx).size();
91 
92  libmesh_assert_equal_to( n_var_dofs, context.get_dof_indices(_product_var_idx).size() );
93 
94  // Element Jacobian * quadrature weight for side integration.
95  const std::vector<libMesh::Real> &JxW_side = side_fe->get_JxW();
96 
97  // The var shape functions at side quadrature points.
98  const std::vector<std::vector<libMesh::Real> >& var_phi_side = side_fe->get_phi();
99 
100  // Physical location of the quadrature points
101  const std::vector<libMesh::Point>& var_qpoint = side_fe->get_xyz();
102 
103  // reactant residual
104  libMesh::DenseSubVector<libMesh::Number> &F_r_var = context.get_elem_residual(_reactant_var_idx);
105 
106  // product residual
107  libMesh::DenseSubVector<libMesh::Number> &F_p_var = context.get_elem_residual(_product_var_idx);
108 
109  unsigned int n_qpoints = context.get_side_qrule().n_points();
110 
111  for (unsigned int qp=0; qp != n_qpoints; qp++)
112  {
113  libMesh::Real jac = JxW_side[qp];
114 
115  if(is_axisymmetric)
116  {
117  const libMesh::Number r = var_qpoint[qp](0);
118  jac *= r;
119  }
120 
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);
124 
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 );
129 
130  const libMesh::Real r_value = this->compute_reactant_mass_flux(rho, Y_r, T);
131 
132  const libMesh::Real p_value = -r_value;
133 
134  for (unsigned int i=0; i != n_var_dofs; i++)
135  {
136  F_r_var(i) += sign*r_value*var_phi_side[i][qp]*jac;
137 
138  F_p_var(i) += sign*p_value*var_phi_side[i][qp]*jac;
139 
140  if( compute_jacobian )
141  libmesh_not_implemented();
142  }
143  }
144 
145  // We're not computing the Jacobian yet
146  return false;
147  }
148 
149  template<typename Chemistry>
151  const CachedValues& cache,
152  const bool request_jacobian )
153  {
154  libmesh_do_once(libmesh_deprecated());
155 
156  libMesh::FEGenericBase<libMesh::Real>* side_fe = NULL;
157  context.get_side_fe( _reactant_var_idx, side_fe );
158 
159  // The number of local degrees of freedom in each variable.
160  const unsigned int n_var_dofs = context.get_dof_indices(_reactant_var_idx).size();
161 
162  libmesh_assert_equal_to( n_var_dofs, context.get_dof_indices(_product_var_idx).size() );
163 
164  // Element Jacobian * quadrature weight for side integration.
165  const std::vector<libMesh::Real> &JxW_side = side_fe->get_JxW();
166 
167  // The var shape functions at side quadrature points.
168  const std::vector<std::vector<libMesh::Real> >& var_phi_side = side_fe->get_phi();
169 
170  // Physical location of the quadrature points
171  const std::vector<libMesh::Point>& var_qpoint = side_fe->get_xyz();
172 
173  // reactant residual
174  libMesh::DenseSubVector<libMesh::Number> &F_r_var = context.get_elem_residual(_reactant_var_idx);
175 
176  // product residual
177  libMesh::DenseSubVector<libMesh::Number> &F_p_var = context.get_elem_residual(_product_var_idx);
178 
179  unsigned int n_qpoints = context.get_side_qrule().n_points();
180 
181  for (unsigned int qp=0; qp != n_qpoints; qp++)
182  {
183  libMesh::Real jac = JxW_side[qp];
184 
186  {
187  const libMesh::Number r = var_qpoint[qp](0);
188  jac *= r;
189  }
190 
191  const libMesh::Real rho = cache.get_cached_values(Cache::MIXTURE_DENSITY)[qp];
192 
193  const libMesh::Real Y_r = cache.get_cached_vector_values(Cache::MASS_FRACTIONS)[qp][this->_reactant_species_idx];
194 
195  const libMesh::Real T = cache.get_cached_values(Cache::TEMPERATURE)[qp];
196 
197  const libMesh::Real r_value = this->compute_reactant_mass_flux(rho, Y_r, T);
198 
199  const libMesh::Real p_value = -r_value;
200 
201  for (unsigned int i=0; i != n_var_dofs; i++)
202  {
203  F_r_var(i) += r_value*var_phi_side[i][qp]*jac;
204 
205  F_p_var(i) += p_value*var_phi_side[i][qp]*jac;
206 
207  if( request_jacobian )
208  {
209  libmesh_not_implemented();
210  }
211  }
212  }
213  }
214 
215 } // end namespace GRINS
static bool is_axisymmetric()
Definition: physics.h:132
unsigned int VariableIndex
More descriptive name of the type used for variable indices.
Definition: var_typedefs.h:42
virtual void init(const libMesh::FEMSystem &system)
Deprecated.
virtual void apply_fluxes(AssemblyContext &context, const CachedValues &cache, const bool request_jacobian)
Deprecated.
GRINS namespace.
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
Definition: cached_values.C:99

Generated on Tue Dec 19 2017 12:47:27 for GRINS-0.8.0 by  doxygen 1.8.9.1