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

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