GRINS-0.8.0
spalart_allmaras.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 
26 // This class
27 #include "grins/spalart_allmaras.h"
28 
29 // GRINS
30 #include "grins/assembly_context.h"
33 
35 #include "grins/parsed_viscosity.h"
39 
40 // libMesh
41 #include "libmesh/quadrature.h"
42 #include "libmesh/elem.h"
43 
44 namespace GRINS
45 {
46 
47  template<class Mu>
48  SpalartAllmaras<Mu>::SpalartAllmaras(const std::string& physics_name, const GetPot& input )
49  : TurbulenceModelsBase<Mu>(physics_name, input), // Define class variables
50  _flow_vars(GRINSPrivate::VariableWarehouse::get_variable_subclass<VelocityVariable>(VariablesParsing::velocity_variable_name(input,physics_name,VariablesParsing::PHYSICS))),
51  _press_var(GRINSPrivate::VariableWarehouse::get_variable_subclass<PressureFEVariable>(VariablesParsing::press_variable_name(input,physics_name,VariablesParsing::PHYSICS))),
52  _turbulence_vars(GRINSPrivate::VariableWarehouse::get_variable_subclass<TurbulenceFEVariables>(VariablesParsing::turb_variable_name(input,physics_name,VariablesParsing::PHYSICS))),
53  _spalart_allmaras_helper(input),
54  _sa_params(input),
55  _no_of_walls(input("Physics/"+PhysicsNaming::spalart_allmaras()+"/no_of_walls", 0)),
56  _infinite_distance(input("Physics/"+PhysicsNaming::spalart_allmaras()+"/infinite_distance", false))
57  {
59 
60  // Loop over the _no_of_walls and fill the wall_ids set
61  for(unsigned int i = 0; i != _no_of_walls; i++)
62  _wall_ids.insert(input("Physics/"+PhysicsNaming::spalart_allmaras()+"/wall_ids", 0, i ));
63 
64  this->_ic_handler = new GenericICHandler( physics_name, input );
65 
69  }
70 
71  template<class Mu>
72  void SpalartAllmaras<Mu>::init_variables( libMesh::FEMSystem* system )
73  {
74  // Init base class.
76 
77  // Initialize Boundary Mesh
78  this->boundary_mesh.reset(new libMesh::SerialMesh(system->get_mesh().comm() , system->get_mesh().mesh_dimension()) );
79 
80  // Use the _wall_ids set to build the boundary mesh object
81  (system->get_mesh()).boundary_info->sync(_wall_ids, *boundary_mesh);
82 
83  //this->distance_function.reset(new DistanceFunction(system->get_equation_systems(), dynamic_cast<libMesh::UnstructuredMesh&>(system->get_mesh()) ));
84  this->distance_function.reset(new DistanceFunction(system->get_equation_systems(), *boundary_mesh));
85 
86  // For now, we are hacking this. Without this initialize function being called
87  // the distance variable will just be zero. For the channel flow, we are just
88  // going to analytically compute the wall distance
89  //this->distance_function->initialize();
90  }
91 
92  template<class Mu>
94  {
95  // We should prerequest all the data
96  // we will need to build the linear system
97  // or evaluate a quantity of interest.
98  context.get_element_fe(_turbulence_vars.nu())->get_JxW();
99  context.get_element_fe(_turbulence_vars.nu())->get_phi();
100  context.get_element_fe(_turbulence_vars.nu())->get_dphi();
101  context.get_element_fe(_turbulence_vars.nu())->get_xyz();
102 
103  context.get_element_fe(_turbulence_vars.nu())->get_phi();
104  context.get_element_fe(_turbulence_vars.nu())->get_xyz();
105 
106  context.get_side_fe(_turbulence_vars.nu())->get_JxW();
107  context.get_side_fe(_turbulence_vars.nu())->get_phi();
108  context.get_side_fe(_turbulence_vars.nu())->get_dphi();
109  context.get_side_fe(_turbulence_vars.nu())->get_xyz();
110 
111  return;
112  }
113 
114  template<class Mu>
115  void SpalartAllmaras<Mu>::set_time_evolving_vars( libMesh::FEMSystem* system )
116  {
117  // Tell the system to march velocity forward in time, but
118  // leave p as a constraint only
119  system->time_evolving(this->_turbulence_vars.nu(), 1);
120  }
121 
122  template<class Mu>
124  ( bool compute_jacobian, AssemblyContext & context )
125  {
126  // Get a pointer to the current element, we need this for computing
127  // the distance to wall for the quadrature points
128  libMesh::Elem &elem_pointer = context.get_elem();
129 
130  // We get some references to cell-specific data that
131  // will be used to assemble the linear system.
132 
133  // Element Jacobian * quadrature weights for interior integration.
134  const std::vector<libMesh::Real> &JxW =
135  context.get_element_fe(this->_turbulence_vars.nu())->get_JxW();
136 
137  // The viscosity shape functions at interior quadrature points.
138  const std::vector<std::vector<libMesh::Real> >& nu_phi =
139  context.get_element_fe(this->_turbulence_vars.nu())->get_phi();
140 
141  // The viscosity shape function gradients (in global coords.)
142  // at interior quadrature points.
143  const std::vector<std::vector<libMesh::RealGradient> >& nu_gradphi =
144  context.get_element_fe(this->_turbulence_vars.nu())->get_dphi();
145 
146  // The number of local degrees of freedom in each variable.
147  const unsigned int n_nu_dofs = context.get_dof_indices(this->_turbulence_vars.nu()).size();
148 
149  // The subvectors and submatrices we need to fill:
150  //
151  // K_{\alpha \beta} = R_{\alpha},{\beta} = \partial{ R_{\alpha} } / \partial{ {\beta} } (where R denotes residual)
152  // e.g., for \alpha = v and \beta = u we get: K{vu} = R_{v},{u}
153  // Note that Kpu, Kpv, Kpw and Fp comes as constraint.
154 
155  //libMesh::DenseSubMatrix<libMesh::Number> &Knunu = context.get_elem_jacobian(this->_turbulence_vars.nu(), this->_turbulence_vars.nu()); // R_{nu},{nu}
156 
157  libMesh::DenseSubVector<libMesh::Number> &Fnu = context.get_elem_residual(this->_turbulence_vars.nu()); // R_{nu}
158 
159  // Now we will build the element Jacobian and residual.
160  // Constructing the residual requires the solution and its
161  // gradient from the previous timestep. This must be
162  // calculated at each quadrature point by summing the
163  // solution degree-of-freedom values by the appropriate
164  // weight functions.
165  unsigned int n_qpoints = context.get_element_qrule().n_points();
166 
167  // Auto pointer to distance fcn evaluated at quad points
168  libMesh::UniquePtr< libMesh::DenseVector<libMesh::Real> > distance_qp;
169 
170  // Fill the vector of distances to quadrature points
171  distance_qp = this->distance_function->interpolate(&elem_pointer, context.get_element_qrule().get_points());
172 
173  for (unsigned int qp=0; qp != n_qpoints; qp++)
174  {
175  // Compute the solution & its gradient at the old Newton iterate.
176  libMesh::Number nu;
177  nu = context.interior_value(this->_turbulence_vars.nu(), qp);
178 
179  libMesh::Gradient grad_nu;
180  grad_nu = context.interior_gradient(this->_turbulence_vars.nu(), qp);
181 
182  libMesh::Real jac = JxW[qp];
183 
184  // The physical viscosity
185  libMesh::Real mu_qp = this->_mu(context, qp);
186 
187  // The vorticity value
188  libMesh::Real vorticity_value_qp = this->_spalart_allmaras_helper.vorticity(context, qp);
189 
190  // The flow velocity
191  libMesh::Number u,v;
192  u = context.interior_value(this->_flow_vars.u(), qp);
193  v = context.interior_value(this->_flow_vars.v(), qp);
194 
195  libMesh::NumberVectorValue U(u,v);
196  if (this->_flow_vars.dim() == 3)
197  U(2) = context.interior_value(this->_flow_vars.w(), qp);
198 
199  //The source term
200  libMesh::Real S_tilde = this->_sa_params.source_fn(nu, mu_qp, (*distance_qp)(qp), vorticity_value_qp, _infinite_distance);
201 
202  // The ft2 function needed for the negative S-A model
203  libMesh::Real chi = nu/mu_qp;
204  libMesh::Real f_t2 = this->_sa_params.get_c_t3()*exp(-this->_sa_params.get_c_t4()*chi*chi);
205 
206  libMesh::Real source_term = this->_sa_params.get_cb1()*(1 - f_t2)*S_tilde*nu;
207  // For a negative turbulent viscosity nu < 0.0 we need to use a different production function
208  if(nu < 0.0)
209  {
210  source_term = this->_sa_params.get_cb1()*(1 - this->_sa_params.get_c_t3())*vorticity_value_qp*nu;
211  }
212 
213  // The wall destruction term
214  libMesh::Real fw = this->_sa_params.destruction_fn(nu, (*distance_qp)(qp), S_tilde, _infinite_distance);
215 
216  libMesh::Real nud = 0.0;
217  if(_infinite_distance)
218  {
219  nud = 0.0;
220  }
221  else
222  {
223  nud = nu/(*distance_qp)(qp);
224  }
225  libMesh::Real nud2 = nud*nud;
226  libMesh::Real kappa2 = (this->_sa_params.get_kappa())*(this->_sa_params.get_kappa());
227  libMesh::Real cw1 = this->_sa_params.get_cb1()/kappa2 + (1.0 + this->_sa_params.get_cb2())/this->_sa_params.get_sigma();
228  libMesh::Real destruction_term = (cw1*fw - (this->_sa_params.get_cb1()/kappa2)*f_t2)*nud2;
229 
230  // For a negative turbulent viscosity nu < 0.0 we need to use a different production function
231  if(nu < 0.0)
232  {
233  destruction_term = -cw1*nud2;
234  }
235 
236  libMesh::Real fn1 = 1.0;
237  // For a negative turbulent viscosity, fn1 needs to be calculated
238  if(nu < 0.0)
239  {
240  libMesh::Real chi3 = chi*chi*chi;
241  fn1 = (this->_sa_params.get_c_n1() + chi3)/(this->_sa_params.get_c_n1() - chi3);
242  }
243 
244  // First, an i-loop over the viscosity degrees of freedom.
245  for (unsigned int i=0; i != n_nu_dofs; i++)
246  {
247  Fnu(i) += jac *
248  ( -this->_rho*(U*grad_nu)*nu_phi[i][qp] // convection term (assumes incompressibility)
249  +source_term*nu_phi[i][qp] // source term
250  + (1./this->_sa_params.get_sigma())*(-(mu_qp+(fn1*nu))*grad_nu*nu_gradphi[i][qp] + this->_sa_params.get_cb2()*grad_nu*grad_nu*nu_phi[i][qp]) // diffusion term
251  - destruction_term*nu_phi[i][qp]); // destruction term
252 
253  // Compute the jacobian if not using numerical jacobians
254  if (compute_jacobian)
255  {
256  libmesh_not_implemented();
257  } // end - if (compute_jacobian)
258 
259  } // end of the outer dof (i) loop
260  } // end of the quadrature point (qp) loop
261  }
262 
263  template<class K>
265  ( bool compute_jacobian, AssemblyContext & context )
266  {
267  // First we get some references to cell-specific data that
268  // will be used to assemble the linear system.
269 
270  // Element Jacobian * quadrature weights for interior integration.
271  const std::vector<libMesh::Real> &JxW =
272  context.get_element_fe(this->_turbulence_vars.nu())->get_JxW();
273 
274  // The viscosity shape functions at interior quadrature points.
275  const std::vector<std::vector<libMesh::Real> >& nu_phi =
276  context.get_element_fe(this->_turbulence_vars.nu())->get_phi();
277 
278  // The number of local degrees of freedom in each variable.
279  const unsigned int n_nu_dofs = context.get_dof_indices(this->_turbulence_vars.nu()).size();
280 
281  // The subvectors and submatrices we need to fill:
282  libMesh::DenseSubVector<libMesh::Real> &F = context.get_elem_residual(this->_turbulence_vars.nu());
283 
284  //libMesh::DenseSubMatrix<libMesh::Real> &M = context.get_elem_jacobian(this->_turbulence_vars.nu(), this->_turbulence_vars.nu());
285 
286  unsigned int n_qpoints = context.get_element_qrule().n_points();
287 
288  for (unsigned int qp = 0; qp != n_qpoints; ++qp)
289  {
290  // For the mass residual, we need to be a little careful.
291  // The time integrator is handling the time-discretization
292  // for us so we need to supply M(u_fixed)*u' for the residual.
293  // u_fixed will be given by the fixed_interior_value function
294  // while u' will be given by the interior_rate function.
295  libMesh::Real nu_dot;
296  context.interior_rate(this->_turbulence_vars.nu(), qp, nu_dot);
297 
298  for (unsigned int i = 0; i != n_nu_dofs; ++i)
299  {
300  F(i) += -JxW[qp]*this->_rho*nu_dot*nu_phi[i][qp];
301 
302  if( compute_jacobian )
303  {
304  libmesh_not_implemented();
305  }// End of check on Jacobian
306 
307  } // End of element dof loop
308 
309  } // End of the quadrature point loop
310  }
311 
312  template<class Mu>
314  ( const std::string & param_name,
316  const
317  {
318  ParameterUser::register_parameter(param_name, param_pointer);
319  this->_mu.register_parameter(param_name, param_pointer);
320  this->_sa_params.register_parameter(param_name, param_pointer);
321  }
322 
323 } // namespace GRINS
324 
325 // Instantiate
virtual void init_context(AssemblyContext &context)
Initialize context for added physics variables.
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:258
TurbulenceFEVariables & _turbulence_vars
static PhysicsName spalart_allmaras()
void check_var_subdomain_consistency(const FEVariablesBase &var) const
Check that var is enabled on at least the subdomains this Physics is.
Definition: physics.C:174
PressureFEVariable & _press_var
Base class for reading and handling initial conditions for physics classes.
Physics class for Turbulence Models.
GRINS namespace.
virtual void init_variables(libMesh::FEMSystem *system)
Initialize variables for this physics.
INSTANTIATE_TURBULENCE_MODELS_SUBCLASS(SpalartAllmaras)
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context)
Time dependent part(s) of physics for element interiors.
std::set< libMesh::boundary_id_type > _wall_ids
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.
virtual void init_variables(libMesh::FEMSystem *)
Initialize variables for this physics.
Definition: physics.h:115
VelocityVariable & _flow_vars
void set_is_constraint_var(bool is_constraint_var)
Set whether or not this is a "constraint" variable.
virtual void mass_residual(bool compute_jacobian, AssemblyContext &context)
Mass matrix part(s) for element interiors. All boundary terms lie within the time_derivative part...
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.
virtual void set_time_evolving_vars(libMesh::FEMSystem *system)
Sets velocity variables to be time-evolving.

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