GRINS-0.6.0
gas_solid_catalytic_wall.h
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 #ifndef GRINS_GAS_SOLID_CATALYTIC_WALL_H
26 #define GRINS_GAS_SOLID_CATALYTIC_WALL_H
27 
28 // GRINS
29 #include "grins/var_typedefs.h"
31 
32 namespace GRINS
33 {
34 
35  template<typename Chemistry>
36  class GasSolidCatalyticWall : public CatalyticWallBase<Chemistry>
37  {
38  public:
39 
40  GasSolidCatalyticWall( const Chemistry& chem_mixture,
41  CatalycityBase& gamma,
42  const unsigned int reactant_gas_species_idx,
43  const unsigned int reactant_solid_species_idx,
44  const unsigned int product_species_idx );
45 
46  virtual ~GasSolidCatalyticWall();
47 
48  virtual void init( const libMesh::FEMSystem& system );
49 
50  virtual void apply_fluxes( AssemblyContext& context,
51  const CachedValues& cache,
52  const bool request_jacobian );
53 
54  libMesh::Real compute_reactant_gas_mass_flux( const libMesh::Real rho,
55  const libMesh::Real Y_r,
56  const libMesh::Real T );
57 
58  // kg/m^2-s
59  libMesh::Real compute_reactant_solid_mass_consumption( const libMesh::Real rho,
60  const libMesh::Real Y_r,
61  const libMesh::Real T );
62 
63  libMesh::Real compute_product_mass_flux( const libMesh::Real rho,
64  const libMesh::Real Y_r,
65  const libMesh::Real T );
66 
67  libMesh::Real compute_reactant_solid_mass_consumption_dT( const libMesh::Real rho,
68  const libMesh::Real Y_r,
69  const libMesh::Real T );
70 
71  libMesh::Real compute_reactant_solid_mass_consumption_dYs( const libMesh::Real rho,
72  const std::vector<libMesh::Real> Y,
73  const libMesh::Real T );
74 
75  protected:
76 
77  const unsigned int _reactant_gas_species_idx;
79 
80  const unsigned int _reactant_solid_species_idx;
81 
82  const unsigned int _product_species_idx;
84 
85  };
86 
87  template<typename Chemistry>
88  inline
90  const libMesh::Real Y_r,
91  const libMesh::Real T )
92  {
93  const libMesh::Real rho_r = rho*Y_r;
94 
95  const libMesh::Real omega_dot = this->omega_dot( rho_r, T );
96 
97  return -omega_dot;
98  }
99 
100  template<typename Chemistry>
101  inline
103  const libMesh::Real Y_r,
104  const libMesh::Real T )
105  {
106  const libMesh::Real rho_r = rho*Y_r;
107 
108  const libMesh::Real M_r = this->_chemistry.M(_reactant_gas_species_idx);
109 
110  const libMesh::Real M_solid = this->_chemistry.M(_reactant_solid_species_idx);
111 
112  const libMesh::Real omega_dot = this->omega_dot( rho_r, T )*M_solid/M_r;
113 
114  return -omega_dot;
115  }
116 
117  template<typename Chemistry>
118  inline
119  libMesh::Real GasSolidCatalyticWall<Chemistry>::compute_product_mass_flux( const libMesh::Real rho,
120  const libMesh::Real Y_r,
121  const libMesh::Real T )
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_p = this->_chemistry.M(_product_species_idx);
128 
129  const libMesh::Real omega_dot = this->omega_dot( rho_r, T )*M_p/M_r;
130 
131  return omega_dot;
132  }
133 
134  template<typename Chemistry>
135  inline
137  const libMesh::Real Y_r,
138  const libMesh::Real T )
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_solid = this->_chemistry.M(_reactant_solid_species_idx);
145 
146  const libMesh::Real domega_dot_dT = this->domega_dot_dT( rho_r, T )*M_solid/M_r;
147 
148  return -domega_dot_dT;
149  }
150 
151  template<typename Chemistry>
152  inline
154  const std::vector<libMesh::Real> Y,
155  const libMesh::Real T )
156  {
157  const libMesh::Real Y_r = Y[_reactant_gas_species_idx];
158 
159  const libMesh::Real rho_r = rho*Y_r;
160 
161  const libMesh::Real M_r = this->_chemistry.M(_reactant_gas_species_idx);
162 
163  const libMesh::Real M_solid = this->_chemistry.M(_reactant_solid_species_idx);
164 
165  const libMesh::Real R = this->_chemistry.R_mix(Y);
166 
167  const libMesh::Real domega_dot_dYs = this->domega_dot_dws( rho_r, Y_r, T, R )*M_solid/M_r;
168 
169  return -domega_dot_dYs;
170  }
171 
172 } // end namespace GRINS
173 
174 #endif // GRINS_GAS_SOLID_CATALYTIC_WALL_H
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)
unsigned int VariableIndex
More descriptive name of the type used for variable indices.
Definition: var_typedefs.h:40
const unsigned int _reactant_gas_species_idx
virtual void init(const libMesh::FEMSystem &system)
libMesh::Real compute_reactant_gas_mass_flux(const libMesh::Real rho, const libMesh::Real Y_r, const libMesh::Real T)
const unsigned int _reactant_solid_species_idx
libMesh::Real compute_reactant_solid_mass_consumption(const libMesh::Real rho, const libMesh::Real Y_r, const libMesh::Real T)
GRINS namespace.
libMesh::Real compute_reactant_solid_mass_consumption_dT(const libMesh::Real rho, const libMesh::Real Y_r, const libMesh::Real T)
virtual void apply_fluxes(AssemblyContext &context, const CachedValues &cache, const bool request_jacobian)
libMesh::Real compute_reactant_solid_mass_consumption_dYs(const libMesh::Real rho, const std::vector< libMesh::Real > Y, const libMesh::Real T)
libMesh::Real compute_product_mass_flux(const libMesh::Real rho, const libMesh::Real Y_r, const libMesh::Real T)

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