GRINS-0.6.0
reacting_low_mach_navier_stokes_base.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 
26 #ifndef GRINS_REACTING_LOW_MACH_NAVIER_STOKES_BASE_H
27 #define GRINS_REACTING_LOW_MACH_NAVIER_STOKES_BASE_H
28 
29 // GRINS
30 #include "grins_config.h"
31 #include "grins/grins_enums.h"
32 #include "grins/physics.h"
33 #include "grins/pressure_pinning.h"
34 #include "grins/assembly_context.h"
35 
36 namespace GRINS
37 {
38  template<typename Mixture, typename Evaluator>
40  {
41  public:
42 
43  ReactingLowMachNavierStokesBase(const PhysicsName& physics_name, const GetPot& input);
45 
47  virtual void read_input_options( const GetPot& input );
48 
49  virtual void init_variables( libMesh::FEMSystem* system );
50 
52  virtual void set_time_evolving_vars( libMesh::FEMSystem* system );
53 
54  // Context initialization
55  virtual void init_context( AssemblyContext& context );
56 
57  unsigned int n_species() const;
58 
59  libMesh::Real T( const libMesh::Point& p, const AssemblyContext& c ) const;
60 
61  void mass_fractions( const libMesh::Point& p, const AssemblyContext& c,
62  std::vector<libMesh::Real>& mass_fracs ) const;
63 
64  libMesh::Real rho( libMesh::Real T, libMesh::Real p0, libMesh::Real R_mix) const;
65 
66  libMesh::Real get_p0_steady( const AssemblyContext& c, unsigned int qp ) const;
67 
68  libMesh::Real get_p0_steady_side( const AssemblyContext& c, unsigned int qp ) const;
69 
70  libMesh::Real get_p0_steady( const AssemblyContext& c, const libMesh::Point& p ) const;
71 
72  libMesh::Real get_p0_transient( const AssemblyContext& c, unsigned int qp ) const;
73 
74  const Mixture& gas_mixture() const;
75 
76  protected:
77 
78  Mixture _gas_mixture;
79 
80  libMesh::Number _p0;
81 
83  unsigned int _dim;
84 
86  unsigned int _n_species;
87 
89  std::vector<VariableIndex> _species_vars; /* Indicies for species densities */
90  VariableIndex _u_var; /* Index for x-velocity field */
91  VariableIndex _v_var; /* Index for y-velocity field */
92  VariableIndex _w_var; /* Index for z-velocity field */
93  VariableIndex _p_var; /* Index for pressure field */
94  VariableIndex _T_var; /* Index for pressure field */
95  VariableIndex _p0_var; /* Index for thermodynamic pressure */
96 
98  std::vector<std::string> _species_var_names;
100 
103 
105  GRINSEnums::Order _species_order, _V_order, _P_order, _T_order;
106 
108  libMesh::Point _g;
109 
112 
114 
115  libMesh::Real _fixed_rho_value;
116 
117  private:
118 
120 
121  }; // class ReactingLowMachNavierStokesBase
122 
123  template<typename Mixture, typename Evaluator>
124  inline
126  { return _n_species; }
127 
128 
129  template<typename Mixture, typename Evaluator>
130  inline
131  libMesh::Real ReactingLowMachNavierStokesBase<Mixture,Evaluator>::T( const libMesh::Point& p,
132  const AssemblyContext& c ) const
133  { return c.point_value(_T_var,p); }
134 
135  template<typename Mixture, typename Evaluator>
136  inline
138  const AssemblyContext& c,
139  std::vector<libMesh::Real>& mass_fracs ) const
140  {
141  libmesh_assert_equal_to(mass_fracs.size(), this->_n_species);
142 
143  for( unsigned int var = 0; var < this->_n_species; var++ )
144  {
145  mass_fracs[var] = c.point_value(_species_vars[var],p);
146  }
147 
148  return;
149  }
150 
151  template<typename Mixture, typename Evaluator>
152  inline
154  libMesh::Real p0,
155  libMesh::Real R_mix) const
156  {
157  libMesh::Real value = 0;
158  if( this->_fixed_density )
159  value = this->_fixed_rho_value;
160  else
161  value = p0/(R_mix*T);
162 
163  return value;
164  }
165 
166  template<typename Mixture, typename Evaluator>
167  inline
169  unsigned int qp ) const
170  {
171  libMesh::Real p0;
172  if( this->_enable_thermo_press_calc )
173  {
174  p0 = c.interior_value( _p0_var, qp );
175  }
176  else
177  {
178  p0 = _p0;
179  }
180  return p0;
181  }
182 
183  template<typename Mixture, typename Evaluator>
184  inline
186  unsigned int qp ) const
187  {
188  libMesh::Real p0;
189  if( this->_enable_thermo_press_calc )
190  {
191  p0 = c.side_value( _p0_var, qp );
192  }
193  else
194  {
195  p0 = _p0;
196  }
197  return p0;
198  }
199 
200  template<typename Mixture, typename Evaluator>
201  inline
203  const libMesh::Point& p ) const
204  {
205  libMesh::Real p0;
206  if( this->_enable_thermo_press_calc )
207  {
208  p0 = c.point_value( _p0_var, p );
209  }
210  else
211  {
212  p0 = _p0;
213  }
214  return p0;
215  }
216 
217  template<typename Mixture, typename Evaluator>
218  inline
220  unsigned int qp ) const
221  {
222  libMesh::Real p0;
223  if( this->_enable_thermo_press_calc )
224  {
225  p0 = c.fixed_interior_value( _p0_var, qp );
226  }
227  else
228  {
229  p0 = _p0;
230  }
231  return p0;
232  }
233 
234  template<typename Mixture, typename Evaluator>
235  inline
237  {
238  return _gas_mixture;
239  }
240 
241 } // namespace GRINS
242 
243 #endif //GRINS_REACTING_LOW_MACH_NAVIER_STOKES_BASE_H
libMesh::Real rho(libMesh::Real T, libMesh::Real p0, libMesh::Real R_mix) const
unsigned int VariableIndex
More descriptive name of the type used for variable indices.
Definition: var_typedefs.h:40
virtual void set_time_evolving_vars(libMesh::FEMSystem *system)
Sets velocity variables to be time-evolving.
virtual void init_context(AssemblyContext &context)
Initialize context for added physics variables.
std::vector< std::string > _species_var_names
Names of each (owned) variable in the system.
unsigned int _dim
Physical dimension of problem.
Physics abstract base class. Defines API for physics to be added to MultiphysicsSystem.
Definition: physics.h:106
virtual void read_input_options(const GetPot &input)
Read options from GetPot input file.
libMesh::Real get_p0_steady_side(const AssemblyContext &c, unsigned int qp) const
void mass_fractions(const libMesh::Point &p, const AssemblyContext &c, std::vector< libMesh::Real > &mass_fracs) const
libMesh::Real get_p0_transient(const AssemblyContext &c, unsigned int qp) const
GRINS namespace.
libMesh::Real get_p0_steady(const AssemblyContext &c, unsigned int qp) const
GRINSEnums::Order _species_order
Element orders, read from input.
GRINSEnums::FEFamily _species_FE_family
Element type, read from input.
std::string PhysicsName
libMesh::Real T(const libMesh::Point &p, const AssemblyContext &c) const
bool _enable_thermo_press_calc
Flag to enable thermodynamic pressure calculation.
virtual void init_variables(libMesh::FEMSystem *system)
Initialize variables for this physics.
std::vector< VariableIndex > _species_vars
Indices for each (owned) variable;.

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