libMesh
miscellaneous_ex17.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 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 // <h1>Miscennaleous Example 17 - Demonstrating mix of preallocated/hash-table matrix
19 // assemblies</h1> \author Alexander D. Lindsay \date 2024
20 
21 // Basic include files needed for the mesh functionality.
22 #include "libmesh/libmesh.h"
23 #include "libmesh/mesh.h"
24 #include "libmesh/mesh_generation.h"
25 #include "libmesh/linear_implicit_system.h"
26 #include "libmesh/equation_systems.h"
27 
28 // Define the Finite Element object.
29 #include "libmesh/fe.h"
30 
31 // Define Gauss quadrature rules.
32 #include "libmesh/quadrature_gauss.h"
33 
34 // Define useful datatypes for finite element
35 // matrix and vector components.
36 #include "libmesh/numeric_vector.h"
37 #include "libmesh/sparse_matrix.h"
38 #include "libmesh/dense_matrix.h"
39 #include "libmesh/dense_vector.h"
40 #include "libmesh/elem.h"
41 #include "libmesh/enum_solver_package.h"
42 
43 // Define the DofMap, which handles degree of freedom
44 // indexing.
45 #include "libmesh/dof_map.h"
46 
47 #ifdef LIBMESH_HAVE_PETSC
48 // include PETSc headers
49 #include "libmesh/petsc_matrix.h"
50 #include "libmesh/petsc_vector.h"
51 #include "petscksp.h"
52 #endif
53 
54 #include <iostream>
55 
56 // Bring in everything from the libMesh namespace
57 using namespace libMesh;
58 
60  const std::string & system_name);
61 
62 // Function prototype for the exact solution.
63 Real exact_solution (const Real x,
64  const Real y,
65  const Real z = 0.);
66 
67 int main (int argc, char ** argv)
68 {
69  // Initialize libraries, like in example 2.
70  LibMeshInit init (argc, argv);
71 
72  // Brief message to the user regarding the program name
73  // and command line arguments.
74  libMesh::out << "Running " << argv[0];
75 
76  for (int i=1; i<argc; i++)
77  libMesh::out << " " << argv[i];
78 
79  libMesh::out << std::endl << std::endl;
80 
81  // Skip this 2D example if libMesh was compiled as 1D-only.
82  libmesh_example_requires(2 <= LIBMESH_DIM, "2D support");
83 
84  // This example is meant to test a PETSc-specific feature, so let's just skip it if
85  // libmesh is not built with Petsc support.
86  libmesh_example_requires(libMesh::default_solver_package() == PETSC_SOLVERS, "--enable-petsc");
87 
88  // Create a mesh, with dimension to be overridden later, distributed
89  // across the default MPI communicator.
90  Mesh mesh(init.comm());
91 
92  // Use the MeshTools::Generation mesh generator to create a uniform
93  // 2D grid on the square [-1,1]^2. We instruct the mesh generator
94  // to build a mesh of 15x15 QUAD9 elements. Building QUAD9
95  // elements instead of the default QUAD4's we used in example 2
96  // allow us to use higher-order approximation.
97  MeshTools::Generation::build_square(mesh, 2, 2, -1., 1., -1., 1., QUAD4);
98 
99  // Print information about the mesh to the screen.
100  // Note that 5x5 QUAD9 elements actually has 11x11 nodes,
101  // so this mesh is significantly larger than the one in example 2.
102  mesh.print_info();
103 
104  // Create an equation systems object.
105  EquationSystems equation_systems (mesh);
106 
107  // Declare the Poisson system and its variables.
108  // The Poisson system is another example of a steady system.
109  auto & system = equation_systems.add_system<LinearImplicitSystem>("Poisson");
110 
111  // Adds the variable "u" to "Poisson". "u"
112  // will be approximated using second-order approximation.
113  system.add_variable("u", FIRST);
114 
115  // Give the system a pointer to the matrix assembly
116  // function. This will be called when needed by the
117  // library.
118  system.attach_assemble_function(assemble_poisson);
119 
120  // Add the preconditioner matrix
121  system.add_matrix("preconditioner");
122 #ifdef LIBMESH_HAVE_PETSC
123 #if PETSC_RELEASE_GREATER_EQUALS(3, 19, 0)
124  system.get_matrix("preconditioner").use_hash_table(true);
125 #endif
126 #endif
127 
128  // Initialize the data structures for the equation system.
129  equation_systems.init();
130 
131  // Prints information about the system to the screen.
132  equation_systems.print_info();
133 
134  // assemble the operators and RHS
135  system.assemble();
136 
137 #ifdef LIBMESH_HAVE_PETSC
138  auto & sys_matrix = cast_ref<PetscMatrix<Number> &>(system.get_system_matrix());
139  auto & pre_matrix = cast_ref<PetscMatrix<Number> &>(system.get_matrix("preconditioner"));
140  LibmeshPetscCall2(system.comm(), PetscOptionsSetValue(NULL, "-ksp_monitor", NULL));
141 
142  auto solve = [&sys_matrix, &pre_matrix, &system]() {
143  // Make sure our matrices are closed
144  sys_matrix.close();
145  pre_matrix.close();
146  system.rhs->close();
147  KSP ksp;
148  LibmeshPetscCall2(system.comm(), KSPCreate(system.comm().get(), &ksp));
149  LibmeshPetscCall2(system.comm(), KSPSetOperators(ksp, sys_matrix.mat(), pre_matrix.mat()));
150  LibmeshPetscCall2(system.comm(), KSPSetFromOptions(ksp));
151  LibmeshPetscCall2(system.comm(),
152  KSPSolve(ksp,
153  cast_ptr<PetscVector<Number> *>(system.rhs)->vec(),
154  cast_ptr<PetscVector<Number> *>(system.solution.get())->vec()));
155  };
156 
157  // solve
158  solve();
159 
160  // MatResetHash added in PETSc version 3.23
161 #if !PETSC_VERSION_LESS_THAN(3, 23, 0)
162  // reset the memory
163  // sys_matrix.restore_original_nonzero_pattern(); # See https://gitlab.com/petsc/petsc/-/merge_requests/8063
164  pre_matrix.restore_original_nonzero_pattern();
165  // zero
166  sys_matrix.zero();
167  pre_matrix.zero();
168  system.rhs->zero();
169  system.solution->zero();
170  system.update();
171  // re-assemble
172  system.assemble();
173  // resolve
174  solve();
175 #endif
176 #endif
177 
178  // All done.
179  return 0;
180 }
181 
182 
183 
184 // We now define the matrix assembly function for the
185 // Poisson system. We need to first compute element
186 // matrices and right-hand sides, and then take into
187 // account the boundary conditions, which will be handled
188 // via a penalty method.
190  const std::string & libmesh_dbg_var(system_name))
191 {
192 
193  // It is a good idea to make sure we are assembling
194  // the proper system.
195  libmesh_assert_equal_to (system_name, "Poisson");
196 
197  // Get a constant reference to the mesh object.
198  const MeshBase & mesh = es.get_mesh();
199 
200  // The dimension that we are running
201  const unsigned int dim = mesh.mesh_dimension();
202 
203  // Get a reference to the LinearImplicitSystem we are solving
204  LinearImplicitSystem & system = es.get_system<LinearImplicitSystem> ("Poisson");
205 
206  // A reference to the DofMap object for this system. The DofMap
207  // object handles the index translation from node and element numbers
208  // to degree of freedom numbers. We will talk more about the DofMap
209  // in future examples.
210  const DofMap & dof_map = system.get_dof_map();
211 
212  // Get a constant reference to the Finite Element type
213  // for the first (and only) variable in the system.
214  FEType fe_type = dof_map.variable_type(0);
215 
216  // Build a Finite Element object of the specified type. Since the
217  // FEBase::build() member dynamically creates memory we will
218  // store the object as a std::unique_ptr<FEBase>. This can be thought
219  // of as a pointer that will clean up after itself. Introduction Example 4
220  // describes some advantages of std::unique_ptr's in the context of
221  // quadrature rules.
222  std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
223 
224  // A 5th order Gauss quadrature rule for numerical integration.
225  QGauss qrule (dim, FIFTH);
226 
227  // Tell the finite element object to use our quadrature rule.
228  fe->attach_quadrature_rule (&qrule);
229 
230  // Declare a special finite element object for
231  // boundary integration.
232  std::unique_ptr<FEBase> fe_face (FEBase::build(dim, fe_type));
233 
234  // Boundary integration requires one quadrature rule,
235  // with dimensionality one less than the dimensionality
236  // of the element.
237  QGauss qface(dim-1, FIFTH);
238 
239  // Tell the finite element object to use our
240  // quadrature rule.
241  fe_face->attach_quadrature_rule (&qface);
242 
243  // Here we define some references to cell-specific data that
244  // will be used to assemble the linear system.
245  //
246  // The element Jacobian * quadrature weight at each integration point.
247  const std::vector<Real> & JxW = fe->get_JxW();
248 
249  // The physical XY locations of the quadrature points on the element.
250  // These might be useful for evaluating spatially varying material
251  // properties at the quadrature points.
252  const std::vector<Point> & q_point = fe->get_xyz();
253 
254  // The element shape functions evaluated at the quadrature points.
255  const std::vector<std::vector<Real>> & phi = fe->get_phi();
256 
257  // The element shape function gradients evaluated at the quadrature
258  // points.
259  const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
260 
261  // Define data structures to contain the element matrix
262  // and right-hand-side vector contribution. Following
263  // basic finite element terminology we will denote these
264  // "Ke" and "Fe". These datatypes are templated on
265  // Number, which allows the same code to work for real
266  // or complex numbers.
269 
270  // This vector will hold the degree of freedom indices for
271  // the element. These define where in the global system
272  // the element degrees of freedom get mapped.
273  std::vector<dof_id_type> dof_indices;
274 
275  // The global system matrix
276  SparseMatrix<Number> & matrix = system.get_system_matrix();
277  // The preconditioning matrix
278  auto & pre_matrix = system.get_matrix("preconditioner");
279 
280  // Now we will loop over all the elements in the mesh.
281  // We will compute the element matrix and right-hand-side
282  // contribution.
283  //
284  // Element ranges are a nice way to iterate through all the
285  // elements, or all the elements that have some property. The
286  // range will iterate from the first to the last element on
287  // the local processor.
288  // It is smart to make this one const so that we don't accidentally
289  // mess it up! In case users later modify this program to include
290  // refinement, we will be safe and will only consider the active
291  // elements; hence we use a variant of the
292  // active_local_element_ptr_range.
293  for (const auto & elem : mesh.active_local_element_ptr_range())
294  {
295  // Get the degree of freedom indices for the
296  // current element. These define where in the global
297  // matrix and right-hand-side this element will
298  // contribute to.
299  dof_map.dof_indices (elem, dof_indices);
300 
301  // Cache the number of degrees of freedom on this element, for
302  // use as a loop bound later. We use cast_int to explicitly
303  // convert from size() (which may be 64-bit) to unsigned int
304  // (which may be 32-bit but which is definitely enough to count
305  // *local* degrees of freedom.
306  const unsigned int n_dofs =
307  cast_int<unsigned int>(dof_indices.size());
308 
309  // Compute the element-specific data for the current
310  // element. This involves computing the location of the
311  // quadrature points (q_point) and the shape functions
312  // (phi, dphi) for the current element.
313  fe->reinit (elem);
314 
315  // With one variable, we should have the same number of degrees
316  // of freedom as shape functions.
317  libmesh_assert_equal_to (n_dofs, phi.size());
318 
319  // Zero the element matrix and right-hand side before
320  // summing them. We use the resize member here because
321  // the number of degrees of freedom might have changed from
322  // the last element. Note that this will be the case if the
323  // element type is different (i.e. the last element was a
324  // triangle, now we are on a quadrilateral).
325 
326  // The DenseMatrix::resize() and the DenseVector::resize()
327  // members will automatically zero out the matrix and vector.
328  Ke.resize (n_dofs, n_dofs);
329 
330  Fe.resize (n_dofs);
331 
332  // Now loop over the quadrature points. This handles
333  // the numeric integration.
334  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
335  {
336 
337  // Now we will build the element matrix. This involves
338  // a double loop to integrate the test functions (i) against
339  // the trial functions (j).
340  for (unsigned int i=0; i != n_dofs; i++)
341  for (unsigned int j=0; j != n_dofs; j++)
342  {
343  Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
344  }
345 
346  // This is the end of the matrix summation loop
347  // Now we build the element right-hand-side contribution.
348  // This involves a single loop in which we integrate the
349  // "forcing function" in the PDE against the test functions.
350  {
351  const Real x = q_point[qp](0);
352  const Real y = q_point[qp](1);
353  const Real eps = 1.e-3;
354 
355 
356  // "fxy" is the forcing function for the Poisson equation.
357  // In this case we set fxy to be a finite difference
358  // Laplacian approximation to the (known) exact solution.
359  //
360  // We will use the second-order accurate FD Laplacian
361  // approximation, which in 2D is
362  //
363  // u_xx + u_yy = (u(i,j-1) + u(i,j+1) +
364  // u(i-1,j) + u(i+1,j) +
365  // -4*u(i,j))/h^2
366  //
367  // Since the value of the forcing function depends only
368  // on the location of the quadrature point (q_point[qp])
369  // we will compute it here, outside of the i-loop
370  const Real fxy = -(exact_solution(x, y-eps) +
371  exact_solution(x, y+eps) +
372  exact_solution(x-eps, y) +
373  exact_solution(x+eps, y) -
374  4.*exact_solution(x, y))/eps/eps;
375 
376  for (unsigned int i=0; i != n_dofs; i++)
377  Fe(i) += JxW[qp]*fxy*phi[i][qp];
378  }
379  }
380 
381  // We have now reached the end of the RHS summation,
382  // and the end of quadrature point loop, so
383  // the interior element integration has
384  // been completed. However, we have not yet addressed
385  // boundary conditions. For this example we will only
386  // consider simple Dirichlet boundary conditions.
387  //
388  // There are several ways Dirichlet boundary conditions
389  // can be imposed. A simple approach, which works for
390  // interpolary bases like the standard Lagrange polynomials,
391  // is to assign function values to the
392  // degrees of freedom living on the domain boundary. This
393  // works well for interpolary bases, but is more difficult
394  // when non-interpolary (e.g Legendre or Hierarchic) bases
395  // are used.
396  //
397  // Dirichlet boundary conditions can also be imposed with a
398  // "penalty" method. In this case essentially the L2 projection
399  // of the boundary values are added to the matrix. The
400  // projection is multiplied by some large factor so that, in
401  // floating point arithmetic, the existing (smaller) entries
402  // in the matrix and right-hand-side are effectively ignored.
403  //
404  // This amounts to adding a term of the form (in latex notation)
405  //
406  // \frac{1}{\epsilon} \int_{\delta \Omega} \phi_i \phi_j = \frac{1}{\epsilon} \int_{\delta \Omega} u \phi_i
407  //
408  // where
409  //
410  // \frac{1}{\epsilon} is the penalty parameter, defined such that \epsilon << 1
411  {
412 
413  // The following loop is over the sides of the element.
414  // If the element has no neighbor on a side then that
415  // side MUST live on a boundary of the domain.
416  for (auto side : elem->side_index_range())
417  if (elem->neighbor_ptr(side) == nullptr)
418  {
419  // The value of the shape functions at the quadrature
420  // points.
421  const std::vector<std::vector<Real>> & phi_face = fe_face->get_phi();
422 
423  // The Jacobian * Quadrature Weight at the quadrature
424  // points on the face.
425  const std::vector<Real> & JxW_face = fe_face->get_JxW();
426 
427  // The XYZ locations (in physical space) of the
428  // quadrature points on the face. This is where
429  // we will interpolate the boundary value function.
430  const std::vector<Point> & qface_point = fe_face->get_xyz();
431 
432  // Compute the shape function values on the element
433  // face.
434  fe_face->reinit(elem, side);
435 
436  // Some shape functions will be 0 on the face, but for
437  // ease of indexing and generality of code we loop over
438  // them anyway
439  libmesh_assert_equal_to (n_dofs, phi_face.size());
440 
441  // Loop over the face quadrature points for integration.
442  for (unsigned int qp=0; qp<qface.n_points(); qp++)
443  {
444  // The location on the boundary of the current
445  // face quadrature point.
446  const Real xf = qface_point[qp](0);
447  const Real yf = qface_point[qp](1);
448 
449  // The penalty value. \frac{1}{\epsilon}
450  // in the discussion above.
451  const Real penalty = 1.e10;
452 
453  // The boundary value.
454  const Real value = exact_solution(xf, yf);
455 
456  // Matrix contribution of the L2 projection.
457  for (unsigned int i=0; i != n_dofs; i++)
458  for (unsigned int j=0; j != n_dofs; j++)
459  Ke(i,j) += JxW_face[qp]*penalty*phi_face[i][qp]*phi_face[j][qp];
460 
461  // Right-hand-side contribution of the L2
462  // projection.
463  for (unsigned int i=0; i != n_dofs; i++)
464  Fe(i) += JxW_face[qp]*penalty*value*phi_face[i][qp];
465  }
466  }
467  }
468 
469  // We have now finished the quadrature point loop,
470  // and have therefore applied all the boundary conditions.
471 
472  // If this assembly program were to be used on an adaptive mesh,
473  // we would have to apply any hanging node constraint equations
474  dof_map.constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
475 
476  // The element matrix and right-hand-side are now built
477  // for this element. Add them to the global matrix and
478  // right-hand-side vector. The SparseMatrix::add_matrix()
479  // and NumericVector::add_vector() members do this for us.
480  matrix.add_matrix(Ke, dof_indices);
481  pre_matrix.add_matrix(Ke, dof_indices);
482  system.rhs->add_vector (Fe, dof_indices);
483  }
484 
485  // All done!
486 }
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:196
This is the EquationSystems class.
This class provides a nice interface to PETSc&#39;s Vec object.
Definition: petsc_vector.h:73
unsigned int dim
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...
void resize(const unsigned int n)
Resize the vector.
Definition: dense_vector.h:396
virtual void add_vector(const T *v, const std::vector< numeric_index_type > &dof_indices)
Computes , where v is a pointer and each dof_indices[i] specifies where to add value v[i]...
Tnew cast_ptr(Told *oldvar)
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
MeshBase & mesh
NumericVector< Number > * rhs
The system matrix.
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.
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:90
The libMesh namespace provides an interface to certain functionality in the library.
const T_sys & get_system(std::string_view name) const
This is the MeshBase class.
Definition: mesh_base.h:75
SolverPackage default_solver_package()
Definition: libmesh.C:1117
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:179
virtual void add_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols)=0
Add the full matrix dm to the SparseMatrix.
void print_info(std::ostream &os=libMesh::out, const unsigned int verbosity=0, const bool global=true) const
Prints relevant information about the mesh.
Definition: mesh_base.C:1562
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
unsigned int add_variable(std::string_view var, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=nullptr)
Adds the variable var to the list of variables for this system.
Definition: system.C:1357
unsigned int n_points() const
Definition: quadrature.h:131
void assemble_poisson(EquationSystems &es, const std::string &system_name)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
OStreamProxy out
Real exact_solution(const Real x, const Real y, const Real z=0.)
This is the exact solution that we are trying to obtain.
const MeshBase & get_mesh() const
static const bool value
Definition: xdr_io.C:54
void resize(const unsigned int new_m, const unsigned int new_n)
Resizes the matrix to the specified size and calls zero().
Definition: dense_matrix.h:895
This class implements specific orders of Gauss quadrature.
unsigned int mesh_dimension() const
Definition: mesh_base.C:372
int main(int argc, char **argv)
virtual void init()
Initialize all the systems.
virtual System & add_system(std::string_view system_type, std::string_view name)
Add the system of type system_type named name to the systems array.
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50
const DofMap & get_dof_map() const
Definition: system.h:2374
const SparseMatrix< Number > & get_matrix(std::string_view mat_name) const
Definition: system.C:1124
const SparseMatrix< Number > & get_system_matrix() const