GRINS-0.7.0
inc_navier_stokes_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 
32 //libMesh
33 #include "libmesh/quadrature.h"
34 
35 namespace GRINS
36 {
37 
38  template<class Mu>
40  const GetPot& input )
41  : IncompressibleNavierStokesStabilizationBase<Mu>(physics_name,input)
42  {}
43 
44  template<class Mu>
46  AssemblyContext& context,
47  CachedValues& /*cache*/ )
48  {
49 #ifdef GRINS_USE_GRVY_TIMERS
50  this->_timer->BeginTimer("IncompressibleNavierStokesSPGSMStabilization::element_time_derivative");
51 #endif
52 
53  // The number of local degrees of freedom in each variable.
54  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
55 
56  // Check number of dofs is same for _flow_vars.u(), v_var and w_var.
57  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.v()).size());
58  if (this->_dim == 3)
59  {
60  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.w()).size());
61  }
62 
63  // Element Jacobian * quadrature weights for interior integration.
64  const std::vector<libMesh::Real> &JxW =
65  context.get_element_fe(this->_flow_vars.u())->get_JxW();
66 
67  // The velocity shape function gradients at interior quadrature points.
68  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
69  context.get_element_fe(this->_flow_vars.u())->get_dphi();
70 
71  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{u}
72  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v()); // R_{v}
73  libMesh::DenseSubVector<libMesh::Number> *Fw = NULL;
74  if(this->_dim == 3)
75  {
76  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
77  }
78 
79  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
80 
81  unsigned int n_qpoints = context.get_element_qrule().n_points();
82 
83  for (unsigned int qp=0; qp != n_qpoints; qp++)
84  {
85  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
86  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
87 
88  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
89  context.interior_value( this->_flow_vars.v(), qp ) );
90  if( this->_dim == 3 )
91  {
92  U(2) = context.interior_value( this->_flow_vars.w(), qp );
93  }
94 
95  // Compute the viscosity at this qp
96  libMesh::Real _mu_qp = this->_mu(context, qp);
97 
98  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, this->_rho, U, _mu_qp, this->_is_steady );
99  libMesh::Real tau_C = this->_stab_helper.compute_tau_continuity( tau_M, g );
100 
101  libMesh::RealGradient RM_s = this->_stab_helper.compute_res_momentum_steady( context, qp, this->_rho, _mu_qp );
102  libMesh::Real RC = this->_stab_helper.compute_res_continuity( context, qp );
103 
104  for (unsigned int i=0; i != n_u_dofs; i++)
105  {
106  Fu(i) += ( - tau_C*RC*u_gradphi[i][qp](0)
107  - tau_M*RM_s(0)*this->_rho*U*u_gradphi[i][qp] )*JxW[qp];
108 
109  Fv(i) += ( - tau_C*RC*u_gradphi[i][qp](1)
110  - tau_M*RM_s(1)*this->_rho*U*u_gradphi[i][qp] )*JxW[qp];
111 
112  if( this->_dim == 3 )
113  {
114  (*Fw)(i) += ( - tau_C*RC*u_gradphi[i][qp](2)
115  - tau_M*RM_s(2)*this->_rho*U*u_gradphi[i][qp] )*JxW[qp];
116  }
117  }
118 
119  if( compute_jacobian )
120  {
121  libmesh_not_implemented();
122  }
123 
124  }
125 
126 #ifdef GRINS_USE_GRVY_TIMERS
127  this->_timer->EndTimer("IncompressibleNavierStokesSPGSMStabilization::element_time_derivative");
128 #endif
129 
130  return;
131  }
132 
133  template<class Mu>
135  AssemblyContext& context,
136  CachedValues& /*cache*/ )
137  {
138 #ifdef GRINS_USE_GRVY_TIMERS
139  this->_timer->BeginTimer("IncompressibleNavierStokesSPGSMStabilization::element_constraint");
140 #endif
141 
142  // The number of local degrees of freedom in each variable.
143  const unsigned int n_p_dofs = context.get_dof_indices(this->_press_var.p()).size();
144 
145  // Element Jacobian * quadrature weights for interior integration.
146  const std::vector<libMesh::Real> &JxW =
147  context.get_element_fe(this->_flow_vars.u())->get_JxW();
148 
149  // The pressure shape functions at interior quadrature points.
150  const std::vector<std::vector<libMesh::RealGradient> >& p_dphi =
151  context.get_element_fe(this->_press_var.p())->get_dphi();
152 
153  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_press_var.p()); // R_{p}
154 
155  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
156 
157  unsigned int n_qpoints = context.get_element_qrule().n_points();
158 
159  for (unsigned int qp=0; qp != n_qpoints; qp++)
160  {
161  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
162  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
163 
164  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
165  context.interior_value( this->_flow_vars.v(), qp ) );
166  if( this->_dim == 3 )
167  {
168  U(2) = context.interior_value( this->_flow_vars.w(), qp );
169  }
170 
171  // Compute the viscosity at this qp
172  libMesh::Real _mu_qp = this->_mu(context, qp);
173 
174  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, this->_rho, U, _mu_qp, this->_is_steady );
175 
176  libMesh::RealGradient RM_s = this->_stab_helper.compute_res_momentum_steady( context, qp, this->_rho, _mu_qp );
177 
178  for (unsigned int i=0; i != n_p_dofs; i++)
179  {
180  Fp(i) += tau_M*RM_s*p_dphi[i][qp]*JxW[qp];
181  }
182 
183  if( compute_jacobian )
184  {
185  libmesh_not_implemented();
186  }
187 
188  }
189 
190 #ifdef GRINS_USE_GRVY_TIMERS
191  this->_timer->EndTimer("IncompressibleNavierStokesSPGSMStabilization::element_constraint");
192 #endif
193 
194  return;
195  }
196 
197  template<class Mu>
199  AssemblyContext& context,
200  CachedValues& /*cache*/ )
201  {
202 #ifdef GRINS_USE_GRVY_TIMERS
203  this->_timer->BeginTimer("IncompressibleNavierStokesSPGSMStabilization::mass_residual");
204 #endif
205 
206  // The number of local degrees of freedom in each variable.
207  const unsigned int n_p_dofs = context.get_dof_indices(this->_press_var.p()).size();
208  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
209 
210  // Element Jacobian * quadrature weights for interior integration.
211  const std::vector<libMesh::Real> &JxW =
212  context.get_element_fe(this->_flow_vars.u())->get_JxW();
213 
214  // The pressure shape functions at interior quadrature points.
215  const std::vector<std::vector<libMesh::RealGradient> >& p_dphi =
216  context.get_element_fe(this->_press_var.p())->get_dphi();
217 
218  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
219  context.get_element_fe(this->_flow_vars.u())->get_dphi();
220 
221  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{p}
222  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v()); // R_{p}
223  libMesh::DenseSubVector<libMesh::Number> *Fw = NULL;
224  if(this->_dim == 3)
225  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
226 
227  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_press_var.p()); // R_{p}
228 
229  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
230 
231  unsigned int n_qpoints = context.get_element_qrule().n_points();
232 
233  for (unsigned int qp=0; qp != n_qpoints; qp++)
234  {
235  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
236  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
237 
238  libMesh::RealGradient U( context.fixed_interior_value( this->_flow_vars.u(), qp ),
239  context.fixed_interior_value( this->_flow_vars.v(), qp ) );
240  // Compute the viscosity at this qp
241  libMesh::Real _mu_qp = this->_mu(context, qp);
242 
243  if( this->_dim == 3 )
244  {
245  U(2) = context.fixed_interior_value( this->_flow_vars.w(), qp );
246  }
247 
248  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, this->_rho, U, _mu_qp, false );
249 
250  libMesh::RealGradient RM_t = this->_stab_helper.compute_res_momentum_transient( context, qp, this->_rho );
251 
252  for (unsigned int i=0; i != n_p_dofs; i++)
253  {
254  Fp(i) -= tau_M*RM_t*p_dphi[i][qp]*JxW[qp];
255  }
256 
257  for (unsigned int i=0; i != n_u_dofs; i++)
258  {
259  Fu(i) -= tau_M*RM_t(0)*this->_rho*U*u_gradphi[i][qp]*JxW[qp];
260 
261  Fv(i) -= tau_M*RM_t(1)*this->_rho*U*u_gradphi[i][qp]*JxW[qp];
262 
263  if( this->_dim == 3 )
264  {
265  (*Fw)(i) -= tau_M*RM_t(2)*this->_rho*U*u_gradphi[i][qp]*JxW[qp];
266  }
267  }
268 
269  if( compute_jacobian )
270  {
271  libmesh_not_implemented();
272  }
273 
274  }
275 
276 #ifdef GRINS_USE_GRVY_TIMERS
277  this->_timer->EndTimer("IncompressibleNavierStokesSPGSMStabilization::mass_residual");
278 #endif
279 
280  return;
281  }
282 
283 } // end namespace GRINS
284 
285 // Instantiate
286 INSTANTIATE_INC_NS_SUBCLASS(IncompressibleNavierStokesSPGSMStabilization);
virtual void element_constraint(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Constraint part(s) of physics for element interiors.
GRINS namespace.
INSTANTIATE_INC_NS_SUBCLASS(IncompressibleNavierStokesSPGSMStabilization)
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Time dependent part(s) of physics for element interiors.
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