libMesh
systems_of_equations_ex3.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2017 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 // <h1>Systems Example 3 - Navier-Stokes with SCALAR Lagrange Multiplier</h1>
21 // \author David Knezevic
22 // \date 2010
23 //
24 // This example shows how the transient Navier-Stokes problem from
25 // example 13 can be solved using a scalar Lagrange multiplier
26 // formulation to constrain the integral of the pressure variable,
27 // rather than pinning the pressure at a single point.
28 
29 // C++ include files that we need
30 #include <iostream>
31 #include <algorithm>
32 #include <sstream>
33 #include <math.h>
34 
35 // Basic include file needed for the mesh functionality.
36 #include "libmesh/libmesh.h"
37 #include "libmesh/mesh.h"
38 #include "libmesh/mesh_generation.h"
39 #include "libmesh/exodusII_io.h"
40 #include "libmesh/equation_systems.h"
41 #include "libmesh/fe.h"
42 #include "libmesh/quadrature_gauss.h"
43 #include "libmesh/dof_map.h"
44 #include "libmesh/sparse_matrix.h"
45 #include "libmesh/numeric_vector.h"
46 #include "libmesh/dense_matrix.h"
47 #include "libmesh/dense_vector.h"
48 #include "libmesh/linear_implicit_system.h"
49 #include "libmesh/transient_system.h"
50 #include "libmesh/perf_log.h"
51 #include "libmesh/boundary_info.h"
52 #include "libmesh/utility.h"
53 #include "libmesh/dirichlet_boundaries.h"
54 #include "libmesh/zero_function.h"
55 #include "libmesh/const_function.h"
56 
57 // For systems of equations the DenseSubMatrix
58 // and DenseSubVector provide convenient ways for
59 // assembling the element matrix and vector on a
60 // component-by-component basis.
61 #include "libmesh/dense_submatrix.h"
62 #include "libmesh/dense_subvector.h"
63 
64 // The definition of a geometric element
65 #include "libmesh/elem.h"
66 
67 // Bring in everything from the libMesh namespace
68 using namespace libMesh;
69 
70 // Function prototype. This function will assemble the system
71 // matrix and right-hand-side.
73  const std::string & system_name);
74 
75 // Function which sets Dirichlet BCs for the lid-driven cavity.
77 
78 // The main program.
79 int main (int argc, char ** argv)
80 {
81  // Initialize libMesh.
82  LibMeshInit init (argc, argv);
83 
84  // This example requires a linear solver package.
85  libmesh_example_requires(libMesh::default_solver_package() != INVALID_SOLVER_PACKAGE,
86  "--enable-petsc, --enable-trilinos, or --enable-eigen");
87 
88  // Skip this 2D example if libMesh was compiled as 1D-only.
89  libmesh_example_requires(2 <= LIBMESH_DIM, "2D support");
90 
91  // This example NaNs with the Eigen sparse linear solvers and
92  // Trilinos solvers, but should work OK with either PETSc or
93  // Laspack.
94  libmesh_example_requires(libMesh::default_solver_package() != EIGEN_SOLVERS, "--enable-petsc or --enable-laspack");
95  libmesh_example_requires(libMesh::default_solver_package() != TRILINOS_SOLVERS, "--enable-petsc or --enable-laspack");
96 
97  // Create a mesh, with dimension to be overridden later, distributed
98  // across the default MPI communicator.
99  Mesh mesh(init.comm());
100 
101  // Use the MeshTools::Generation mesh generator to create a uniform
102  // 2D grid on the square [-1,1]^2. We instruct the mesh generator
103  // to build a mesh of 8x8 Quad9 elements in 2D. Building these
104  // higher-order elements allows us to use higher-order
105  // approximation, as in example 3.
107  20, 20,
108  0., 1.,
109  0., 1.,
110  QUAD9);
111 
112  // Print information about the mesh to the screen.
113  mesh.print_info();
114 
115  // Create an equation systems object.
116  EquationSystems equation_systems (mesh);
117 
118  // Declare the system and its variables.
119  // Creates a transient system named "Navier-Stokes"
121  equation_systems.add_system<TransientLinearImplicitSystem> ("Navier-Stokes");
122 
123  // Add the variables "vel_x" & "vel_y" to "Navier-Stokes". They
124  // will be approximated using second-order approximation.
125  system.add_variable ("vel_x", SECOND);
126  system.add_variable ("vel_y", SECOND);
127 
128  // Add the variable "p" to "Navier-Stokes". This will
129  // be approximated with a first-order basis,
130  // providing an LBB-stable pressure-velocity pair.
131  system.add_variable ("p", FIRST);
132 
133  // Add a scalar Lagrange multiplier to constrain the
134  // pressure to have zero mean.
135  system.add_variable ("alpha", FIRST, SCALAR);
136 
137  // Give the system a pointer to the matrix assembly
138  // function.
139  system.attach_assemble_function (assemble_stokes);
140 
141  // Set Dirichlet boundary conditions.
142  set_lid_driven_bcs(system);
143 
144  // Initialize the data structures for the equation system.
145  equation_systems.init ();
146 
147  // Prints information about the system to the screen.
148  equation_systems.print_info();
149 
150  // Create a performance-logging object for this example
151  PerfLog perf_log("Systems Example 3");
152 
153  // Get a reference to the Stokes system to use later.
154  TransientLinearImplicitSystem & navier_stokes_system =
155  equation_systems.get_system<TransientLinearImplicitSystem>("Navier-Stokes");
156 
157  // Now we begin the timestep loop to compute the time-accurate
158  // solution of the equations.
159  const Real dt = 0.1;
160  navier_stokes_system.time = 0.0;
161  const unsigned int n_timesteps = 15;
162 
163  // The number of steps and the stopping criterion are also required
164  // for the nonlinear iterations.
165  const unsigned int n_nonlinear_steps = 15;
166  const Real nonlinear_tolerance = TOLERANCE*10;
167 
168  // We also set a standard linear solver flag in the EquationSystems object
169  // which controls the maximum number of linear solver iterations allowed.
170  equation_systems.parameters.set<unsigned int>("linear solver maximum iterations") = 250;
171 
172  // Tell the system of equations what the timestep is by using
173  // the set_parameter function. The matrix assembly routine can
174  // then reference this parameter.
175  equation_systems.parameters.set<Real> ("dt") = dt;
176 
177  // The kinematic viscosity, nu = mu/rho, units of length**2/time.
178  equation_systems.parameters.set<Real> ("nu") = .007;
179 
180  // The first thing to do is to get a copy of the solution at
181  // the current nonlinear iteration. This value will be used to
182  // determine if we can exit the nonlinear loop.
184  last_nonlinear_soln (navier_stokes_system.solution->clone());
185 
186 #ifdef LIBMESH_HAVE_EXODUS_API
187  // Since we are not doing adaptivity, write all solutions to a single Exodus file.
188  ExodusII_IO exo_io(mesh);
189 
190  // Write out the initial condition
191  exo_io.write_equation_systems ("out.e", equation_systems);
192 #endif
193 
194  for (unsigned int t_step=1; t_step<=n_timesteps; ++t_step)
195  {
196  // Increment the time counter, set the time and the
197  // time step size as parameters in the EquationSystem.
198  navier_stokes_system.time += dt;
199 
200  // A pretty update message
201  libMesh::out << "\n\n*** Solving time step "
202  << t_step
203  << ", time = "
204  << navier_stokes_system.time
205  << " ***"
206  << std::endl;
207 
208  // Now we need to update the solution vector from the
209  // previous time step. This is done directly through
210  // the reference to the Stokes system.
211  *navier_stokes_system.old_local_solution = *navier_stokes_system.current_local_solution;
212 
213  // At the beginning of each solve, reset the linear solver tolerance
214  // to a "reasonable" starting value.
215  const Real initial_linear_solver_tol = 1.e-6;
216  equation_systems.parameters.set<Real> ("linear solver tolerance") = initial_linear_solver_tol;
217 
218  // We'll set this flag when convergence is (hopefully) achieved.
219  bool converged = false;
220 
221  // Now we begin the nonlinear loop
222  for (unsigned int l=0; l<n_nonlinear_steps; ++l)
223  {
224  // Update the nonlinear solution.
225  last_nonlinear_soln->zero();
226  last_nonlinear_soln->add(*navier_stokes_system.solution);
227 
228  // Assemble & solve the linear system.
229  perf_log.push("linear solve");
230  equation_systems.get_system("Navier-Stokes").solve();
231  perf_log.pop("linear solve");
232 
233  // Compute the difference between this solution and the last
234  // nonlinear iterate.
235  last_nonlinear_soln->add (-1., *navier_stokes_system.solution);
236 
237  // Close the vector before computing its norm
238  last_nonlinear_soln->close();
239 
240  // Compute the l2 norm of the difference
241  const Real norm_delta = last_nonlinear_soln->l2_norm();
242 
243  // How many iterations were required to solve the linear system?
244  const unsigned int n_linear_iterations = navier_stokes_system.n_linear_iterations();
245 
246  // What was the final residual of the linear system?
247  const Real final_linear_residual = navier_stokes_system.final_linear_residual();
248 
249  // If the solver did no work (sometimes -ksp_converged_reason
250  // says "Linear solve converged due to CONVERGED_RTOL
251  // iterations 0") but the nonlinear residual norm is above
252  // the tolerance, we need to pick an even lower linear
253  // solver tolerance and try again. Note that the tolerance
254  // is relative to the norm of the RHS, which for this
255  // particular problem does not go to zero, since we are
256  // solving for the full solution rather than the update.
257  //
258  // Similarly, if the solver did no work and this is the 0th
259  // nonlinear step, it means that the delta between solutions
260  // is being inaccurately measured as "0" since the solution
261  // did not change. Decrease the tolerance and try again.
262  if (n_linear_iterations == 0 &&
263  (navier_stokes_system.final_linear_residual() >= nonlinear_tolerance || l==0))
264  {
265  Real old_linear_solver_tolerance = equation_systems.parameters.get<Real> ("linear solver tolerance");
266  equation_systems.parameters.set<Real> ("linear solver tolerance") = 1.e-3 * old_linear_solver_tolerance;
267  continue;
268  }
269 
270  // Print out convergence information for the linear and
271  // nonlinear iterations.
272  libMesh::out << "Linear solver converged at step: "
273  << n_linear_iterations
274  << ", final residual: "
275  << final_linear_residual
276  << " Nonlinear convergence: ||u - u_old|| = "
277  << norm_delta
278  << std::endl;
279 
280  // Terminate the solution iteration if the difference between
281  // this nonlinear iterate and the last is sufficiently small, AND
282  // if the most recent linear system was solved to a sufficient tolerance.
283  if ((norm_delta < nonlinear_tolerance) &&
284  (navier_stokes_system.final_linear_residual() < nonlinear_tolerance))
285  {
286  libMesh::out << " Nonlinear solver converged at step "
287  << l
288  << std::endl;
289  converged = true;
290  break;
291  }
292 
293  // Otherwise, decrease the linear system tolerance. For the inexact Newton
294  // method, the linear solver tolerance needs to decrease as we get closer to
295  // the solution to ensure quadratic convergence. The new linear solver tolerance
296  // is chosen (heuristically) as the square of the previous linear system residual norm.
297  //Real flr2 = final_linear_residual*final_linear_residual;
298  Real new_linear_solver_tolerance = std::min(Utility::pow<2>(final_linear_residual), initial_linear_solver_tol);
299  equation_systems.parameters.set<Real> ("linear solver tolerance") = new_linear_solver_tolerance;
300  } // end nonlinear loop
301 
302  // Don't keep going if we failed to converge.
303  if (!converged)
304  libmesh_error_msg("Error: Newton iterations failed to converge!");
305 
306 #ifdef LIBMESH_HAVE_EXODUS_API
307  // Write out every nth timestep to file.
308  const unsigned int write_interval = 1;
309 
310  if ((t_step+1)%write_interval == 0)
311  {
312  exo_io.write_timestep("out.e",
313  equation_systems,
314  t_step+1, // we're off by one since we wrote the IC and the Exodus numbering is 1-based.
315  navier_stokes_system.time);
316  }
317 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
318  } // end timestep loop.
319 
320  // All done.
321  return 0;
322 }
323 
324 
325 
326 
327 
328 
329 // The matrix assembly function to be called at each time step to
330 // prepare for the linear solve.
332  const std::string & libmesh_dbg_var(system_name))
333 {
334  // It is a good idea to make sure we are assembling
335  // the proper system.
336  libmesh_assert_equal_to (system_name, "Navier-Stokes");
337 
338  // Get a constant reference to the mesh object.
339  const MeshBase & mesh = es.get_mesh();
340 
341  // The dimension that we are running
342  const unsigned int dim = mesh.mesh_dimension();
343 
344  // Get a reference to the Stokes system object.
345  TransientLinearImplicitSystem & navier_stokes_system =
346  es.get_system<TransientLinearImplicitSystem> ("Navier-Stokes");
347 
348  // Numeric ids corresponding to each variable in the system
349  const unsigned int u_var = navier_stokes_system.variable_number ("vel_x");
350  const unsigned int v_var = navier_stokes_system.variable_number ("vel_y");
351  const unsigned int p_var = navier_stokes_system.variable_number ("p");
352  const unsigned int alpha_var = navier_stokes_system.variable_number ("alpha");
353 
354  // Get the Finite Element type for "vel_x". Note this will be
355  // the same as the type for "vel_y".
356  FEType fe_vel_type = navier_stokes_system.variable_type(u_var);
357 
358  // Get the Finite Element type for "p".
359  FEType fe_pres_type = navier_stokes_system.variable_type(p_var);
360 
361  // Build a Finite Element object of the specified type for
362  // the velocity variables.
363  UniquePtr<FEBase> fe_vel (FEBase::build(dim, fe_vel_type));
364 
365  // Build a Finite Element object of the specified type for
366  // the pressure variables.
367  UniquePtr<FEBase> fe_pres (FEBase::build(dim, fe_pres_type));
368 
369  // A Gauss quadrature rule for numerical integration.
370  // Let the FEType object decide what order rule is appropriate.
371  QGauss qrule (dim, fe_vel_type.default_quadrature_order());
372 
373  // Tell the finite element objects to use our quadrature rule.
374  fe_vel->attach_quadrature_rule (&qrule);
375  fe_pres->attach_quadrature_rule (&qrule);
376 
377  // Here we define some references to cell-specific data that
378  // will be used to assemble the linear system.
379  //
380  // The element Jacobian * quadrature weight at each integration point.
381  const std::vector<Real> & JxW = fe_vel->get_JxW();
382 
383  // The element shape functions evaluated at the quadrature points.
384  const std::vector<std::vector<Real>> & phi = fe_vel->get_phi();
385 
386  // The element shape function gradients for the velocity
387  // variables evaluated at the quadrature points.
388  const std::vector<std::vector<RealGradient>> & dphi = fe_vel->get_dphi();
389 
390  // The element shape functions for the pressure variable
391  // evaluated at the quadrature points.
392  const std::vector<std::vector<Real>> & psi = fe_pres->get_phi();
393 
394  // The value of the linear shape function gradients at the quadrature points
395  // const std::vector<std::vector<RealGradient>> & dpsi = fe_pres->get_dphi();
396 
397  // A reference to the DofMap object for this system. The DofMap
398  // object handles the index translation from node and element numbers
399  // to degree of freedom numbers. We will talk more about the DofMap
400  // in future examples.
401  const DofMap & dof_map = navier_stokes_system.get_dof_map();
402 
403  // Define data structures to contain the element matrix
404  // and right-hand-side vector contribution. Following
405  // basic finite element terminology we will denote these
406  // "Ke" and "Fe".
409 
411  Kuu(Ke), Kuv(Ke), Kup(Ke),
412  Kvu(Ke), Kvv(Ke), Kvp(Ke),
413  Kpu(Ke), Kpv(Ke), Kpp(Ke);
414  DenseSubMatrix<Number> Kalpha_p(Ke), Kp_alpha(Ke);
415 
417  Fu(Fe),
418  Fv(Fe),
419  Fp(Fe);
420 
421  // This vector will hold the degree of freedom indices for
422  // the element. These define where in the global system
423  // the element degrees of freedom get mapped.
424  std::vector<dof_id_type> dof_indices;
425  std::vector<dof_id_type> dof_indices_u;
426  std::vector<dof_id_type> dof_indices_v;
427  std::vector<dof_id_type> dof_indices_p;
428  std::vector<dof_id_type> dof_indices_alpha;
429 
430  // Find out what the timestep size parameter is from the system, and
431  // the value of theta for the theta method. We use implicit Euler (theta=1)
432  // for this simulation even though it is only first-order accurate in time.
433  // The reason for this decision is that the second-order Crank-Nicolson
434  // method is notoriously oscillatory for problems with discontinuous
435  // initial data such as the lid-driven cavity. Therefore,
436  // we sacrifice accuracy in time for stability, but since the solution
437  // reaches steady state relatively quickly we can afford to take small
438  // timesteps. If you monitor the initial nonlinear residual for this
439  // simulation, you should see that it is monotonically decreasing in time.
440  const Real dt = es.parameters.get<Real>("dt");
441  const Real theta = 1.;
442 
443  // The kinematic viscosity, multiplies the "viscous" terms.
444  const Real nu = es.parameters.get<Real>("nu");
445 
446  // Now we will loop over all the elements in the mesh that
447  // live on the local processor. We will compute the element
448  // matrix and right-hand-side contribution. Since the mesh
449  // will be refined we want to only consider the ACTIVE elements,
450  // hence we use a variant of the active_elem_iterator.
451  for (const auto & elem : mesh.active_local_element_ptr_range())
452  {
453  // Get the degree of freedom indices for the
454  // current element. These define where in the global
455  // matrix and right-hand-side this element will
456  // contribute to.
457  dof_map.dof_indices (elem, dof_indices);
458  dof_map.dof_indices (elem, dof_indices_u, u_var);
459  dof_map.dof_indices (elem, dof_indices_v, v_var);
460  dof_map.dof_indices (elem, dof_indices_p, p_var);
461  dof_map.dof_indices (elem, dof_indices_alpha, alpha_var);
462 
463  const unsigned int n_dofs = dof_indices.size();
464  const unsigned int n_u_dofs = dof_indices_u.size();
465  const unsigned int n_v_dofs = dof_indices_v.size();
466  const unsigned int n_p_dofs = dof_indices_p.size();
467 
468  // Compute the element-specific data for the current
469  // element. This involves computing the location of the
470  // quadrature points (q_point) and the shape functions
471  // (phi, dphi) for the current element.
472  fe_vel->reinit (elem);
473  fe_pres->reinit (elem);
474 
475  // Zero the element matrix and right-hand side before
476  // summing them. We use the resize member here because
477  // the number of degrees of freedom might have changed from
478  // the last element. Note that this will be the case if the
479  // element type is different (i.e. the last element was a
480  // triangle, now we are on a quadrilateral).
481  Ke.resize (n_dofs, n_dofs);
482  Fe.resize (n_dofs);
483 
484  // Reposition the submatrices... The idea is this:
485  //
486  // - - - -
487  // | Kuu Kuv Kup | | Fu |
488  // Ke = | Kvu Kvv Kvp |; Fe = | Fv |
489  // | Kpu Kpv Kpp | | Fp |
490  // - - - -
491  //
492  // The DenseSubMatrix.reposition () member takes the
493  // (row_offset, column_offset, row_size, column_size).
494  //
495  // Similarly, the DenseSubVector.reposition () member
496  // takes the (row_offset, row_size)
497  Kuu.reposition (u_var*n_u_dofs, u_var*n_u_dofs, n_u_dofs, n_u_dofs);
498  Kuv.reposition (u_var*n_u_dofs, v_var*n_u_dofs, n_u_dofs, n_v_dofs);
499  Kup.reposition (u_var*n_u_dofs, p_var*n_u_dofs, n_u_dofs, n_p_dofs);
500 
501  Kvu.reposition (v_var*n_v_dofs, u_var*n_v_dofs, n_v_dofs, n_u_dofs);
502  Kvv.reposition (v_var*n_v_dofs, v_var*n_v_dofs, n_v_dofs, n_v_dofs);
503  Kvp.reposition (v_var*n_v_dofs, p_var*n_v_dofs, n_v_dofs, n_p_dofs);
504 
505  Kpu.reposition (p_var*n_u_dofs, u_var*n_u_dofs, n_p_dofs, n_u_dofs);
506  Kpv.reposition (p_var*n_u_dofs, v_var*n_u_dofs, n_p_dofs, n_v_dofs);
507  Kpp.reposition (p_var*n_u_dofs, p_var*n_u_dofs, n_p_dofs, n_p_dofs);
508 
509  // Also, add a row and a column to constrain the pressure
510  Kp_alpha.reposition (p_var*n_u_dofs, p_var*n_u_dofs+n_p_dofs, n_p_dofs, 1);
511  Kalpha_p.reposition (p_var*n_u_dofs+n_p_dofs, p_var*n_u_dofs, 1, n_p_dofs);
512 
513 
514  Fu.reposition (u_var*n_u_dofs, n_u_dofs);
515  Fv.reposition (v_var*n_u_dofs, n_v_dofs);
516  Fp.reposition (p_var*n_u_dofs, n_p_dofs);
517 
518  // Now we will build the element matrix and right-hand-side.
519  // Constructing the RHS requires the solution and its
520  // gradient from the previous timestep. This must be
521  // calculated at each quadrature point by summing the
522  // solution degree-of-freedom values by the appropriate
523  // weight functions.
524  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
525  {
526  // Values to hold the solution & its gradient at the previous timestep.
527  Number u = 0., u_old = 0.;
528  Number v = 0., v_old = 0.;
529  Number p_old = 0.;
530  Gradient grad_u, grad_u_old;
531  Gradient grad_v, grad_v_old;
532 
533  // Compute the velocity & its gradient from the previous timestep
534  // and the old Newton iterate.
535  for (unsigned int l=0; l<n_u_dofs; l++)
536  {
537  // From the old timestep:
538  u_old += phi[l][qp]*navier_stokes_system.old_solution (dof_indices_u[l]);
539  v_old += phi[l][qp]*navier_stokes_system.old_solution (dof_indices_v[l]);
540  grad_u_old.add_scaled (dphi[l][qp], navier_stokes_system.old_solution (dof_indices_u[l]));
541  grad_v_old.add_scaled (dphi[l][qp], navier_stokes_system.old_solution (dof_indices_v[l]));
542 
543  // From the previous Newton iterate:
544  u += phi[l][qp]*navier_stokes_system.current_solution (dof_indices_u[l]);
545  v += phi[l][qp]*navier_stokes_system.current_solution (dof_indices_v[l]);
546  grad_u.add_scaled (dphi[l][qp], navier_stokes_system.current_solution (dof_indices_u[l]));
547  grad_v.add_scaled (dphi[l][qp], navier_stokes_system.current_solution (dof_indices_v[l]));
548  }
549 
550  // Compute the old pressure value at this quadrature point.
551  for (unsigned int l=0; l<n_p_dofs; l++)
552  p_old += psi[l][qp]*navier_stokes_system.old_solution (dof_indices_p[l]);
553 
554  // Definitions for convenience. It is sometimes simpler to do a
555  // dot product if you have the full vector at your disposal.
556  const NumberVectorValue U_old (u_old, v_old);
557  const NumberVectorValue U (u, v);
558  const Number u_x = grad_u(0);
559  const Number u_y = grad_u(1);
560  const Number v_x = grad_v(0);
561  const Number v_y = grad_v(1);
562 
563  // First, an i-loop over the velocity degrees of freedom.
564  // We know that n_u_dofs == n_v_dofs so we can compute contributions
565  // for both at the same time.
566  for (unsigned int i=0; i<n_u_dofs; i++)
567  {
568  Fu(i) += JxW[qp]*(u_old*phi[i][qp] - // mass-matrix term
569  (1.-theta)*dt*(U_old*grad_u_old)*phi[i][qp] + // convection term
570  (1.-theta)*dt*p_old*dphi[i][qp](0) - // pressure term on rhs
571  (1.-theta)*dt*nu*(grad_u_old*dphi[i][qp]) + // diffusion term on rhs
572  theta*dt*(U*grad_u)*phi[i][qp]); // Newton term
573 
574 
575  Fv(i) += JxW[qp]*(v_old*phi[i][qp] - // mass-matrix term
576  (1.-theta)*dt*(U_old*grad_v_old)*phi[i][qp] + // convection term
577  (1.-theta)*dt*p_old*dphi[i][qp](1) - // pressure term on rhs
578  (1.-theta)*dt*nu*(grad_v_old*dphi[i][qp]) + // diffusion term on rhs
579  theta*dt*(U*grad_v)*phi[i][qp]); // Newton term
580 
581 
582  // Note that the Fp block is identically zero unless we are using
583  // some kind of artificial compressibility scheme...
584 
585  // Matrix contributions for the uu and vv couplings.
586  for (unsigned int j=0; j<n_u_dofs; j++)
587  {
588  Kuu(i,j) += JxW[qp]*(phi[i][qp]*phi[j][qp] + // mass matrix term
589  theta*dt*nu*(dphi[i][qp]*dphi[j][qp]) + // diffusion term
590  theta*dt*(U*dphi[j][qp])*phi[i][qp] + // convection term
591  theta*dt*u_x*phi[i][qp]*phi[j][qp]); // Newton term
592 
593  Kuv(i,j) += JxW[qp]*theta*dt*u_y*phi[i][qp]*phi[j][qp]; // Newton term
594 
595  Kvv(i,j) += JxW[qp]*(phi[i][qp]*phi[j][qp] + // mass matrix term
596  theta*dt*nu*(dphi[i][qp]*dphi[j][qp]) + // diffusion term
597  theta*dt*(U*dphi[j][qp])*phi[i][qp] + // convection term
598  theta*dt*v_y*phi[i][qp]*phi[j][qp]); // Newton term
599 
600  Kvu(i,j) += JxW[qp]*theta*dt*v_x*phi[i][qp]*phi[j][qp]; // Newton term
601  }
602 
603  // Matrix contributions for the up and vp couplings.
604  for (unsigned int j=0; j<n_p_dofs; j++)
605  {
606  Kup(i,j) += JxW[qp]*(-theta*dt*psi[j][qp]*dphi[i][qp](0));
607  Kvp(i,j) += JxW[qp]*(-theta*dt*psi[j][qp]*dphi[i][qp](1));
608  }
609  }
610 
611  // Now an i-loop over the pressure degrees of freedom. This code computes
612  // the matrix entries due to the continuity equation. Note: To maintain a
613  // symmetric matrix, we may (or may not) multiply the continuity equation by
614  // negative one. Here we do not.
615  for (unsigned int i=0; i<n_p_dofs; i++)
616  {
617  Kp_alpha(i,0) += JxW[qp]*psi[i][qp];
618  Kalpha_p(0,i) += JxW[qp]*psi[i][qp];
619  for (unsigned int j=0; j<n_u_dofs; j++)
620  {
621  Kpu(i,j) += JxW[qp]*psi[i][qp]*dphi[j][qp](0);
622  Kpv(i,j) += JxW[qp]*psi[i][qp]*dphi[j][qp](1);
623  }
624  }
625  } // end of the quadrature point qp-loop
626 
627  // Since we're using heterogeneous DirichletBoundary objects for
628  // the boundary conditions, we need to call a specific function
629  // to constrain the element stiffness matrix.
630  dof_map.heterogenously_constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
631 
632  // The element matrix and right-hand-side are now built
633  // for this element. Add them to the global matrix and
634  // right-hand-side vector. The SparseMatrix::add_matrix()
635  // and NumericVector::add_vector() members do this for us.
636  navier_stokes_system.matrix->add_matrix (Ke, dof_indices);
637  navier_stokes_system.rhs->add_vector (Fe, dof_indices);
638  } // end of element loop
639 
640  // We can set the mean of the pressure by setting Falpha. Typically
641  // a value of zero is chosen, but the value should be arbitrary.
642  navier_stokes_system.rhs->add(navier_stokes_system.rhs->size()-1, 10.);
643 }
644 
645 
646 
648 {
649  unsigned short int
650  u_var = system.variable_number("vel_x"),
651  v_var = system.variable_number("vel_y");
652 
653  // Get a convenient reference to the System's DofMap
654  DofMap & dof_map = system.get_dof_map();
655 
656  {
657  // u=v=0 on bottom, left, right
658  std::set<boundary_id_type> boundary_ids;
659  boundary_ids.insert(0);
660  boundary_ids.insert(1);
661  boundary_ids.insert(3);
662 
663  std::vector<unsigned int> variables;
664  variables.push_back(u_var);
665  variables.push_back(v_var);
666 
667  dof_map.add_dirichlet_boundary(DirichletBoundary(boundary_ids,
668  variables,
670  }
671  {
672  // u=1 on top
673  std::set<boundary_id_type> boundary_ids;
674  boundary_ids.insert(2);
675 
676  std::vector<unsigned int> variables;
677  variables.push_back(u_var);
678 
679  dof_map.add_dirichlet_boundary(DirichletBoundary(boundary_ids,
680  variables,
681  ConstFunction<Number>(1.)));
682  }
683  {
684  // v=0 on top
685  std::set<boundary_id_type> boundary_ids;
686  boundary_ids.insert(2);
687 
688  std::vector<unsigned int> variables;
689  variables.push_back(v_var);
690 
691  dof_map.add_dirichlet_boundary(DirichletBoundary(boundary_ids,
692  variables,
694  }
695 }
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
This is the EquationSystems class.
void pop(const char *label, const char *header="")
Pop the event label off the stack, resuming any lower event.
Definition: perf_log.C:162
UniquePtr< NumericVector< Number > > old_local_solution
All the values I need to compute my contribution to the simulation at hand.
ConstFunction that simply returns 0.
Definition: zero_function.h:35
void add_scaled(const TypeVector< T2 > &, const T)
Add a scaled value to this vector without creating a temporary.
Definition: type_vector.h:624
unsigned int dim
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Definition: exodusII_io.h:52
void resize(const unsigned int n)
Resize the vector.
Definition: dense_vector.h:350
MeshBase & mesh
This class allows one to associate Dirichlet boundary values with a given set of mesh boundary ids an...
This class provides a specific system class.
const T & get(const std::string &) const
Definition: parameters.h:430
void build_square(UnstructuredMesh &mesh, const unsigned int nx, const unsigned int ny, const Real xmin=0., const Real xmax=1., const Real ymin=0., const Real ymax=1., const ElemType type=INVALID_ELEM, const bool gauss_lobatto_grid=false)
A specialized build_cube() for 2D meshes.
static const Real TOLERANCE
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:62
The libMesh namespace provides an interface to certain functionality in the library.
This is the MeshBase class.
Definition: mesh_base.h:68
std::unique_ptr< T > UniquePtr
Definition: auto_ptr.h:46
Defines a dense subvector for use in finite element computations.
void set_lid_driven_bcs(TransientLinearImplicitSystem &system)
SolverPackage default_solver_package()
Definition: libmesh.C:995
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
The PerfLog class allows monitoring of specific events.
Definition: perf_log.h:125
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:167
void reposition(const unsigned int ioff, const unsigned int n)
Changes the location of the subvector in the parent vector.
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
Order default_quadrature_order() const
Definition: fe_type.h:332
Defines a dense submatrix for use in Finite Element-type computations.
void init(triangulateio &t)
Initializes the fields of t to NULL/0 as necessary.
virtual System & add_system(const std::string &system_type, const std::string &name)
Add the system of type system_type named name to the systems array.
void push(const char *label, const char *header="")
Push the event label onto the stack, pausing any active event.
Definition: perf_log.C:132
const Parallel::Communicator & comm() const
Definition: libmesh.h:81
PetscErrorCode Vec Mat libmesh_dbg_var(j)
void assemble_stokes(EquationSystems &es, const std::string &system_name)
void heterogenously_constrain_element_matrix_and_vector(DenseMatrix< Number > &matrix, DenseVector< Number > &rhs, std::vector< dof_id_type > &elem_dofs, bool asymmetric_constraint_rows=true, int qoi_index=-1) const
Constrains the element matrix and vector.
int main(int argc, char **argv)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
T & set(const std::string &)
Definition: parameters.h:469
OStreamProxy out
unsigned int mesh_dimension() const
Definition: mesh_base.C:148
void resize(const unsigned int new_m, const unsigned int new_n)
Resize the matrix.
Definition: dense_matrix.h:776
This class implements specific orders of Gauss quadrature.
Parameters parameters
Data structure holding arbitrary parameters.
Function that returns a single value that never changes.
const MeshBase & get_mesh() const
void add_dirichlet_boundary(const DirichletBoundary &dirichlet_boundary)
Adds a copy of the specified Dirichlet boundary to the system.
Number old_solution(const dof_id_type global_dof_number) const
void reposition(const unsigned int ioff, const unsigned int joff, const unsigned int new_m, const unsigned int new_n)
Changes the location of the submatrix in the parent matrix.
const T_sys & get_system(const std::string &name) const
virtual void init()
Initialize all the systems.
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50
long double min(long double a, double b)
void print_info(std::ostream &os=libMesh::out) const
Prints relevant information about the mesh.
Definition: mesh_base.C:448
void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
Fills the vector di with the global degree of freedom indices for the element.
Definition: dof_map.C:1917
static UniquePtr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.