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

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