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