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