GRINS-0.6.0
spalart_allmaras_spgsm_stab.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 // This class
27 
28 // GRINS
29 #include "grins/assembly_context.h"
31 #include "grins/parsed_viscosity.h"
34 
35 //libMesh
36 #include "libmesh/quadrature.h"
37 
38 namespace GRINS
39 {
40 
41  template<class Mu>
43  const GetPot& input )
44  : SpalartAllmarasStabilizationBase<Mu>(physics_name,input)
45  {
46  this->read_input_options(input);
47 
48  return;
49  }
50 
51  template<class Mu>
53  {
54  return;
55  }
56 
57  template<class Mu>
58  void SpalartAllmarasSPGSMStabilization<Mu>::init_variables( libMesh::FEMSystem* system )
59  {
60  // Init base class variables for stab_helper and distance function initialization
62 
63  return;
64  }
65 
66  template<class Mu>
68  AssemblyContext& context,
69  CachedValues& /*cache*/ )
70  {
71 #ifdef GRINS_USE_GRVY_TIMERS
72  this->_timer->BeginTimer("SpalartAllmarasSPGSMStabilization::element_time_derivative");
73 #endif
74 
75  // Get a pointer to the current element, we need this for computing the distance to wall for the
76  // quadrature points
77  libMesh::Elem &elem_pointer = context.get_elem();
78 
79  // The number of local degrees of freedom in each variable.
80  const unsigned int n_nu_dofs = context.get_dof_indices(this->_turbulence_vars.nu_var()).size();
81 
82  // Element Jacobian * quadrature weights for interior integration.
83  const std::vector<libMesh::Real> &JxW =
84  context.get_element_fe(this->_turbulence_vars.nu_var())->get_JxW();
85 
86  // The viscosity shape function gradients (in global coords.)
87  // at interior quadrature points.
88  const std::vector<std::vector<libMesh::RealGradient> >& nu_gradphi =
89  context.get_element_fe(this->_turbulence_vars.nu_var())->get_dphi();
90 
91  // Quadrature point locations
92  //const std::vector<libMesh::Point>& nu_qpoint =
93  //context.get_element_fe(this->_turbulence_vars.nu_var())->get_xyz();
94 
95  //libMesh::DenseSubMatrix<libMesh::Number> &Knunu = context.get_elem_jacobian(this->_turbulence_vars.nu_var(), this->_turbulence_vars.nu_var()); // R_{nu},{nu}
96 
97  libMesh::DenseSubVector<libMesh::Number> &Fnu = context.get_elem_residual(this->_turbulence_vars.nu_var()); // R_{nu}
98 
99  libMesh::FEBase* fe = context.get_element_fe(this->_turbulence_vars.nu_var());
100 
101  unsigned int n_qpoints = context.get_element_qrule().n_points();
102 
103  // Auto pointer to distance fcn evaluated at quad points
104  libMesh::AutoPtr< libMesh::DenseVector<libMesh::Real> > distance_qp;
105 
106  // Fill the vector of distances to quadrature points
107  distance_qp = this->distance_function->interpolate(&elem_pointer, context.get_element_qrule().get_points());
108 
109  for (unsigned int qp=0; qp != n_qpoints; qp++)
110  {
111  libMesh::Gradient grad_nu;
112  grad_nu = context.interior_gradient(this->_turbulence_vars.nu_var(), qp);
113 
114  libMesh::Real jac = JxW[qp];
115 
116  // The physical viscosity
117  libMesh::Real _mu_qp = this->_mu(context, qp);
118 
119  // To be fixed
120  // For the channel flow we will just set the distance function analytically
121  //(*distance_qp)(qp) = std::min(fabs(y),fabs(1 - y));
122 
123  // The flow velocity
124  libMesh::Number u,v;
125  u = context.interior_value(this->_flow_vars.u_var(), qp);
126  v = context.interior_value(this->_flow_vars.v_var(), qp);
127 
128  libMesh::NumberVectorValue U(u,v);
129  if (this->_dim == 3)
130  U(2) = context.interior_value(this->_flow_vars.w_var(), qp);
131 
132  // Stabilization terms
133 
134  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
135  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
136 
137  libMesh::Real tau_spalart = this->_stab_helper.compute_tau_spalart( context, qp, g, G, this->_rho, U, _mu_qp, this->_is_steady );
138 
139  libMesh::Number RM_spalart = this->_stab_helper.compute_res_spalart_steady( context, qp, this->_rho, _mu_qp, (*distance_qp)(qp) );
140 
141  for (unsigned int i=0; i != n_nu_dofs; i++)
142  {
143  Fnu(i) += jac*( -tau_spalart*RM_spalart*this->_rho*(U*nu_gradphi[i][qp]) );
144  }
145 
146  if( compute_jacobian )
147  {
148  libmesh_not_implemented();
149  }
150 
151  }
152 
153 #ifdef GRINS_USE_GRVY_TIMERS
154  this->_timer->EndTimer("SpalartAllmarasSPGSMStabilization::element_time_derivative");
155 #endif
156 
157  return;
158  }
159 
160  template<class Mu>
162  AssemblyContext& context,
163  CachedValues& /*cache*/ )
164  {
165 #ifdef GRINS_USE_GRVY_TIMERS
166  this->_timer->BeginTimer("SpalartAllmarasSPGSMStabilization::mass_residual");
167 #endif
168 
169  // Get a pointer to the current element, we need this for computing the distance to wall for the
170  // quadrature points
171  libMesh::Elem &elem_pointer = context.get_elem();
172 
173  // The number of local degrees of freedom in each variable.
174  const unsigned int n_nu_dofs = context.get_dof_indices(this->_turbulence_vars.nu_var()).size();
175 
176  // Element Jacobian * quadrature weights for interior integration.
177  const std::vector<libMesh::Real> &JxW =
178  context.get_element_fe(this->_turbulence_vars.nu_var())->get_JxW();
179 
180  // The pressure shape functions at interior quadrature points.
181  const std::vector<std::vector<libMesh::RealGradient> >& nu_gradphi =
182  context.get_element_fe(this->_turbulence_vars.nu_var())->get_dphi();
183 
184  libMesh::DenseSubVector<libMesh::Number> &Fnu = context.get_elem_residual(this->_turbulence_vars.nu_var()); // R_{nu}
185 
186  libMesh::FEBase* fe = context.get_element_fe(this->_turbulence_vars.nu_var());
187 
188  unsigned int n_qpoints = context.get_element_qrule().n_points();
189 
190  // Auto pointer to distance fcn evaluated at quad points
191  libMesh::AutoPtr< libMesh::DenseVector<libMesh::Real> > distance_qp;
192 
193  // Fill the vector of distances to quadrature points
194  distance_qp = this->distance_function->interpolate(&elem_pointer, context.get_element_qrule().get_points());
195 
196  for (unsigned int qp=0; qp != n_qpoints; qp++)
197  {
198  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
199  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
200 
201  libMesh::RealGradient U( context.fixed_interior_value( this->_flow_vars.u_var(), qp ),
202  context.fixed_interior_value( this->_flow_vars.v_var(), qp ) );
203  // Compute the viscosity at this qp
204  libMesh::Real _mu_qp = this->_mu(context, qp);
205 
206  if( this->_dim == 3 )
207  {
208  U(2) = context.fixed_interior_value( this->_flow_vars.w_var(), qp );
209  }
210 
211  libMesh::Real tau_spalart = this->_stab_helper.compute_tau_spalart( context, qp, g, G, this->_rho, U, _mu_qp, this->_is_steady );
212 
213  libMesh::Real RM_spalart = this->_stab_helper.compute_res_spalart_transient( context, qp, this->_rho );
214 
215  for (unsigned int i=0; i != n_nu_dofs; i++)
216  {
217  Fnu(i) += -JxW[qp]*tau_spalart*RM_spalart*this->_rho*(U*nu_gradphi[i][qp]);
218  }
219 
220  if( compute_jacobian )
221  {
222  libmesh_not_implemented();
223  }
224 
225  }
226 
227 #ifdef GRINS_USE_GRVY_TIMERS
228  this->_timer->EndTimer("SpalartAllmarasSPGSMStabilization::mass_residual");
229 #endif
230 
231  return;
232  }
233 
234 } // end namespace GRINS
235 
236 // Instantiate
237 INSTANTIATE_TURBULENCE_MODELS_SUBCLASS(SpalartAllmarasSPGSMStabilization);
INSTANTIATE_TURBULENCE_MODELS_SUBCLASS(SpalartAllmarasSPGSMStabilization)
virtual void init_variables(libMesh::FEMSystem *system)
Initialize variables for this physics.
virtual void init_variables(libMesh::FEMSystem *system)
Initialize variables for this physics.
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Time dependent part(s) of physics for element interiors.
GRINS namespace.
virtual void read_input_options(const GetPot &input)
Read options from GetPot input file. By default, nothing is read.
Definition: physics.C:70
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...

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