libMesh
Functions
systems_of_equations_ex3.C File Reference

Go to the source code of this file.

Functions

void assemble_stokes (EquationSystems &es, const std::string &system_name)
 
void set_lid_driven_bcs (TransientLinearImplicitSystem &system)
 
int main (int argc, char **argv)
 
void assemble_stokes (EquationSystems &es, const std::string &libmesh_dbg_var(system_name))
 

Function Documentation

◆ assemble_stokes() [1/2]

void assemble_stokes ( EquationSystems es,
const std::string &  libmesh_dbg_varsystem_name 
)

Definition at line 337 of file systems_of_equations_ex3.C.

339 {
340  // It is a good idea to make sure we are assembling
341  // the proper system.
342  libmesh_assert_equal_to (system_name, "Navier-Stokes");
343 
344 #if LIBMESH_DIM > 1
345  // Get a constant reference to the mesh object.
346  const MeshBase & mesh = es.get_mesh();
347 
348  // The dimension that we are running
349  const unsigned int dim = mesh.mesh_dimension();
350 
351  // Get a reference to the Stokes system object.
352  TransientLinearImplicitSystem & navier_stokes_system =
353  es.get_system<TransientLinearImplicitSystem> ("Navier-Stokes");
354 
355  // Numeric ids corresponding to each variable in the system
356  const unsigned int u_var = navier_stokes_system.variable_number ("vel_x");
357  const unsigned int v_var = navier_stokes_system.variable_number ("vel_y");
358  const unsigned int p_var = navier_stokes_system.variable_number ("p");
359  const unsigned int alpha_var = navier_stokes_system.variable_number ("alpha");
360 
361  // Get the Finite Element type for "vel_x". Note this will be
362  // the same as the type for "vel_y".
363  FEType fe_vel_type = navier_stokes_system.variable_type(u_var);
364 
365  // Get the Finite Element type for "p".
366  FEType fe_pres_type = navier_stokes_system.variable_type(p_var);
367 
368  // Build a Finite Element object of the specified type for
369  // the velocity variables.
370  std::unique_ptr<FEBase> fe_vel (FEBase::build(dim, fe_vel_type));
371 
372  // Build a Finite Element object of the specified type for
373  // the pressure variables.
374  std::unique_ptr<FEBase> fe_pres (FEBase::build(dim, fe_pres_type));
375 
376  // A Gauss quadrature rule for numerical integration.
377  // Let the FEType object decide what order rule is appropriate.
378  QGauss qrule (dim, fe_vel_type.default_quadrature_order());
379 
380  // Tell the finite element objects to use our quadrature rule.
381  fe_vel->attach_quadrature_rule (&qrule);
382  fe_pres->attach_quadrature_rule (&qrule);
383 
384  // Here we define some references to cell-specific data that
385  // will be used to assemble the linear system.
386  //
387  // The element Jacobian * quadrature weight at each integration point.
388  const std::vector<Real> & JxW = fe_vel->get_JxW();
389 
390  // The element shape functions evaluated at the quadrature points.
391  const std::vector<std::vector<Real>> & phi = fe_vel->get_phi();
392 
393  // The element shape function gradients for the velocity
394  // variables evaluated at the quadrature points.
395  const std::vector<std::vector<RealGradient>> & dphi = fe_vel->get_dphi();
396 
397  // The element shape functions for the pressure variable
398  // evaluated at the quadrature points.
399  const std::vector<std::vector<Real>> & psi = fe_pres->get_phi();
400 
401  // The value of the linear shape function gradients at the quadrature points
402  // const std::vector<std::vector<RealGradient>> & dpsi = fe_pres->get_dphi();
403 
404  // A reference to the DofMap object for this system. The DofMap
405  // object handles the index translation from node and element numbers
406  // to degree of freedom numbers. We will talk more about the DofMap
407  // in future examples.
408  const DofMap & dof_map = navier_stokes_system.get_dof_map();
409 
410  // Define data structures to contain the element matrix
411  // and right-hand-side vector contribution. Following
412  // basic finite element terminology we will denote these
413  // "Ke" and "Fe".
416 
418  Kuu(Ke), Kuv(Ke), Kup(Ke),
419  Kvu(Ke), Kvv(Ke), Kvp(Ke),
420  Kpu(Ke), Kpv(Ke), Kpp(Ke);
421  DenseSubMatrix<Number> Kalpha_p(Ke), Kp_alpha(Ke);
422 
424  Fu(Fe),
425  Fv(Fe),
426  Fp(Fe);
427 
428  // This vector will hold the degree of freedom indices for
429  // the element. These define where in the global system
430  // the element degrees of freedom get mapped.
431  std::vector<dof_id_type> dof_indices;
432  std::vector<dof_id_type> dof_indices_u;
433  std::vector<dof_id_type> dof_indices_v;
434  std::vector<dof_id_type> dof_indices_p;
435  std::vector<dof_id_type> dof_indices_alpha;
436 
437  // Find out what the timestep size parameter is from the system, and
438  // the value of theta for the theta method. We use implicit Euler (theta=1)
439  // for this simulation even though it is only first-order accurate in time.
440  // The reason for this decision is that the second-order Crank-Nicolson
441  // method is notoriously oscillatory for problems with discontinuous
442  // initial data such as the lid-driven cavity. Therefore,
443  // we sacrifice accuracy in time for stability, but since the solution
444  // reaches steady state relatively quickly we can afford to take small
445  // timesteps. If you monitor the initial nonlinear residual for this
446  // simulation, you should see that it is monotonically decreasing in time.
447  const Real dt = es.parameters.get<Real>("dt");
448  const Real theta = 1.;
449 
450  // The kinematic viscosity, multiplies the "viscous" terms.
451  const Real nu = es.parameters.get<Real>("nu");
452 
453  // Now we will loop over all the elements in the mesh that
454  // live on the local processor. We will compute the element
455  // matrix and right-hand-side contribution. Since the mesh
456  // will be refined we want to only consider the ACTIVE elements,
457  // hence we use a variant of the active_elem_iterator.
458  for (const auto & elem : mesh.active_local_element_ptr_range())
459  {
460  // Get the degree of freedom indices for the
461  // current element. These define where in the global
462  // matrix and right-hand-side this element will
463  // contribute to.
464  dof_map.dof_indices (elem, dof_indices);
465  dof_map.dof_indices (elem, dof_indices_u, u_var);
466  dof_map.dof_indices (elem, dof_indices_v, v_var);
467  dof_map.dof_indices (elem, dof_indices_p, p_var);
468  dof_map.dof_indices (elem, dof_indices_alpha, alpha_var);
469 
470  const unsigned int n_dofs = dof_indices.size();
471  const unsigned int n_u_dofs = dof_indices_u.size();
472  const unsigned int n_v_dofs = dof_indices_v.size();
473  const unsigned int n_p_dofs = dof_indices_p.size();
474 
475  // Compute the element-specific data for the current
476  // element. This involves computing the location of the
477  // quadrature points (q_point) and the shape functions
478  // (phi, dphi) for the current element.
479  fe_vel->reinit (elem);
480  fe_pres->reinit (elem);
481 
482  // Zero the element matrix and right-hand side before
483  // summing them. We use the resize member here because
484  // the number of degrees of freedom might have changed from
485  // the last element. Note that this will be the case if the
486  // element type is different (i.e. the last element was a
487  // triangle, now we are on a quadrilateral).
488  Ke.resize (n_dofs, n_dofs);
489  Fe.resize (n_dofs);
490 
491  // Reposition the submatrices... The idea is this:
492  //
493  // - - - -
494  // | Kuu Kuv Kup | | Fu |
495  // Ke = | Kvu Kvv Kvp |; Fe = | Fv |
496  // | Kpu Kpv Kpp | | Fp |
497  // - - - -
498  //
499  // The DenseSubMatrix.reposition () member takes the
500  // (row_offset, column_offset, row_size, column_size).
501  //
502  // Similarly, the DenseSubVector.reposition () member
503  // takes the (row_offset, row_size)
504  Kuu.reposition (u_var*n_u_dofs, u_var*n_u_dofs, n_u_dofs, n_u_dofs);
505  Kuv.reposition (u_var*n_u_dofs, v_var*n_u_dofs, n_u_dofs, n_v_dofs);
506  Kup.reposition (u_var*n_u_dofs, p_var*n_u_dofs, n_u_dofs, n_p_dofs);
507 
508  Kvu.reposition (v_var*n_v_dofs, u_var*n_v_dofs, n_v_dofs, n_u_dofs);
509  Kvv.reposition (v_var*n_v_dofs, v_var*n_v_dofs, n_v_dofs, n_v_dofs);
510  Kvp.reposition (v_var*n_v_dofs, p_var*n_v_dofs, n_v_dofs, n_p_dofs);
511 
512  Kpu.reposition (p_var*n_u_dofs, u_var*n_u_dofs, n_p_dofs, n_u_dofs);
513  Kpv.reposition (p_var*n_u_dofs, v_var*n_u_dofs, n_p_dofs, n_v_dofs);
514  Kpp.reposition (p_var*n_u_dofs, p_var*n_u_dofs, n_p_dofs, n_p_dofs);
515 
516  // Also, add a row and a column to constrain the pressure
517  Kp_alpha.reposition (p_var*n_u_dofs, p_var*n_u_dofs+n_p_dofs, n_p_dofs, 1);
518  Kalpha_p.reposition (p_var*n_u_dofs+n_p_dofs, p_var*n_u_dofs, 1, n_p_dofs);
519 
520 
521  Fu.reposition (u_var*n_u_dofs, n_u_dofs);
522  Fv.reposition (v_var*n_u_dofs, n_v_dofs);
523  Fp.reposition (p_var*n_u_dofs, n_p_dofs);
524 
525  // Now we will build the element matrix and right-hand-side.
526  // Constructing the RHS requires the solution and its
527  // gradient from the previous timestep. This must be
528  // calculated at each quadrature point by summing the
529  // solution degree-of-freedom values by the appropriate
530  // weight functions.
531  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
532  {
533  // Values to hold the solution & its gradient at the previous timestep.
534  Number u = 0., u_old = 0.;
535  Number v = 0., v_old = 0.;
536  Number p_old = 0.;
537  Gradient grad_u, grad_u_old;
538  Gradient grad_v, grad_v_old;
539 
540  // Compute the velocity & its gradient from the previous timestep
541  // and the old Newton iterate.
542  for (unsigned int l=0; l<n_u_dofs; l++)
543  {
544  // From the old timestep:
545  u_old += phi[l][qp]*navier_stokes_system.old_solution (dof_indices_u[l]);
546  v_old += phi[l][qp]*navier_stokes_system.old_solution (dof_indices_v[l]);
547  grad_u_old.add_scaled (dphi[l][qp], navier_stokes_system.old_solution (dof_indices_u[l]));
548  grad_v_old.add_scaled (dphi[l][qp], navier_stokes_system.old_solution (dof_indices_v[l]));
549 
550  // From the previous Newton iterate:
551  u += phi[l][qp]*navier_stokes_system.current_solution (dof_indices_u[l]);
552  v += phi[l][qp]*navier_stokes_system.current_solution (dof_indices_v[l]);
553  grad_u.add_scaled (dphi[l][qp], navier_stokes_system.current_solution (dof_indices_u[l]));
554  grad_v.add_scaled (dphi[l][qp], navier_stokes_system.current_solution (dof_indices_v[l]));
555  }
556 
557  // Compute the old pressure value at this quadrature point.
558  for (unsigned int l=0; l<n_p_dofs; l++)
559  p_old += psi[l][qp]*navier_stokes_system.old_solution (dof_indices_p[l]);
560 
561  // Definitions for convenience. It is sometimes simpler to do a
562  // dot product if you have the full vector at your disposal.
563  const NumberVectorValue U_old (u_old, v_old);
564  const NumberVectorValue U (u, v);
565  const Number u_x = grad_u(0);
566  const Number u_y = grad_u(1);
567  const Number v_x = grad_v(0);
568  const Number v_y = grad_v(1);
569 
570  // First, an i-loop over the velocity degrees of freedom.
571  // We know that n_u_dofs == n_v_dofs so we can compute contributions
572  // for both at the same time.
573  for (unsigned int i=0; i<n_u_dofs; i++)
574  {
575  Fu(i) += JxW[qp]*(u_old*phi[i][qp] - // mass-matrix term
576  (1.-theta)*dt*(U_old*grad_u_old)*phi[i][qp] + // convection term
577  (1.-theta)*dt*p_old*dphi[i][qp](0) - // pressure term on rhs
578  (1.-theta)*dt*nu*(grad_u_old*dphi[i][qp]) + // diffusion term on rhs
579  theta*dt*(U*grad_u)*phi[i][qp]); // Newton term
580 
581 
582  Fv(i) += JxW[qp]*(v_old*phi[i][qp] - // mass-matrix term
583  (1.-theta)*dt*(U_old*grad_v_old)*phi[i][qp] + // convection term
584  (1.-theta)*dt*p_old*dphi[i][qp](1) - // pressure term on rhs
585  (1.-theta)*dt*nu*(grad_v_old*dphi[i][qp]) + // diffusion term on rhs
586  theta*dt*(U*grad_v)*phi[i][qp]); // Newton term
587 
588 
589  // Note that the Fp block is identically zero unless we are using
590  // some kind of artificial compressibility scheme...
591 
592  // Matrix contributions for the uu and vv couplings.
593  for (unsigned int j=0; j<n_u_dofs; j++)
594  {
595  Kuu(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*u_x*phi[i][qp]*phi[j][qp]); // Newton term
599 
600  Kuv(i,j) += JxW[qp]*theta*dt*u_y*phi[i][qp]*phi[j][qp]; // Newton term
601 
602  Kvv(i,j) += JxW[qp]*(phi[i][qp]*phi[j][qp] + // mass matrix term
603  theta*dt*nu*(dphi[i][qp]*dphi[j][qp]) + // diffusion term
604  theta*dt*(U*dphi[j][qp])*phi[i][qp] + // convection term
605  theta*dt*v_y*phi[i][qp]*phi[j][qp]); // Newton term
606 
607  Kvu(i,j) += JxW[qp]*theta*dt*v_x*phi[i][qp]*phi[j][qp]; // Newton term
608  }
609 
610  // Matrix contributions for the up and vp couplings.
611  for (unsigned int j=0; j<n_p_dofs; j++)
612  {
613  Kup(i,j) += JxW[qp]*(-theta*dt*psi[j][qp]*dphi[i][qp](0));
614  Kvp(i,j) += JxW[qp]*(-theta*dt*psi[j][qp]*dphi[i][qp](1));
615  }
616  }
617 
618  // Now an i-loop over the pressure degrees of freedom. This code computes
619  // the matrix entries due to the continuity equation. Note: To maintain a
620  // symmetric matrix, we may (or may not) multiply the continuity equation by
621  // negative one. Here we do not.
622  for (unsigned int i=0; i<n_p_dofs; i++)
623  {
624  Kp_alpha(i,0) += JxW[qp]*psi[i][qp];
625  Kalpha_p(0,i) += JxW[qp]*psi[i][qp];
626  for (unsigned int j=0; j<n_u_dofs; j++)
627  {
628  Kpu(i,j) += JxW[qp]*psi[i][qp]*dphi[j][qp](0);
629  Kpv(i,j) += JxW[qp]*psi[i][qp]*dphi[j][qp](1);
630  }
631  }
632  } // end of the quadrature point qp-loop
633 
634  // Since we're using heterogeneous DirichletBoundary objects for
635  // the boundary conditions, we need to call a specific function
636  // to constrain the element stiffness matrix.
637  dof_map.heterogenously_constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
638 
639  // The element matrix and right-hand-side are now built
640  // for this element. Add them to the global matrix and
641  // right-hand-side vector. The SparseMatrix::add_matrix()
642  // and NumericVector::add_vector() members do this for us.
643  navier_stokes_system.matrix->add_matrix (Ke, dof_indices);
644  navier_stokes_system.rhs->add_vector (Fe, dof_indices);
645  } // end of element loop
646 
647  // We can set the mean of the pressure by setting Falpha. Typically
648  // a value of zero is chosen, but the value should be arbitrary.
649  navier_stokes_system.rhs->add(navier_stokes_system.rhs->size()-1, 10.);
650 #else
651  libmesh_ignore(es);
652 #endif
653 }

References libMesh::MeshBase::active_local_element_ptr_range(), libMesh::TypeVector< T >::add_scaled(), libMesh::FEGenericBase< OutputType >::build(), libMesh::FEType::default_quadrature_order(), dim, libMesh::DofMap::dof_indices(), libMesh::Parameters::get(), libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), libMesh::DofMap::heterogenously_constrain_element_matrix_and_vector(), libMesh::libmesh_ignore(), mesh, libMesh::MeshBase::mesh_dimension(), libMesh::TransientSystem< Base >::old_solution(), libMesh::EquationSystems::parameters, libMesh::Real, libMesh::DenseSubVector< T >::reposition(), libMesh::DenseSubMatrix< T >::reposition(), libMesh::DenseVector< T >::resize(), and libMesh::DenseMatrix< T >::resize().

◆ assemble_stokes() [2/2]

void assemble_stokes ( EquationSystems es,
const std::string &  system_name 
)

Referenced by main().

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 80 of file systems_of_equations_ex3.C.

81 {
82  // Initialize libMesh.
83  LibMeshInit init (argc, argv);
84 
85  // This example requires a linear solver package.
86  libmesh_example_requires(libMesh::default_solver_package() != INVALID_SOLVER_PACKAGE,
87  "--enable-petsc, --enable-trilinos, or --enable-eigen");
88 
89  // Skip this 2D example if libMesh was compiled as 1D-only.
90  libmesh_example_requires(2 <= LIBMESH_DIM, "2D support");
91 
92  // We use Dirichlet boundary conditions here
93 #ifndef LIBMESH_ENABLE_DIRICHLET
94  libmesh_example_requires(false, "--enable-dirichlet");
95 #endif
96 
97  // This example NaNs with the Eigen sparse linear solvers and
98  // Trilinos solvers, but should work OK with either PETSc or
99  // Laspack.
100  libmesh_example_requires(libMesh::default_solver_package() != EIGEN_SOLVERS, "--enable-petsc or --enable-laspack");
101  libmesh_example_requires(libMesh::default_solver_package() != TRILINOS_SOLVERS, "--enable-petsc or --enable-laspack");
102 
103  // Create a mesh, with dimension to be overridden later, distributed
104  // across the default MPI communicator.
105  Mesh mesh(init.comm());
106 
107  // Use the MeshTools::Generation mesh generator to create a uniform
108  // 2D grid on the square [-1,1]^2. We instruct the mesh generator
109  // to build a mesh of 8x8 Quad9 elements in 2D. Building these
110  // higher-order elements allows us to use higher-order
111  // approximation, as in example 3.
113  20, 20,
114  0., 1.,
115  0., 1.,
116  QUAD9);
117 
118  // Print information about the mesh to the screen.
119  mesh.print_info();
120 
121  // Create an equation systems object.
122  EquationSystems equation_systems (mesh);
123 
124  // Declare the system and its variables.
125  // Creates a transient system named "Navier-Stokes"
127  equation_systems.add_system<TransientLinearImplicitSystem> ("Navier-Stokes");
128 
129  // Add the variables "vel_x" & "vel_y" to "Navier-Stokes". They
130  // will be approximated using second-order approximation.
131  system.add_variable ("vel_x", SECOND);
132  system.add_variable ("vel_y", SECOND);
133 
134  // Add the variable "p" to "Navier-Stokes". This will
135  // be approximated with a first-order basis,
136  // providing an LBB-stable pressure-velocity pair.
137  system.add_variable ("p", FIRST);
138 
139  // Add a scalar Lagrange multiplier to constrain the
140  // pressure to have zero mean.
141  system.add_variable ("alpha", FIRST, SCALAR);
142 
143  // Give the system a pointer to the matrix assembly
144  // function.
145  system.attach_assemble_function (assemble_stokes);
146 
147  // Set Dirichlet boundary conditions.
148  set_lid_driven_bcs(system);
149 
150  // Initialize the data structures for the equation system.
151  equation_systems.init ();
152 
153  // Prints information about the system to the screen.
154  equation_systems.print_info();
155 
156  // Create a performance-logging object for this example
157  PerfLog perf_log("Systems Example 3");
158 
159  // Get a reference to the Stokes system to use later.
160  TransientLinearImplicitSystem & navier_stokes_system =
161  equation_systems.get_system<TransientLinearImplicitSystem>("Navier-Stokes");
162 
163  // Now we begin the timestep loop to compute the time-accurate
164  // solution of the equations.
165  const Real dt = 0.1;
166  navier_stokes_system.time = 0.0;
167  const unsigned int n_timesteps = 15;
168 
169  // The number of steps and the stopping criterion are also required
170  // for the nonlinear iterations.
171  const unsigned int n_nonlinear_steps = 15;
172  const Real nonlinear_tolerance = TOLERANCE*10;
173 
174  // We also set a standard linear solver flag in the EquationSystems object
175  // which controls the maximum number of linear solver iterations allowed.
176  equation_systems.parameters.set<unsigned int>("linear solver maximum iterations") = 250;
177 
178  // Tell the system of equations what the timestep is by using
179  // the set_parameter function. The matrix assembly routine can
180  // then reference this parameter.
181  equation_systems.parameters.set<Real> ("dt") = dt;
182 
183  // The kinematic viscosity, nu = mu/rho, units of length**2/time.
184  equation_systems.parameters.set<Real> ("nu") = .007;
185 
186  // The first thing to do is to get a copy of the solution at
187  // the current nonlinear iteration. This value will be used to
188  // determine if we can exit the nonlinear loop.
189  std::unique_ptr<NumericVector<Number>>
190  last_nonlinear_soln (navier_stokes_system.solution->clone());
191 
192 #ifdef LIBMESH_HAVE_EXODUS_API
193  // Since we are not doing adaptivity, write all solutions to a single Exodus file.
194  ExodusII_IO exo_io(mesh);
195 
196  // Write out the initial condition
197  exo_io.write_equation_systems ("out.e", equation_systems);
198 #endif
199 
200  for (unsigned int t_step=1; t_step<=n_timesteps; ++t_step)
201  {
202  // Increment the time counter, set the time and the
203  // time step size as parameters in the EquationSystem.
204  navier_stokes_system.time += dt;
205 
206  // A pretty update message
207  libMesh::out << "\n\n*** Solving time step "
208  << t_step
209  << ", time = "
210  << navier_stokes_system.time
211  << " ***"
212  << std::endl;
213 
214  // Now we need to update the solution vector from the
215  // previous time step. This is done directly through
216  // the reference to the Stokes system.
217  *navier_stokes_system.old_local_solution = *navier_stokes_system.current_local_solution;
218 
219  // At the beginning of each solve, reset the linear solver tolerance
220  // to a "reasonable" starting value.
221  const Real initial_linear_solver_tol = 1.e-6;
222  equation_systems.parameters.set<Real> ("linear solver tolerance") = initial_linear_solver_tol;
223 
224  // We'll set this flag when convergence is (hopefully) achieved.
225  bool converged = false;
226 
227  // Now we begin the nonlinear loop
228  for (unsigned int l=0; l<n_nonlinear_steps; ++l)
229  {
230  // Update the nonlinear solution.
231  last_nonlinear_soln->zero();
232  last_nonlinear_soln->add(*navier_stokes_system.solution);
233 
234  // Assemble & solve the linear system.
235  perf_log.push("linear solve");
236  equation_systems.get_system("Navier-Stokes").solve();
237  perf_log.pop("linear solve");
238 
239  // Compute the difference between this solution and the last
240  // nonlinear iterate.
241  last_nonlinear_soln->add (-1., *navier_stokes_system.solution);
242 
243  // Close the vector before computing its norm
244  last_nonlinear_soln->close();
245 
246  // Compute the l2 norm of the difference
247  const Real norm_delta = last_nonlinear_soln->l2_norm();
248 
249  // How many iterations were required to solve the linear system?
250  const unsigned int n_linear_iterations = navier_stokes_system.n_linear_iterations();
251 
252  // What was the final residual of the linear system?
253  const Real final_linear_residual = navier_stokes_system.final_linear_residual();
254 
255  // If the solver did no work (sometimes -ksp_converged_reason
256  // says "Linear solve converged due to CONVERGED_RTOL
257  // iterations 0") but the nonlinear residual norm is above
258  // the tolerance, we need to pick an even lower linear
259  // solver tolerance and try again. Note that the tolerance
260  // is relative to the norm of the RHS, which for this
261  // particular problem does not go to zero, since we are
262  // solving for the full solution rather than the update.
263  //
264  // Similarly, if the solver did no work and this is the 0th
265  // nonlinear step, it means that the delta between solutions
266  // is being inaccurately measured as "0" since the solution
267  // did not change. Decrease the tolerance and try again.
268  if (n_linear_iterations == 0 &&
269  (navier_stokes_system.final_linear_residual() >= nonlinear_tolerance || l==0))
270  {
271  Real old_linear_solver_tolerance = equation_systems.parameters.get<Real> ("linear solver tolerance");
272  equation_systems.parameters.set<Real> ("linear solver tolerance") = 1.e-3 * old_linear_solver_tolerance;
273  continue;
274  }
275 
276  // Print out convergence information for the linear and
277  // nonlinear iterations.
278  libMesh::out << "Linear solver converged at step: "
279  << n_linear_iterations
280  << ", final residual: "
281  << final_linear_residual
282  << " Nonlinear convergence: ||u - u_old|| = "
283  << norm_delta
284  << std::endl;
285 
286  // Terminate the solution iteration if the difference between
287  // this nonlinear iterate and the last is sufficiently small, AND
288  // if the most recent linear system was solved to a sufficient tolerance.
289  if ((norm_delta < nonlinear_tolerance) &&
290  (navier_stokes_system.final_linear_residual() < nonlinear_tolerance))
291  {
292  libMesh::out << " Nonlinear solver converged at step "
293  << l
294  << std::endl;
295  converged = true;
296  break;
297  }
298 
299  // Otherwise, decrease the linear system tolerance. For the inexact Newton
300  // method, the linear solver tolerance needs to decrease as we get closer to
301  // the solution to ensure quadratic convergence. The new linear solver tolerance
302  // is chosen (heuristically) as the square of the previous linear system residual norm.
303  //Real flr2 = final_linear_residual*final_linear_residual;
304  Real new_linear_solver_tolerance = std::min(Utility::pow<2>(final_linear_residual), initial_linear_solver_tol);
305  equation_systems.parameters.set<Real> ("linear solver tolerance") = new_linear_solver_tolerance;
306  } // end nonlinear loop
307 
308  // Don't keep going if we failed to converge.
309  if (!converged)
310  libmesh_error_msg("Error: Newton iterations failed to converge!");
311 
312 #ifdef LIBMESH_HAVE_EXODUS_API
313  // Write out every nth timestep to file.
314  const unsigned int write_interval = 1;
315 
316  if ((t_step+1)%write_interval == 0)
317  {
318  exo_io.write_timestep("out.e",
319  equation_systems,
320  t_step+1, // we're off by one since we wrote the IC and the Exodus numbering is 1-based.
321  navier_stokes_system.time);
322  }
323 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
324  } // end timestep loop.
325 
326  // All done.
327  return 0;
328 }

References libMesh::EquationSystems::add_system(), assemble_stokes(), libMesh::MeshTools::Generation::build_square(), libMesh::default_solver_package(), libMesh::EIGEN_SOLVERS, libMesh::FIRST, libMesh::Parameters::get(), libMesh::EquationSystems::get_system(), libMesh::TriangleWrapper::init(), libMesh::EquationSystems::init(), libMesh::INVALID_SOLVER_PACKAGE, mesh, libMesh::TransientSystem< Base >::old_local_solution, libMesh::out, libMesh::EquationSystems::parameters, libMesh::PerfLog::pop(), libMesh::EquationSystems::print_info(), libMesh::MeshBase::print_info(), libMesh::PerfLog::push(), libMesh::QUAD9, libMesh::Real, libMesh::SCALAR, libMesh::SECOND, libMesh::Parameters::set(), set_lid_driven_bcs(), libMesh::TOLERANCE, and libMesh::TRILINOS_SOLVERS.

◆ set_lid_driven_bcs()

void set_lid_driven_bcs ( TransientLinearImplicitSystem system)

Definition at line 657 of file systems_of_equations_ex3.C.

658 {
659 #ifdef LIBMESH_ENABLE_DIRICHLET
660  unsigned short int
661  u_var = system.variable_number("vel_x"),
662  v_var = system.variable_number("vel_y");
663 
664  // Get a convenient reference to the System's DofMap
665  DofMap & dof_map = system.get_dof_map();
666 
667  {
668  // u=v=0 on bottom, left, right
669  std::set<boundary_id_type> boundary_ids;
670  boundary_ids.insert(0);
671  boundary_ids.insert(1);
672  boundary_ids.insert(3);
673 
674  std::vector<unsigned int> variables;
675  variables.push_back(u_var);
676  variables.push_back(v_var);
677 
678  dof_map.add_dirichlet_boundary(DirichletBoundary(boundary_ids,
679  variables,
681  }
682  {
683  // u=1 on top
684  std::set<boundary_id_type> boundary_ids;
685  boundary_ids.insert(2);
686 
687  std::vector<unsigned int> variables;
688  variables.push_back(u_var);
689 
690  dof_map.add_dirichlet_boundary(DirichletBoundary(boundary_ids,
691  variables,
692  ConstFunction<Number>(1.)));
693  }
694  {
695  // v=0 on top
696  std::set<boundary_id_type> boundary_ids;
697  boundary_ids.insert(2);
698 
699  std::vector<unsigned int> variables;
700  variables.push_back(v_var);
701 
702  dof_map.add_dirichlet_boundary(DirichletBoundary(boundary_ids,
703  variables,
705  }
706 #else
707  libmesh_ignore(system);
708 #endif
709 }

References libMesh::DofMap::add_dirichlet_boundary(), and libMesh::libmesh_ignore().

Referenced by main().

libMesh::Number
Real Number
Definition: libmesh_common.h:195
libMesh::Mesh
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50
libMesh::EquationSystems::get_mesh
const MeshBase & get_mesh() const
Definition: equation_systems.h:637
libMesh::TransientSystem
Manages storage and variables for transient systems.
Definition: transient_system.h:57
libMesh::QGauss
This class implements specific orders of Gauss quadrature.
Definition: quadrature_gauss.h:39
libMesh::DofMap::dof_indices
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:1967
libMesh::DofMap::add_dirichlet_boundary
void add_dirichlet_boundary(const DirichletBoundary &dirichlet_boundary)
Adds a copy of the specified Dirichlet boundary to the system.
Definition: dof_map_constraints.C:4390
assemble_stokes
void assemble_stokes(EquationSystems &es, const std::string &system_name)
libMesh::MeshBase::active_local_element_ptr_range
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
libMesh::DenseSubMatrix
Defines a dense submatrix for use in Finite Element-type computations.
Definition: dense_submatrix.h:45
libMesh::EquationSystems::get_system
const T_sys & get_system(const std::string &name) const
Definition: equation_systems.h:757
libMesh::TOLERANCE
static const Real TOLERANCE
Definition: libmesh_common.h:128
libMesh::DenseMatrix< Number >
libMesh::default_solver_package
SolverPackage default_solver_package()
Definition: libmesh.C:993
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::MeshBase::mesh_dimension
unsigned int mesh_dimension() const
Definition: mesh_base.C:135
libMesh::ExodusII_IO
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Definition: exodusII_io.h:51
libMesh::SECOND
Definition: enum_order.h:43
dim
unsigned int dim
Definition: adaptivity_ex3.C:113
libMesh::DenseMatrix::resize
void resize(const unsigned int new_m, const unsigned int new_n)
Resize the matrix.
Definition: dense_matrix.h:822
libMesh::VectorValue< Number >
libMesh::TriangleWrapper::init
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
libMesh::MeshTools::Generation::build_square
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.
Definition: mesh_generation.C:1501
libMesh::MeshBase
This is the MeshBase class.
Definition: mesh_base.h:78
libMesh::INVALID_SOLVER_PACKAGE
Definition: enum_solver_package.h:43
libMesh::libmesh_ignore
void libmesh_ignore(const Args &...)
Definition: libmesh_common.h:526
libMesh::PerfLog
The PerfLog class allows monitoring of specific events.
Definition: perf_log.h:125
libMesh::TypeVector::add_scaled
void add_scaled(const TypeVector< T2 > &, const T &)
Add a scaled value to this vector without creating a temporary.
Definition: type_vector.h:665
libMesh::LibMeshInit
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:83
libMesh::ZeroFunction
ConstFunction that simply returns 0.
Definition: zero_function.h:36
libMesh::EquationSystems
This is the EquationSystems class.
Definition: equation_systems.h:74
libMesh::DenseSubVector< Number >
libMesh::TRILINOS_SOLVERS
Definition: enum_solver_package.h:37
libMesh::DenseVector::resize
void resize(const unsigned int n)
Resize the vector.
Definition: dense_vector.h:355
libMesh::FEType
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
libMesh::DofMap::heterogenously_constrain_element_matrix_and_vector
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.
Definition: dof_map.h:2040
set_lid_driven_bcs
void set_lid_driven_bcs(TransientLinearImplicitSystem &system)
Definition: systems_of_equations_ex3.C:657
libMesh::DofMap
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:176
libMesh::MeshBase::print_info
void print_info(std::ostream &os=libMesh::out) const
Prints relevant information about the mesh.
Definition: mesh_base.C:585
libMesh::DirichletBoundary
This class allows one to associate Dirichlet boundary values with a given set of mesh boundary ids an...
Definition: dirichlet_boundaries.h:88
libMesh::TransientSystem::old_solution
Number old_solution(const dof_id_type global_dof_number) const
Definition: transient_system.C:122
libMesh::QUAD9
Definition: enum_elem_type.h:43
libMesh::SCALAR
Definition: enum_fe_family.h:58
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::FEType::default_quadrature_order
Order default_quadrature_order() const
Definition: fe_type.h:353
libMesh::TransientSystem::old_local_solution
std::unique_ptr< NumericVector< Number > > old_local_solution
All the values I need to compute my contribution to the simulation at hand.
Definition: transient_system.h:119
libMesh::Parameters::get
const T & get(const std::string &) const
Definition: parameters.h:421
libMesh::out
OStreamProxy out
libMesh::FIRST
Definition: enum_order.h:42
libMesh::ConstFunction
Function that returns a single value that never changes.
Definition: const_function.h:43
libMesh::EIGEN_SOLVERS
Definition: enum_solver_package.h:40
libMesh::DenseVector< Number >
libMesh::EquationSystems::parameters
Parameters parameters
Data structure holding arbitrary parameters.
Definition: equation_systems.h:557