libMesh
eigenproblems_ex1.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 
19 
20 // <h1>Eigenproblems Example 1 - Solving an Eigen Problem</h1>
21 // \author Steffen Petersen
22 // \date 2005
23 //
24 // This example introduces the EigenSystem and shows
25 // how libMesh can be used for eigenvalue analysis.
26 //
27 // For solving eigen problems, libMesh interfaces
28 // SLEPc (www.grycap.upv.es/slepc/) which again is based on PETSc.
29 // Hence, this example will only work if the library is compiled
30 // with SLEPc support enabled.
31 //
32 // In this example some eigenvalues for a standard symmetric eigenvalue
33 // problem A*x=lambda*x are computed, where the matrix A
34 // is assembled according to a mass matrix.
35 
36 
37 // libMesh include files.
38 #include "libmesh/libmesh.h"
39 #include "libmesh/mesh.h"
40 #include "libmesh/mesh_generation.h"
41 #include "libmesh/exodusII_io.h"
42 #include "libmesh/eigen_system.h"
43 #include "libmesh/equation_systems.h"
44 #include "libmesh/fe.h"
45 #include "libmesh/quadrature_gauss.h"
46 #include "libmesh/dense_matrix.h"
47 #include "libmesh/sparse_matrix.h"
48 #include "libmesh/numeric_vector.h"
49 #include "libmesh/dof_map.h"
50 #include "libmesh/getpot.h"
51 
52 // Bring in everything from the libMesh namespace
53 using namespace libMesh;
54 
55 // Function prototype. This is the function that will assemble
56 // the eigen system. Here, we will simply assemble a mass matrix.
58  const std::string & system_name);
59 
60 int main (int argc, char ** argv)
61 {
62  // Initialize libMesh and the dependent libraries.
63  LibMeshInit init (argc, argv);
64 
65 #ifdef LIBMESH_DEFAULT_SINGLE_PRECISION
66  // SLEPc currently gives us a nasty crash with Real==float
67  libmesh_example_requires(false, "--disable-singleprecision");
68 #endif
69 
70 #ifndef LIBMESH_HAVE_SLEPC
71  libmesh_example_requires(false, "--enable-slepc");
72 #else
73  // Check for proper usage.
74  libmesh_error_msg_if(argc < 3, "\nUsage: " << argv[0] << " -n <number of eigen values>");
75 
76  // Tell the user what we are doing.
77  libMesh::out << "Running " << argv[0];
78 
79  for (int i=1; i<argc; i++)
80  libMesh::out << " " << argv[i];
81 
82  libMesh::out << std::endl << std::endl;
83 
84  // Get the number of eigen values to be computed from "-n", and
85  // possibly the mesh size from -nx and -ny
86 
87  const int nev = libMesh::command_line_next("-n", 5),
88  nx = libMesh::command_line_next("-nx", 20),
89  ny = libMesh::command_line_next("-ny", 20);
90 
91  // Skip this 2D example if libMesh was compiled as 1D-only.
92  libmesh_example_requires(2 <= LIBMESH_DIM, "2D support");
93 
94  // Create a mesh, with dimension to be overridden later, on the
95  // default MPI communicator.
96  Mesh mesh(init.comm());
97 
98  // Use the internal mesh generator to create a uniform
99  // 2D grid on a square.
101  nx, ny,
102  -1., 1.,
103  -1., 1.,
104  QUAD4);
105 
106  // Print information about the mesh to the screen.
107  mesh.print_info();
108 
109  // Create an equation systems object.
110  EquationSystems equation_systems (mesh);
111 
112  // Create a EigenSystem named "Eigensystem" and (for convenience)
113  // use a reference to the system we create.
114  EigenSystem & eigen_system =
115  equation_systems.add_system<EigenSystem> ("Eigensystem");
116 
117  // Declare the system variables.
118  // Adds the variable "p" to "Eigensystem". "p"
119  // will be approximated using second-order approximation.
120  eigen_system.add_variable("p", FIRST);
121 
122  // Give the system a pointer to the matrix assembly
123  // function defined below.
125 
126  // Set necessary parameters used in EigenSystem::solve(),
127  // i.e. the number of requested eigenpairs nev and the number
128  // of basis vectors ncv used in the solution algorithm. Note that
129  // ncv >= nev must hold and ncv >= 2*nev is recommended.
130  equation_systems.parameters.set<unsigned int>("eigenpairs") = nev;
131  equation_systems.parameters.set<unsigned int>("basis vectors") = nev*3;
132 
133  // You may optionally change the default eigensolver used by SLEPc.
134  // The Krylov-Schur method is mathematically equivalent to implicitly
135  // restarted Arnoldi, the method of Arpack, so there is currently no
136  // point in using SLEPc with Arpack.
137  // ARNOLDI = default in SLEPc 2.3.1 and earlier
138  // KRYLOVSCHUR default in SLEPc 2.3.2 and later
139  // eigen_system.get_eigen_solver().set_eigensolver_type(KRYLOVSCHUR);
140 
141  // Set the solver tolerance and the maximum number of iterations.
142  equation_systems.parameters.set<Real>
143  ("linear solver tolerance") = pow(TOLERANCE, 5./3.);
144  equation_systems.parameters.set<unsigned int>
145  ("linear solver maximum iterations") = 1000;
146 
147  // Initialize the data structures for the equation system.
148  equation_systems.init();
149 
150  // Prints information about the system to the screen.
151  equation_systems.print_info();
152 
153  // Solve the system "Eigensystem".
154  eigen_system.solve();
155 
156  // Get the number of converged eigen pairs.
157  unsigned int nconv = eigen_system.get_n_converged();
158 
159  libMesh::out << "Number of converged eigenpairs: " << nconv
160  << "\n" << std::endl;
161 
162  // Get the last converged eigenpair
163  if (nconv != 0)
164  {
165  eigen_system.get_eigenpair(nconv-1);
166 
167 #ifdef LIBMESH_HAVE_EXODUS_API
168  // Write the eigen vector to file.
169  ExodusII_IO (mesh).write_equation_systems ("out.e", equation_systems);
170 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
171  }
172  else
173  {
174  libMesh::out << "WARNING: Solver did not converge!\n" << nconv << std::endl;
175  }
176 
177  // All done.
178  return 0;
179 #endif // LIBMESH_HAVE_SLEPC
180 }
181 
182 
183 
185  const std::string & libmesh_dbg_var(system_name))
186 {
187  // It is a good idea to make sure we are assembling
188  // the proper system.
189  libmesh_assert_equal_to (system_name, "Eigensystem");
190 
191 #ifdef LIBMESH_HAVE_SLEPC
192 
193  // Get a constant reference to the mesh object.
194  const MeshBase & mesh = es.get_mesh();
195 
196  // The dimension that we are running.
197  const unsigned int dim = mesh.mesh_dimension();
198 
199  // Get a reference to our system.
200  EigenSystem & eigen_system = es.get_system<EigenSystem> ("Eigensystem");
201 
202  // Get a constant reference to the Finite Element type
203  // for the first (and only) variable in the system.
204  FEType fe_type = eigen_system.get_dof_map().variable_type(0);
205 
206  // A reference to the system matrix
207  SparseMatrix<Number> & matrix_A = eigen_system.get_matrix_A();
208 
209  // Build a Finite Element object of the specified type. Since the
210  // FEBase::build() member dynamically creates memory we will
211  // store the object as a std::unique_ptr<FEBase>. This can be thought
212  // of as a pointer that will clean up after itself.
213  std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
214 
215  // A Gauss quadrature rule for numerical integration.
216  // Use the default quadrature order.
217  QGauss qrule (dim, fe_type.default_quadrature_order());
218 
219  // Tell the finite element object to use our quadrature rule.
220  fe->attach_quadrature_rule (&qrule);
221 
222  // The element Jacobian * quadrature weight at each integration point.
223  const std::vector<Real> & JxW = fe->get_JxW();
224 
225  // The element shape functions evaluated at the quadrature points.
226  const std::vector<std::vector<Real>> & phi = fe->get_phi();
227 
228  // A reference to the DofMap object for this system. The DofMap
229  // object handles the index translation from node and element numbers
230  // to degree of freedom numbers.
231  const DofMap & dof_map = eigen_system.get_dof_map();
232 
233  // The element mass matrix.
235 
236  // This vector will hold the degree of freedom indices for
237  // the element. These define where in the global system
238  // the element degrees of freedom get mapped.
239  std::vector<dof_id_type> dof_indices;
240 
241  // Now we will loop over all the elements in the mesh that
242  // live on the local processor. We will compute the element
243  // matrix and right-hand-side contribution. In case users
244  // later modify this program to include refinement, we will
245  // be safe and will only consider the active elements;
246  // hence we use a variant of the active_elem_iterator.
247  for (const auto & elem : mesh.active_local_element_ptr_range())
248  {
249  // Get the degree of freedom indices for the
250  // current element. These define where in the global
251  // matrix and right-hand-side this element will
252  // contribute to.
253  dof_map.dof_indices (elem, dof_indices);
254 
255  // Compute the element-specific data for the current
256  // element. This involves computing the location of the
257  // quadrature points (q_point) and the shape functions
258  // (phi, dphi) for the current element.
259  fe->reinit (elem);
260 
261  // Zero the element matrices and rhs before
262  // summing them. We use the resize member here because
263  // the number of degrees of freedom might have changed from
264  // the last element. Note that this will be the case if the
265  // element type is different (i.e. the last element was a
266  // triangle, now we are on a quadrilateral).
267  const unsigned int n_dofs =
268  cast_int<unsigned int>(dof_indices.size());
269  Me.resize (n_dofs, n_dofs);
270 
271  // Now loop over the quadrature points. This handles
272  // the numeric integration.
273  //
274  // We will build the element matrix. This involves
275  // a double loop to integrate the test functions (i) against
276  // the trial functions (j).
277  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
278  for (unsigned int i=0; i != n_dofs; i++)
279  for (unsigned int j=0; j != n_dofs; j++)
280  Me(i,j) += JxW[qp]*phi[i][qp]*phi[j][qp];
281 
282  // On an unrefined mesh, constrain_element_matrix does
283  // nothing. If this assembly function is ever repurposed to
284  // run on a refined mesh, getting the hanging node constraints
285  // right will be important. Note that, even with
286  // asymmetric_constraint_rows = false, the constrained dof
287  // diagonals still exist in the matrix, with diagonal entries
288  // that are there to ensure non-singular matrices for linear
289  // solves but which would generate positive non-physical
290  // eigenvalues for eigensolves.
291  dof_map.constrain_element_matrix(Me, dof_indices, false);
292 
293  // Finally, simply add the element contribution to the
294  // overall matrix.
295  matrix_A.add_matrix (Me, dof_indices);
296  } // end of element loop
297 
298 #else
299  // Avoid compiler warnings
300  libmesh_ignore(es);
301 #endif // LIBMESH_HAVE_SLEPC
302 }
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:196
virtual std::pair< Real, Real > get_eigenpair(dof_id_type i)
Definition: eigen_system.C:308
T command_line_next(std::string name, T default_value)
Use GetPot&#39;s search()/next() functions to get following arguments from the command line...
Definition: libmesh.C:1078
This is the EquationSystems class.
const SparseMatrix< Number > & get_matrix_A() const
Definition: eigen_system.C:333
void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
Definition: dof_map.C:2164
static constexpr Real TOLERANCE
unsigned int dim
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Definition: exodusII_io.h:52
const FEType & variable_type(const unsigned int c) const
Definition: dof_map.h:2220
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
MeshBase & mesh
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
int main(int argc, char **argv)
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:179
void libmesh_ignore(const Args &...)
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.
T pow(const T &x)
Definition: utility.h:328
virtual void write_equation_systems(const std::string &fname, const EquationSystems &es, const std::set< std::string > *system_names=nullptr) override
Writes out the solution for no specific time or timestep.
Definition: exodusII_io.C:2033
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
void attach_assemble_function(void fptr(EquationSystems &es, const std::string &name))
Register a user function to use in assembling the system matrix and RHS.
Definition: system.C:2161
virtual void solve() override
Assembles & solves the eigen system.
Definition: eigen_system.C:287
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
T & set(const std::string &)
Definition: parameters.h:469
OStreamProxy out
const MeshBase & get_mesh() const
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
Parameters parameters
Data structure holding arbitrary parameters.
unsigned int get_n_converged() const
Definition: eigen_system.h:130
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.
Manages consistently variables, degrees of freedom, and coefficient vectors for eigenvalue problems...
Definition: eigen_system.h:55
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
void assemble_mass(EquationSystems &es, const std::string &system_name)
void constrain_element_matrix(DenseMatrix< Number > &matrix, std::vector< dof_id_type > &elem_dofs, bool asymmetric_constraint_rows=true) const
Constrains the element matrix.
Definition: dof_map.h:2317