libMesh
Functions
systems_of_equations_ex2.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)
 
void set_stagnation_bcs (TransientLinearImplicitSystem &system)
 
void set_poiseuille_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 342 of file systems_of_equations_ex2.C.

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

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 92 of file systems_of_equations_ex2.C.

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

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::SECOND, libMesh::Parameters::set(), set_lid_driven_bcs(), and libMesh::TRILINOS_SOLVERS.

◆ set_lid_driven_bcs()

void set_lid_driven_bcs ( TransientLinearImplicitSystem system)

Definition at line 672 of file systems_of_equations_ex2.C.

673 {
674  // This problem *does* require a pressure pin, there are Dirichlet
675  // boundary conditions for u and v on the entire boundary.
676  system.get_equation_systems().parameters.set<bool>("pin_pressure") = true;
677 
678 #ifdef LIBMESH_ENABLE_DIRICHLET
679  unsigned short int
680  u_var = system.variable_number("vel_x"),
681  v_var = system.variable_number("vel_y");
682 
683  // Get a convenient reference to the System's DofMap
684  DofMap & dof_map = system.get_dof_map();
685 
686  {
687  // u=v=0 on bottom, left, right
688  std::set<boundary_id_type> boundary_ids;
689  boundary_ids.insert(0);
690  boundary_ids.insert(1);
691  boundary_ids.insert(3);
692 
693  std::vector<unsigned int> variables;
694  variables.push_back(u_var);
695  variables.push_back(v_var);
696 
697  dof_map.add_dirichlet_boundary(DirichletBoundary(boundary_ids,
698  variables,
700  }
701  {
702  // u=1 on top
703  std::set<boundary_id_type> boundary_ids;
704  boundary_ids.insert(2);
705 
706  std::vector<unsigned int> variables;
707  variables.push_back(u_var);
708 
709  dof_map.add_dirichlet_boundary(DirichletBoundary(boundary_ids,
710  variables,
711  ConstFunction<Number>(1.)));
712  }
713  {
714  // v=0 on top
715  std::set<boundary_id_type> boundary_ids;
716  boundary_ids.insert(2);
717 
718  std::vector<unsigned int> variables;
719  variables.push_back(v_var);
720 
721  dof_map.add_dirichlet_boundary(DirichletBoundary(boundary_ids,
722  variables,
724  }
725 #endif // LIBMESH_ENABLE_DIRICHLET
726 }

References libMesh::DofMap::add_dirichlet_boundary().

Referenced by main().

◆ set_poiseuille_bcs()

void set_poiseuille_bcs ( TransientLinearImplicitSystem system)

Definition at line 818 of file systems_of_equations_ex2.C.

819 {
820  // This problem does not require a pressure pin, the Neumann outlet
821  // BCs are sufficient to set the value of the pressure.
822  system.get_equation_systems().parameters.set<bool>("pin_pressure") = false;
823 
824 #ifdef LIBMESH_ENABLE_DIRICHLET
825  unsigned short int
826  u_var = system.variable_number("vel_x"),
827  v_var = system.variable_number("vel_y");
828 
829  // Get a convenient reference to the System's DofMap
830  DofMap & dof_map = system.get_dof_map();
831 
832  {
833  // u=v=0 on top, bottom
834  std::set<boundary_id_type> boundary_ids;
835  boundary_ids.insert(0);
836  boundary_ids.insert(2);
837 
838  std::vector<unsigned int> variables;
839  variables.push_back(u_var);
840  variables.push_back(v_var);
841 
842  dof_map.add_dirichlet_boundary(DirichletBoundary(boundary_ids,
843  variables,
845  }
846  {
847  // u=quadratic on left
848  std::set<boundary_id_type> boundary_ids;
849  boundary_ids.insert(3);
850 
851  std::vector<unsigned int> variables;
852  variables.push_back(u_var);
853 
854  dof_map.add_dirichlet_boundary(DirichletBoundary(boundary_ids,
855  variables,
856  ParsedFunction<Number>("4*y*(1-y)")));
857  }
858  {
859  // v=0 on left
860  std::set<boundary_id_type> boundary_ids;
861  boundary_ids.insert(3);
862 
863  std::vector<unsigned int> variables;
864  variables.push_back(v_var);
865 
866  dof_map.add_dirichlet_boundary(DirichletBoundary(boundary_ids,
867  variables,
869  }
870 #endif // LIBMESH_ENABLE_DIRICHLET
871 }

References libMesh::DofMap::add_dirichlet_boundary().

◆ set_stagnation_bcs()

void set_stagnation_bcs ( TransientLinearImplicitSystem system)

Definition at line 730 of file systems_of_equations_ex2.C.

731 {
732  // This problem does not require a pressure pin, the Neumann outlet
733  // BCs are sufficient to set the value of the pressure.
734  system.get_equation_systems().parameters.set<bool>("pin_pressure") = false;
735 
736 #ifdef LIBMESH_ENABLE_DIRICHLET
737  unsigned short int
738  u_var = system.variable_number("vel_x"),
739  v_var = system.variable_number("vel_y");
740 
741  // Get a convenient reference to the System's DofMap
742  DofMap & dof_map = system.get_dof_map();
743 
744  {
745  // u=v=0 on bottom
746  std::set<boundary_id_type> boundary_ids;
747  boundary_ids.insert(0);
748 
749  std::vector<unsigned int> variables;
750  variables.push_back(u_var);
751  variables.push_back(v_var);
752 
753  dof_map.add_dirichlet_boundary(DirichletBoundary(boundary_ids,
754  variables,
756  }
757  {
758  // u=0 on left (symmetry)
759  std::set<boundary_id_type> boundary_ids;
760  boundary_ids.insert(3);
761 
762  std::vector<unsigned int> variables;
763  variables.push_back(u_var);
764 
765  dof_map.add_dirichlet_boundary(DirichletBoundary(boundary_ids,
766  variables,
768  }
769  {
770  // u = k*x on top
771  std::set<boundary_id_type> boundary_ids;
772  boundary_ids.insert(2);
773 
774  std::vector<unsigned int> variables;
775  variables.push_back(u_var);
776 
777  // Set up ParsedFunction parameters
778  std::vector<std::string> additional_vars;
779  additional_vars.push_back("k");
780  std::vector<Number> initial_vals;
781  initial_vals.push_back(1.);
782 
783  dof_map.add_dirichlet_boundary(DirichletBoundary(boundary_ids,
784  variables,
786  &additional_vars,
787  &initial_vals)));
788  }
789  {
790  // v = -k*y on top
791  std::set<boundary_id_type> boundary_ids;
792  boundary_ids.insert(2);
793 
794  std::vector<unsigned int> variables;
795  variables.push_back(v_var);
796 
797  // Set up ParsedFunction parameters
798  std::vector<std::string> additional_vars;
799  additional_vars.push_back("k");
800  std::vector<Number> initial_vals;
801  initial_vals.push_back(1.);
802 
803  // Note: we have to specify LOCAL_VARIABLE_ORDER here, since we're
804  // using a ParsedFunction to set the value of v_var, which is
805  // actually the second variable in the system.
806  dof_map.add_dirichlet_boundary(DirichletBoundary(boundary_ids,
807  variables,
808  ParsedFunction<Number>("-k*y",
809  &additional_vars,
810  &initial_vals),
812  }
813 #endif // LIBMESH_ENABLE_DIRICHLET
814 }

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

libMesh::Number
Real Number
Definition: libmesh_common.h:195
libMesh::ParsedFunction
A Function generated (via FParser) by parsing a mathematical expression.
Definition: parsed_function.h:60
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
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::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::LOCAL_VARIABLE_ORDER
Definition: dirichlet_boundaries.h:62
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
set_lid_driven_bcs
void set_lid_driven_bcs(TransientLinearImplicitSystem &system)
Definition: systems_of_equations_ex2.C:672
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
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::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 >
assemble_stokes
void assemble_stokes(EquationSystems &es, const std::string &system_name)
libMesh::EquationSystems::parameters
Parameters parameters
Data structure holding arbitrary parameters.
Definition: equation_systems.h:557