GRINS-0.7.0
stokes.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 
26 // This class
27 #include "grins/stokes.h"
28 
29 // GRINS
30 #include "grins_config.h"
32 #include "grins/assembly_context.h"
34 #include "grins/multiphysics_sys.h"
35 
36 // libMesh
37 #include "libmesh/fem_context.h"
38 #include "libmesh/quadrature.h"
39 
40 namespace GRINS
41 {
42 
43  template<class Mu>
44  Stokes<Mu>::Stokes(const std::string& physics_name, const GetPot& input )
45  : IncompressibleNavierStokesBase<Mu>(physics_name,
46  PhysicsNaming::stokes(), /* "core" Physics name */
47  input),
48  _p_pinning(input,physics_name),
49  _pin_pressure( input("Physics/"+PhysicsNaming::stokes()+"/pin_pressure", false ) )
50  {
51  this->_ic_handler = new GenericICHandler( physics_name, input );
52  }
53 
54  template<class Mu>
56  {
57  return;
58  }
59 
60  template<class Mu>
62  {
63  if( _pin_pressure )
64  _p_pinning.check_pin_location(system.get_mesh());
65  }
66 
67  template<class Mu>
68  void Stokes<Mu>::element_time_derivative( bool compute_jacobian,
69  AssemblyContext& context,
70  CachedValues& /*cache*/ )
71  {
72 #ifdef GRINS_USE_GRVY_TIMERS
73  this->_timer->BeginTimer("Stokes::element_time_derivative");
74 #endif
75 
76  // The number of local degrees of freedom in each variable.
77  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
78  const unsigned int n_p_dofs = context.get_dof_indices(this->_press_var.p()).size();
79 
80  // Check number of dofs is same for this->_flow_vars.u(), v_var and w_var.
81  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.v()).size());
82  if (this->_dim == 3)
83  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.w()).size());
84 
85  // We get some references to cell-specific data that
86  // will be used to assemble the linear system.
87 
88  // Element Jacobian * quadrature weights for interior integration.
89  const std::vector<libMesh::Real> &JxW =
90  context.get_element_fe(this->_flow_vars.u())->get_JxW();
91 
92  // The velocity shape function gradients (in global coords.)
93  // at interior quadrature points.
94  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
95  context.get_element_fe(this->_flow_vars.u())->get_dphi();
96 
97  // The pressure shape functions at interior quadrature points.
98  const std::vector<std::vector<libMesh::Real> >& p_phi =
99  context.get_element_fe(this->_press_var.p())->get_phi();
100 
101  // The subvectors and submatrices we need to fill:
102  //
103  // K_{\alpha \beta} = R_{\alpha},{\beta} = \partial{ R_{\alpha} } / \partial{ {\beta} } (where R denotes residual)
104  // e.g., for \alpha = v and \beta = u we get: K{vu} = R_{v},{u}
105  // Note that Kpu, Kpv, Kpw and Fp comes as constraint.
106 
107  libMesh::DenseSubMatrix<libMesh::Number> &Kuu = context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.u()); // R_{u},{u}
108  libMesh::DenseSubMatrix<libMesh::Number> &Kvv = context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.v()); // R_{v},{v}
109  libMesh::DenseSubMatrix<libMesh::Number>* Kww = NULL;
110 
111  libMesh::DenseSubMatrix<libMesh::Number> &Kup = context.get_elem_jacobian(this->_flow_vars.u(), this->_press_var.p()); // R_{u},{p}
112  libMesh::DenseSubMatrix<libMesh::Number> &Kvp = context.get_elem_jacobian(this->_flow_vars.v(), this->_press_var.p()); // R_{v},{p}
113  libMesh::DenseSubMatrix<libMesh::Number>* Kwp = NULL;
114 
115  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{u}
116  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v()); // R_{v}
117  libMesh::DenseSubVector<libMesh::Number>* Fw = NULL;
118 
119  if( this->_dim == 3 )
120  {
121  Kww = &context.get_elem_jacobian(this->_flow_vars.w(), this->_flow_vars.w()); // R_{w},{w}
122  Kwp = &context.get_elem_jacobian(this->_flow_vars.w(), this->_press_var.p()); // R_{w},{p}
123  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
124  }
125 
126  // Now we will build the element Jacobian and residual.
127  // Constructing the residual requires the solution and its
128  // gradient from the previous timestep. This must be
129  // calculated at each quadrature point by summing the
130  // solution degree-of-freedom values by the appropriate
131  // weight functions.
132  unsigned int n_qpoints = context.get_element_qrule().n_points();
133 
134  for (unsigned int qp=0; qp != n_qpoints; qp++)
135  {
136  // Compute the solution & its gradient at the old Newton iterate.
137  libMesh::Number p, u, v, w;
138  p = context.interior_value(this->_press_var.p(), qp);
139  u = context.interior_value(this->_flow_vars.u(), qp);
140  v = context.interior_value(this->_flow_vars.v(), qp);
141  if (this->_dim == 3)
142  w = context.interior_value(this->_flow_vars.w(), qp);
143 
144  libMesh::Gradient grad_u, grad_v, grad_w;
145  grad_u = context.interior_gradient(this->_flow_vars.u(), qp);
146  grad_v = context.interior_gradient(this->_flow_vars.v(), qp);
147  if (this->_dim == 3)
148  grad_w = context.interior_gradient(this->_flow_vars.w(), qp);
149 
150  libMesh::NumberVectorValue Uvec (u,v);
151  if (this->_dim == 3)
152  Uvec(2) = w;
153 
154  // Compute the viscosity at this qp
155  libMesh::Real _mu_qp = this->_mu(context, qp);
156 
157  // First, an i-loop over the velocity degrees of freedom.
158  // We know that n_u_dofs == n_v_dofs so we can compute contributions
159  // for both at the same time.
160  for (unsigned int i=0; i != n_u_dofs; i++)
161  {
162  Fu(i) += JxW[qp] *
163  ( p*u_gradphi[i][qp](0) // pressure term
164  -_mu_qp*(u_gradphi[i][qp]*grad_u) ); // diffusion term
165 
166  Fv(i) += JxW[qp] *
167  ( p*u_gradphi[i][qp](1) // pressure term
168  -_mu_qp*(u_gradphi[i][qp]*grad_v) ); // diffusion term
169  if (this->_dim == 3)
170  {
171  (*Fw)(i) += JxW[qp] *
172  ( p*u_gradphi[i][qp](2) // pressure term
173  -_mu_qp*(u_gradphi[i][qp]*grad_w) ); // diffusion term
174  }
175 
176  if (compute_jacobian)
177  {
178  for (unsigned int j=0; j != n_u_dofs; j++)
179  {
180  // TODO: precompute some terms like:
181  // (Uvec*vel_gblgradphivec[j][qp]),
182  // vel_phi[i][qp]*vel_phi[j][qp],
183  // (vel_gblgradphivec[i][qp]*vel_gblgradphivec[j][qp])
184 
185  Kuu(i,j) += JxW[qp] * context.get_elem_solution_derivative() *
186  (-_mu_qp*(u_gradphi[i][qp]*u_gradphi[j][qp])); // diffusion term
187 
188  Kvv(i,j) += JxW[qp] * context.get_elem_solution_derivative() *
189  (-_mu_qp*(u_gradphi[i][qp]*u_gradphi[j][qp])); // diffusion term
190 
191  if (this->_dim == 3)
192  {
193  (*Kww)(i,j) += JxW[qp] * context.get_elem_solution_derivative() *
194  (-_mu_qp*(u_gradphi[i][qp]*u_gradphi[j][qp])); // diffusion term
195  }
196  } // end of the inner dof (j) loop
197 
198  // Matrix contributions for the up, vp and wp couplings
199  for (unsigned int j=0; j != n_p_dofs; j++)
200  {
201  Kup(i,j) += context.get_elem_solution_derivative() * JxW[qp]*u_gradphi[i][qp](0)*p_phi[j][qp];
202  Kvp(i,j) += context.get_elem_solution_derivative() * JxW[qp]*u_gradphi[i][qp](1)*p_phi[j][qp];
203  if (this->_dim == 3)
204  (*Kwp)(i,j) += context.get_elem_solution_derivative() * JxW[qp]*u_gradphi[i][qp](2)*p_phi[j][qp];
205  } // end of the inner dof (j) loop
206 
207  } // end - if (compute_jacobian && context.get_elem_solution_derivative())
208 
209  } // end of the outer dof (i) loop
210  } // end of the quadrature point (qp) loop
211 
212 #ifdef GRINS_USE_GRVY_TIMERS
213  this->_timer->EndTimer("Stokes::element_time_derivative");
214 #endif
215 
216  return;
217  }
218 
219  template<class Mu>
220  void Stokes<Mu>::element_constraint( bool compute_jacobian,
221  AssemblyContext& context,
222  CachedValues& /*cache*/ )
223  {
224 #ifdef GRINS_USE_GRVY_TIMERS
225  this->_timer->BeginTimer("Stokes::element_constraint");
226 #endif
227 
228  // The number of local degrees of freedom in each variable.
229  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
230  const unsigned int n_p_dofs = context.get_dof_indices(this->_press_var.p()).size();
231 
232  // We get some references to cell-specific data that
233  // will be used to assemble the linear system.
234 
235  // Element Jacobian * quadrature weights for interior integration.
236  const std::vector<libMesh::Real> &JxW =
237  context.get_element_fe(this->_flow_vars.u())->get_JxW();
238 
239  // The velocity shape function gradients (in global coords.)
240  // at interior quadrature points.
241  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
242  context.get_element_fe(this->_flow_vars.u())->get_dphi();
243 
244  // The pressure shape functions at interior quadrature points.
245  const std::vector<std::vector<libMesh::Real> >& p_phi =
246  context.get_element_fe(this->_press_var.p())->get_phi();
247 
248  // The subvectors and submatrices we need to fill:
249  //
250  // Kpu, Kpv, Kpw, Fp
251 
252  libMesh::DenseSubMatrix<libMesh::Number> &Kpu = context.get_elem_jacobian(this->_press_var.p(), this->_flow_vars.u()); // R_{p},{u}
253  libMesh::DenseSubMatrix<libMesh::Number> &Kpv = context.get_elem_jacobian(this->_press_var.p(), this->_flow_vars.v()); // R_{p},{v}
254  libMesh::DenseSubMatrix<libMesh::Number>* Kpw = NULL;
255 
256  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_press_var.p()); // R_{p}
257 
258  if( this->_dim == 3 )
259  {
260  Kpw = &context.get_elem_jacobian(this->_press_var.p(), this->_flow_vars.w()); // R_{p},{w}
261  }
262 
263  // Add the constraint given by the continuity equation.
264  unsigned int n_qpoints = context.get_element_qrule().n_points();
265  for (unsigned int qp=0; qp != n_qpoints; qp++)
266  {
267  // Compute the velocity gradient at the old Newton iterate.
268  libMesh::Gradient grad_u, grad_v, grad_w;
269  grad_u = context.interior_gradient(this->_flow_vars.u(), qp);
270  grad_v = context.interior_gradient(this->_flow_vars.v(), qp);
271  if (this->_dim == 3)
272  grad_w = context.interior_gradient(this->_flow_vars.w(), qp);
273 
274  // Now a loop over the pressure degrees of freedom. This
275  // computes the contributions of the continuity equation.
276  for (unsigned int i=0; i != n_p_dofs; i++)
277  {
278  Fp(i) += JxW[qp] * p_phi[i][qp] *
279  (grad_u(0) + grad_v(1));
280  if (this->_dim == 3)
281  Fp(i) += JxW[qp] * p_phi[i][qp] *
282  (grad_w(2));
283 
284  if (compute_jacobian)
285  {
286  for (unsigned int j=0; j != n_u_dofs; j++)
287  {
288  Kpu(i,j) += context.get_elem_solution_derivative() * JxW[qp]*p_phi[i][qp]*u_gradphi[j][qp](0);
289  Kpv(i,j) += context.get_elem_solution_derivative() * JxW[qp]*p_phi[i][qp]*u_gradphi[j][qp](1);
290  if (this->_dim == 3)
291  (*Kpw)(i,j) += context.get_elem_solution_derivative() * JxW[qp]*p_phi[i][qp]*u_gradphi[j][qp](2);
292  } // end of the inner dof (j) loop
293 
294  } // end - if (compute_jacobian && context.get_elem_solution_derivative())
295 
296  } // end of the outer dof (i) loop
297  } // end of the quadrature point (qp) loop
298 
299 
300  // Pin p = p_value at p_point
301  if( _pin_pressure )
302  {
303  _p_pinning.pin_value( context, compute_jacobian, this->_press_var.p() );
304  }
305 
306 
307 #ifdef GRINS_USE_GRVY_TIMERS
308  this->_timer->EndTimer("Stokes::element_constraint");
309 #endif
310 
311  return;
312  }
313 
314  template<class Mu>
315  void Stokes<Mu>::mass_residual( bool compute_jacobian,
316  AssemblyContext& context,
317  CachedValues& /*cache*/)
318  {
319  // Element Jacobian * quadrature weights for interior integration
320  // We assume the same for each flow variable
321  const std::vector<libMesh::Real> &JxW =
322  context.get_element_fe(this->_flow_vars.u())->get_JxW();
323 
324  // The shape functions at interior quadrature points.
325  // We assume the same for each flow variable
326  const std::vector<std::vector<libMesh::Real> >& u_phi =
327  context.get_element_fe(this->_flow_vars.u())->get_phi();
328 
329  // The number of local degrees of freedom in each variable
330  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
331 
332  // The subvectors and submatrices we need to fill:
333  libMesh::DenseSubVector<libMesh::Real> &F_u = context.get_elem_residual(this->_flow_vars.u());
334  libMesh::DenseSubVector<libMesh::Real> &F_v = context.get_elem_residual(this->_flow_vars.v());
335  libMesh::DenseSubVector<libMesh::Real>* F_w = NULL;
336 
337  libMesh::DenseSubMatrix<libMesh::Real> &M_uu = context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.u());
338  libMesh::DenseSubMatrix<libMesh::Real> &M_vv = context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.v());
339  libMesh::DenseSubMatrix<libMesh::Real>* M_ww = NULL;
340 
341  if( this->_dim == 3 )
342  {
343  F_w = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
344  M_ww = &context.get_elem_jacobian(this->_flow_vars.w(), this->_flow_vars.w());
345  }
346 
347  unsigned int n_qpoints = context.get_element_qrule().n_points();
348 
349  for (unsigned int qp = 0; qp != n_qpoints; ++qp)
350  {
351  // For the mass residual, we need to be a little careful.
352  // The time integrator is handling the time-discretization
353  // for us so we need to supply M(u_fixed)*u' for the residual.
354  // u_fixed will be given by the fixed_interior_value function
355  // while u' will be given by the interior_rate function.
356  libMesh::Real u_dot, v_dot, w_dot = 0.0;
357  context.interior_rate(this->_flow_vars.u(), qp, u_dot);
358  context.interior_rate(this->_flow_vars.v(), qp, v_dot);
359 
360  if( this->_dim == 3 )
361  context.interior_rate(this->_flow_vars.w(), qp, w_dot);
362 
363  for (unsigned int i = 0; i != n_u_dofs; ++i)
364  {
365  F_u(i) -= JxW[qp]*this->_rho*u_dot*u_phi[i][qp];
366  F_v(i) -= JxW[qp]*this->_rho*v_dot*u_phi[i][qp];
367 
368  if( this->_dim == 3 )
369  (*F_w)(i) -= JxW[qp]*this->_rho*w_dot*u_phi[i][qp];
370 
371  if( compute_jacobian )
372  {
373  for (unsigned int j=0; j != n_u_dofs; j++)
374  {
375  // Assuming rho is constant w.r.t. u, v, w
376  // and T (if Boussinesq added).
377  libMesh::Real value = context.get_elem_solution_derivative() * JxW[qp]*this->_rho*u_phi[i][qp]*u_phi[j][qp];
378 
379  M_uu(i,j) -= value;
380  M_vv(i,j) -= value;
381 
382  if( this->_dim == 3)
383  {
384  (*M_ww)(i,j) -= value;
385  }
386 
387  } // End dof loop
388  } // End Jacobian check
389  } // End dof loop
390  } // End quadrature loop
391 
392  return;
393  }
394 
395 } // namespace GRINS
396 
397 // Instantiate
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:269
Physics class for Incompressible Navier-Stokes.
INSTANTIATE_INC_NS_SUBCLASS(Stokes)
Base class for reading and handling initial conditions for physics classes.
GRINS namespace.
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...
Definition: stokes.C:315
virtual void auxiliary_init(MultiphysicsSystem &system)
Any auxillary initialization a Physics class may need.
Definition: stokes.C:61
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Time dependent part(s) of physics for element interiors.
Definition: stokes.C:68
Interface with libMesh for solving Multiphysics problems.
virtual void element_constraint(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Constraint part(s) of physics for element interiors.
Definition: stokes.C:220

Generated on Thu Jun 2 2016 21:52:28 for GRINS-0.7.0 by  doxygen 1.8.10