GRINS-0.8.0
List of all members | Public Member Functions | Protected Attributes
GRINS::GasSolidCatalyticWall< Chemistry > Class Template Reference

#include <gas_solid_catalytic_wall.h>

Inheritance diagram for GRINS::GasSolidCatalyticWall< Chemistry >:
Inheritance graph
[legend]
Collaboration diagram for GRINS::GasSolidCatalyticWall< Chemistry >:
Collaboration graph
[legend]

Public Member Functions

 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)
 
 GasSolidCatalyticWall (const Chemistry &chem_mixture, CatalycityBase &gamma, const unsigned int reactant_gas_species_idx, const unsigned int reactant_solid_species_idx, const unsigned int product_species_idx)
 Deprecated. More...
 
virtual ~GasSolidCatalyticWall ()
 
virtual void init (const libMesh::FEMSystem &system)
 Deprecated. More...
 
virtual bool eval_flux (bool compute_jacobian, AssemblyContext &context, libMesh::Real sign, bool is_axisymmetric)
 
virtual void apply_fluxes (AssemblyContext &context, const CachedValues &cache, const bool request_jacobian)
 Deprecated. More...
 
libMesh::Real compute_reactant_gas_mass_flux (const libMesh::Real rho, const libMesh::Real Y_r, const libMesh::Real T)
 
libMesh::Real compute_reactant_solid_mass_consumption (const libMesh::Real rho, const libMesh::Real Y_r, 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 compute_reactant_solid_mass_consumption_dT (const libMesh::Real rho, const libMesh::Real Y_r, const libMesh::Real T)
 
libMesh::Real compute_reactant_solid_mass_consumption_dYs (const libMesh::Real rho, const std::vector< libMesh::Real > Y, const libMesh::Real T)
 
- Public Member Functions inherited from GRINS::CatalyticWallBase< Chemistry >
 CatalyticWallBase (SharedPtr< Chemistry > &chem, SharedPtr< CatalycityBase > &gamma, const std::vector< VariableIndex > &species_vars, VariableIndex T_var, libMesh::Real p0, unsigned int reactant_species_idx)
 
 CatalyticWallBase (const Chemistry &chem_mixture, CatalycityBase &gamma, const unsigned int reactant_species_idx)
 Deprecated. More...
 
virtual ~CatalyticWallBase ()
 
libMesh::Real rho (libMesh::Real T, libMesh::Real p0, libMesh::Real R_mix) const
 
libMesh::Real omega_dot (const libMesh::Real rho_s, const libMesh::Real T) const
 $ \rho_s \gamma \sqrt{ \frac{R_s T}{2\pi} } $ More...
 
libMesh::Real domega_dot_dws (const libMesh::Real rho_s, const libMesh::Real w_s, const libMesh::Real T, const libMesh::Real R) const
 
libMesh::Real domega_dot_dT (const libMesh::Real rho_s, const libMesh::Real T) const
 
void set_catalycity_params (const std::vector< libMesh::Real > &params)
 
virtual void register_parameter (const std::string &param_name, libMesh::ParameterMultiAccessor< libMesh::Number > &param_pointer) const
 Each subclass will register its copy of an independent. More...
 
- Public Member Functions inherited from GRINS::NeumannBCAbstract
 NeumannBCAbstract ()
 
virtual ~NeumannBCAbstract ()
 
- Public Member Functions inherited from GRINS::ParameterUser
 ParameterUser (const std::string &user_name)
 
virtual ~ParameterUser ()
 
virtual void set_parameter (libMesh::Number &param_variable, const GetPot &input, const std::string &param_name, libMesh::Number param_default)
 Each subclass can simultaneously read a parameter value from. More...
 
virtual void set_parameter (libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &func, const GetPot &input, const std::string &func_param_name, const std::string &param_default)
 Each subclass can simultaneously read a parsed function from. More...
 
virtual void set_parameter (libMesh::ParsedFEMFunction< libMesh::Number > &func, const GetPot &input, const std::string &func_param_name, const std::string &param_default)
 Each subclass can simultaneously read a parsed function from. More...
 
virtual void move_parameter (const libMesh::Number &old_parameter, libMesh::Number &new_parameter)
 When cloning an object, we need to update parameter pointers. More...
 
virtual void move_parameter (const libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &old_func, libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &new_func)
 When cloning an object, we need to update parameter pointers. More...
 
virtual void move_parameter (const libMesh::ParsedFEMFunction< libMesh::Number > &old_func, libMesh::ParsedFEMFunction< libMesh::Number > &new_func)
 When cloning an object, we need to update parameter pointers. More...
 

Protected Attributes

const unsigned int _reactant_gas_species_idx
 
VariableIndex _reactant_gas_var_idx
 
const unsigned int _reactant_solid_species_idx
 
const unsigned int _product_species_idx
 
VariableIndex _product_var_idx
 
- Protected Attributes inherited from GRINS::CatalyticWallBase< Chemistry >
SharedPtr< Chemistry > _chem_ptr
 
const Chemistry & _chemistry
 Deprecated. More...
 
SharedPtr< CatalycityBase_gamma_ptr
 
libMesh::UniquePtr< CatalycityBase_gamma_s
 Deprecated. More...
 
const libMesh::Real _C
 $ \sqrt{ \frac{R_s}{2\pi} } $ More...
 
std::vector< VariableIndex_species_vars
 
VariableIndex _T_var
 
libMesh::Real _p0
 Thermodynamic pressure. More...
 

Additional Inherited Members

- Static Public Attributes inherited from GRINS::ParameterUser
static std::string zero_vector_function = std::string("{0}")
 A parseable function string with LIBMESH_DIM components, all 0. More...
 
- Protected Member Functions inherited from GRINS::CatalyticWallBase< Chemistry >
libMesh::Real eval_gamma (libMesh::Real T) const
 Temporary helper to deal with intermediate refactoring. More...
 
libMesh::Real eval_gamma_dT (libMesh::Real T) const
 Temporary helper to deal with intermediate refactoring. More...
 

Detailed Description

template<typename Chemistry>
class GRINS::GasSolidCatalyticWall< Chemistry >

Definition at line 36 of file gas_solid_catalytic_wall.h.

Constructor & Destructor Documentation

template<typename Chemistry >
GRINS::GasSolidCatalyticWall< Chemistry >::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 
)

Definition at line 37 of file gas_solid_catalytic_wall.C.

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])
51  {}
const unsigned int _reactant_gas_species_idx
const unsigned int _reactant_solid_species_idx
template<typename Chemistry >
GRINS::GasSolidCatalyticWall< Chemistry >::GasSolidCatalyticWall ( const Chemistry &  chem_mixture,
CatalycityBase gamma,
const unsigned int  reactant_gas_species_idx,
const unsigned int  reactant_solid_species_idx,
const unsigned int  product_species_idx 
)

Deprecated.

Definition at line 54 of file gas_solid_catalytic_wall.C.

59  : CatalyticWallBase<Chemistry>(chem_mixture,gamma,reactant_gas_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)
63  {
64  libmesh_deprecated();
65  }
const unsigned int _reactant_gas_species_idx
const unsigned int _reactant_solid_species_idx
template<typename Chemistry>
virtual GRINS::GasSolidCatalyticWall< Chemistry >::~GasSolidCatalyticWall ( )
inlinevirtual

Definition at line 56 of file gas_solid_catalytic_wall.h.

56 {};

Member Function Documentation

template<typename Chemistry >
void GRINS::GasSolidCatalyticWall< Chemistry >::apply_fluxes ( AssemblyContext context,
const CachedValues cache,
const bool  request_jacobian 
)
virtual

Deprecated.

Implements GRINS::CatalyticWallBase< Chemistry >.

Definition at line 159 of file gas_solid_catalytic_wall.C.

References GRINS::CachedValues::get_cached_values(), GRINS::CachedValues::get_cached_vector_values(), GRINS::Physics::is_axisymmetric(), GRINS::Cache::MASS_FRACTIONS, GRINS::Cache::MIXTURE_DENSITY, and GRINS::Cache::TEMPERATURE.

162  {
163  libmesh_do_once(libmesh_deprecated());
164 
165  libMesh::FEGenericBase<libMesh::Real>* side_fe = NULL;
166  context.get_side_fe( _reactant_gas_var_idx, side_fe );
167 
168  // The number of local degrees of freedom in each variable.
169  const unsigned int n_var_dofs = context.get_dof_indices(_reactant_gas_var_idx).size();
170 
171  libmesh_assert_equal_to( n_var_dofs, context.get_dof_indices(_product_var_idx).size() );
172 
173  // Element Jacobian * quadrature weight for side integration.
174  const std::vector<libMesh::Real> &JxW_side = side_fe->get_JxW();
175 
176  // The var shape functions at side quadrature points.
177  const std::vector<std::vector<libMesh::Real> >& var_phi_side = side_fe->get_phi();
178 
179  // Physical location of the quadrature points
180  const std::vector<libMesh::Point>& var_qpoint = side_fe->get_xyz();
181 
182  // reactant residual
183  libMesh::DenseSubVector<libMesh::Number> &F_r_var = context.get_elem_residual(_reactant_gas_var_idx);
184 
185  // product residual
186  libMesh::DenseSubVector<libMesh::Number> &F_p_var = context.get_elem_residual(_product_var_idx);
187 
188  unsigned int n_qpoints = context.get_side_qrule().n_points();
189 
190  for (unsigned int qp=0; qp != n_qpoints; qp++)
191  {
192  libMesh::Real jac = JxW_side[qp];
193 
195  {
196  const libMesh::Number r = var_qpoint[qp](0);
197  jac *= r;
198  }
199 
200  const libMesh::Real rho = cache.get_cached_values(Cache::MIXTURE_DENSITY)[qp];
201 
202  const libMesh::Real Y_r = cache.get_cached_vector_values(Cache::MASS_FRACTIONS)[qp][this->_reactant_gas_species_idx];
203 
204  const libMesh::Real T = cache.get_cached_values(Cache::TEMPERATURE)[qp];
205 
206  const libMesh::Real r_value = this->compute_reactant_gas_mass_flux(rho, Y_r, T);
207 
208  const libMesh::Real p_value = this->compute_product_mass_flux(rho, Y_r, T);
209 
210  for (unsigned int i=0; i != n_var_dofs; i++)
211  {
212  F_r_var(i) += r_value*var_phi_side[i][qp]*jac;
213 
214  F_p_var(i) += p_value*var_phi_side[i][qp]*jac;
215 
216  if( request_jacobian )
217  {
218  libmesh_not_implemented();
219  }
220  }
221  }
222 
223  return;
224  }
static bool is_axisymmetric()
Definition: physics.h:132
const unsigned int _reactant_gas_species_idx
libMesh::Real compute_reactant_gas_mass_flux(const libMesh::Real rho, const libMesh::Real Y_r, 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
template<typename Chemistry >
libMesh::Real GRINS::GasSolidCatalyticWall< Chemistry >::compute_product_mass_flux ( const libMesh::Real  rho,
const libMesh::Real  Y_r,
const libMesh::Real  T 
)
inline

Definition at line 136 of file gas_solid_catalytic_wall.h.

Referenced by test().

139  {
140  const libMesh::Real rho_r = rho*Y_r;
141 
142  const libMesh::Real M_r = this->_chemistry.M(_reactant_gas_species_idx);
143 
144  const libMesh::Real M_p = this->_chemistry.M(_product_species_idx);
145 
146  const libMesh::Real omega_dot = this->omega_dot( rho_r, T )*M_p/M_r;
147 
148  return omega_dot;
149  }
const unsigned int _reactant_gas_species_idx
libMesh::Real omega_dot(const libMesh::Real rho_s, const libMesh::Real T) const
const Chemistry & _chemistry
Deprecated.
libMesh::Real rho(libMesh::Real T, libMesh::Real p0, libMesh::Real R_mix) const
template<typename Chemistry >
libMesh::Real GRINS::GasSolidCatalyticWall< Chemistry >::compute_reactant_gas_mass_flux ( const libMesh::Real  rho,
const libMesh::Real  Y_r,
const libMesh::Real  T 
)
inline

Definition at line 106 of file gas_solid_catalytic_wall.h.

Referenced by test().

109  {
110  const libMesh::Real rho_r = rho*Y_r;
111 
112  const libMesh::Real omega_dot = this->omega_dot( rho_r, T );
113 
114  return -omega_dot;
115  }
libMesh::Real omega_dot(const libMesh::Real rho_s, const libMesh::Real T) const
libMesh::Real rho(libMesh::Real T, libMesh::Real p0, libMesh::Real R_mix) const
template<typename Chemistry >
libMesh::Real GRINS::GasSolidCatalyticWall< Chemistry >::compute_reactant_solid_mass_consumption ( const libMesh::Real  rho,
const libMesh::Real  Y_r,
const libMesh::Real  T 
)
inline

Definition at line 119 of file gas_solid_catalytic_wall.h.

Referenced by test().

122  {
123  const libMesh::Real rho_r = rho*Y_r;
124 
125  const libMesh::Real M_r = this->_chemistry.M(_reactant_gas_species_idx);
126 
127  const libMesh::Real M_solid = this->_chemistry.M(_reactant_solid_species_idx);
128 
129  const libMesh::Real omega_dot = this->omega_dot( rho_r, T )*M_solid/M_r;
130 
131  return -omega_dot;
132  }
const unsigned int _reactant_gas_species_idx
const unsigned int _reactant_solid_species_idx
libMesh::Real omega_dot(const libMesh::Real rho_s, const libMesh::Real T) const
const Chemistry & _chemistry
Deprecated.
libMesh::Real rho(libMesh::Real T, libMesh::Real p0, libMesh::Real R_mix) const
template<typename Chemistry >
libMesh::Real GRINS::GasSolidCatalyticWall< Chemistry >::compute_reactant_solid_mass_consumption_dT ( const libMesh::Real  rho,
const libMesh::Real  Y_r,
const libMesh::Real  T 
)
inline

Definition at line 153 of file gas_solid_catalytic_wall.h.

156  {
157  const libMesh::Real rho_r = rho*Y_r;
158 
159  const libMesh::Real M_r = this->_chemistry.M(_reactant_gas_species_idx);
160 
161  const libMesh::Real M_solid = this->_chemistry.M(_reactant_solid_species_idx);
162 
163  const libMesh::Real domega_dot_dT = this->domega_dot_dT( rho_r, T )*M_solid/M_r;
164 
165  return -domega_dot_dT;
166  }
const unsigned int _reactant_gas_species_idx
const unsigned int _reactant_solid_species_idx
libMesh::Real domega_dot_dT(const libMesh::Real rho_s, const libMesh::Real T) const
const Chemistry & _chemistry
Deprecated.
libMesh::Real rho(libMesh::Real T, libMesh::Real p0, libMesh::Real R_mix) const
template<typename Chemistry >
libMesh::Real GRINS::GasSolidCatalyticWall< Chemistry >::compute_reactant_solid_mass_consumption_dYs ( const libMesh::Real  rho,
const std::vector< libMesh::Real >  Y,
const libMesh::Real  T 
)
inline

Definition at line 170 of file gas_solid_catalytic_wall.h.

173  {
174  const libMesh::Real Y_r = Y[_reactant_gas_species_idx];
175 
176  const libMesh::Real rho_r = rho*Y_r;
177 
178  const libMesh::Real M_r = this->_chemistry.M(_reactant_gas_species_idx);
179 
180  const libMesh::Real M_solid = this->_chemistry.M(_reactant_solid_species_idx);
181 
182  const libMesh::Real R = this->_chemistry.R_mix(Y);
183 
184  const libMesh::Real domega_dot_dYs = this->domega_dot_dws( rho_r, Y_r, T, R )*M_solid/M_r;
185 
186  return -domega_dot_dYs;
187  }
const unsigned int _reactant_gas_species_idx
const unsigned int _reactant_solid_species_idx
libMesh::Real domega_dot_dws(const libMesh::Real rho_s, const libMesh::Real w_s, const libMesh::Real T, const libMesh::Real R) const
const Chemistry & _chemistry
Deprecated.
libMesh::Real rho(libMesh::Real T, libMesh::Real p0, libMesh::Real R_mix) const
template<typename Chemistry >
bool GRINS::GasSolidCatalyticWall< Chemistry >::eval_flux ( bool  compute_jacobian,
AssemblyContext context,
libMesh::Real  sign,
bool  is_axisymmetric 
)
virtual

Implements GRINS::NeumannBCAbstract.

Definition at line 87 of file gas_solid_catalytic_wall.C.

91  {
92  libMesh::FEGenericBase<libMesh::Real>* side_fe = NULL;
93  context.get_side_fe( _reactant_gas_var_idx, side_fe );
94 
95  // The number of local degrees of freedom in each variable.
96  const unsigned int n_var_dofs = context.get_dof_indices(_reactant_gas_var_idx).size();
97 
98  libmesh_assert_equal_to( n_var_dofs, context.get_dof_indices(_product_var_idx).size() );
99 
100  // Element Jacobian * quadrature weight for side integration.
101  const std::vector<libMesh::Real> &JxW_side = side_fe->get_JxW();
102 
103  // The var shape functions at side quadrature points.
104  const std::vector<std::vector<libMesh::Real> >& var_phi_side = side_fe->get_phi();
105 
106  // Physical location of the quadrature points
107  const std::vector<libMesh::Point>& var_qpoint = side_fe->get_xyz();
108 
109  // reactant residual
110  libMesh::DenseSubVector<libMesh::Number> &F_r_var = context.get_elem_residual(_reactant_gas_var_idx);
111 
112  // product residual
113  libMesh::DenseSubVector<libMesh::Number> &F_p_var = context.get_elem_residual(_product_var_idx);
114 
115  unsigned int n_qpoints = context.get_side_qrule().n_points();
116 
117  for (unsigned int qp=0; qp != n_qpoints; qp++)
118  {
119  libMesh::Real jac = JxW_side[qp];
120 
121  if(is_axisymmetric)
122  {
123  const libMesh::Number r = var_qpoint[qp](0);
124  jac *= r;
125  }
126 
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);
130 
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 );
134 
135  libMesh::Real Y_r = mass_fractions[this->_reactant_gas_species_idx];
136 
137  const libMesh::Real r_value = this->compute_reactant_gas_mass_flux(rho, Y_r, T);
138 
139  const libMesh::Real p_value = this->compute_product_mass_flux(rho, Y_r, T);
140 
141  for (unsigned int i=0; i != n_var_dofs; i++)
142  {
143  F_r_var(i) += sign*r_value*var_phi_side[i][qp]*jac;
144 
145  F_p_var(i) += sign*p_value*var_phi_side[i][qp]*jac;
146 
147  if( compute_jacobian )
148  {
149  libmesh_not_implemented();
150  }
151  }
152  }
153 
154  // We're not computing the Jacobian yet
155  return false;
156  }
const unsigned int _reactant_gas_species_idx
libMesh::Real compute_reactant_gas_mass_flux(const libMesh::Real rho, const libMesh::Real Y_r, const libMesh::Real T)
libMesh::Real _p0
Thermodynamic pressure.
SharedPtr< Chemistry > _chem_ptr
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
template<typename Chemistry >
void GRINS::GasSolidCatalyticWall< Chemistry >::init ( const libMesh::FEMSystem &  system)
virtual

Deprecated.

Reimplemented from GRINS::CatalyticWallBase< Chemistry >.

Definition at line 68 of file gas_solid_catalytic_wall.C.

69  {
70  libmesh_do_once(libmesh_deprecated());
71 
72  const std::string r_var_name = std::string("w_"+this->_chemistry.species_name( this->_reactant_gas_species_idx ) );
73 
74  const std::string p_var_name = std::string("w_"+this->_chemistry.species_name( this->_product_species_idx ) );
75 
76  libmesh_assert( system.has_variable( r_var_name ) );
77  libmesh_assert( system.has_variable( p_var_name ) );
78 
79  this->_reactant_gas_var_idx = system.variable_number( r_var_name );
80 
81  this->_product_var_idx = system.variable_number( p_var_name );
82 
83  return;
84  }
const Chemistry & _chemistry
Deprecated.

Member Data Documentation

template<typename Chemistry>
const unsigned int GRINS::GasSolidCatalyticWall< Chemistry >::_product_species_idx
protected

Definition at line 99 of file gas_solid_catalytic_wall.h.

template<typename Chemistry>
VariableIndex GRINS::GasSolidCatalyticWall< Chemistry >::_product_var_idx
protected

Definition at line 100 of file gas_solid_catalytic_wall.h.

template<typename Chemistry>
const unsigned int GRINS::GasSolidCatalyticWall< Chemistry >::_reactant_gas_species_idx
protected

Definition at line 94 of file gas_solid_catalytic_wall.h.

template<typename Chemistry>
VariableIndex GRINS::GasSolidCatalyticWall< Chemistry >::_reactant_gas_var_idx
protected

Definition at line 95 of file gas_solid_catalytic_wall.h.

template<typename Chemistry>
const unsigned int GRINS::GasSolidCatalyticWall< Chemistry >::_reactant_solid_species_idx
protected

Definition at line 97 of file gas_solid_catalytic_wall.h.


The documentation for this class was generated from the following files:

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