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

Generated on Mon Jun 22 2015 21:32:20 for GRINS-0.6.0 by  doxygen 1.8.9.1