libMesh
eigenproblems_ex2.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 2 - Solving a generalized Eigen Problem</h1>
21 // \author Steffen Petersen
22 // \date 2006
23 //
24 // This example shows how the previous EigenSolver example
25 // can be adapted to solve generalized eigenvalue problems.
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 generalized symmetric
33 // eigenvalue problem A*x=lambda*B*x are computed, where the
34 // matrices A and B are assembled according to stiffness and
35 // mass matrix, respectively.
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/slepc_eigen_solver.h"
45 #include "libmesh/fe.h"
46 #include "libmesh/quadrature_gauss.h"
47 #include "libmesh/dense_matrix.h"
48 #include "libmesh/sparse_matrix.h"
49 #include "libmesh/numeric_vector.h"
50 #include "libmesh/dof_map.h"
51 #include "libmesh/enum_eigen_solver_type.h"
52 #include "libmesh/getpot.h"
53 
54 // Bring in everything from the libMesh namespace
55 using namespace libMesh;
56 
57 
58 // Function prototype. This is the function that will assemble
59 // the eigen system. Here, we will simply assemble a mass matrix.
61  const std::string & system_name);
62 
63 
64 
65 int main (int argc, char ** argv)
66 {
67  // Initialize libMesh and the dependent libraries.
68  LibMeshInit init (argc, argv);
69 
70  // This example is designed for the SLEPc eigen solver interface.
71 #ifndef LIBMESH_HAVE_SLEPC
72  if (init.comm().rank() == 0)
73  libMesh::err << "ERROR: This example requires libMesh to be\n"
74  << "compiled with SLEPc eigen solvers support!"
75  << std::endl;
76 
77  return 0;
78 #else
79 
80 #ifdef LIBMESH_DEFAULT_SINGLE_PRECISION
81  // SLEPc currently gives us a nasty crash with Real==float
82  libmesh_example_requires(false, "--disable-singleprecision");
83 #endif
84 
85  // Check for proper usage.
86  libmesh_error_msg_if(argc < 3, "\nUsage: " << argv[0] << " -n <number of eigen values>");
87 
88  // Tell the user what we are doing.
89  libMesh::out << "Running " << argv[0];
90 
91  for (int i=1; i<argc; i++)
92  libMesh::out << " " << argv[i];
93 
94  libMesh::out << std::endl << std::endl;
95 
96  // Get the number of eigen values to be computed from "-n", and
97  // possibly the mesh size from -nx and -ny
98 
99  const int nev = libMesh::command_line_next("-n", 5),
100  nx = libMesh::command_line_next("-nx", 20),
101  ny = libMesh::command_line_next("-ny", 20);
102 
103 
104  // Skip this 2D example if libMesh was compiled as 1D-only.
105  libmesh_example_requires(2 <= LIBMESH_DIM, "2D support");
106 
107  // Create a mesh, with dimension to be overridden later, on the
108  // default MPI communicator.
109  Mesh mesh(init.comm());
110 
111  // Use the internal mesh generator to create a uniform
112  // 2D grid on a square.
114  nx, ny,
115  -1., 1.,
116  -1., 1.,
117  QUAD4);
118 
119  // Print information about the mesh to the screen.
120  mesh.print_info();
121 
122  // Create an equation systems object.
123  EquationSystems equation_systems (mesh);
124 
125  // Create a EigenSystem named "Eigensystem" and (for convenience)
126  // use a reference to the system we create.
127  EigenSystem & eigen_system =
128  equation_systems.add_system<EigenSystem> ("Eigensystem");
129 
130  // Declare the system variables.
131  // Adds the variable "p" to "Eigensystem". "p"
132  // will be approximated using second-order approximation.
133  eigen_system.add_variable("p", FIRST);
134 
135  // Give the system a pointer to the matrix assembly
136  // function defined below.
138 
139  // Set necessary parameters used in EigenSystem::solve(),
140  // i.e. the number of requested eigenpairs nev and the number
141  // of basis vectors ncv used in the solution algorithm. Note that
142  // ncv >= nev must hold and ncv >= 2*nev is recommended.
143  equation_systems.parameters.set<unsigned int>("eigenpairs") = nev;
144  equation_systems.parameters.set<unsigned int>("basis vectors") = nev*3;
145 
146  // You may optionally change the default eigensolver used by SLEPc.
147  // The Krylov-Schur method is mathematically equivalent to implicitly
148  // restarted Arnoldi, the method of Arpack, so there is currently no
149  // point in using SLEPc with Arpack.
150  // ARNOLDI = default in SLEPc 2.3.1 and earlier
151  // KRYLOVSCHUR default in SLEPc 2.3.2 and later
152  // eigen_system.get_eigen_solver().set_eigensolver_type(KRYLOVSCHUR);
153 
154  // Set the solver tolerance and the maximum number of iterations.
155  equation_systems.parameters.set<Real>("linear solver tolerance") = pow(TOLERANCE, 5./3.);
156  equation_systems.parameters.set<unsigned int>("linear solver maximum iterations") = 1000;
157 
158  // Set the type of the problem, here we deal with
159  // a generalized Hermitian problem.
160  eigen_system.set_eigenproblem_type(GHEP);
161 
162  // Set the eigenvalues to be computed. Note that not
163  // all solvers in SLEPc support this capability.
164  eigen_system.get_eigen_solver().set_position_of_spectrum(2.3);
165 
166  // Initialize the data structures for the equation system.
167  equation_systems.init();
168 
169  // Prints information about the system to the screen.
170  equation_systems.print_info();
171 
172 #if SLEPC_VERSION_LESS_THAN(3,1,0)
173  libmesh_error_msg("SLEPc 3.1 is required to call EigenSolver::set_initial_space()");
174 #else
175  // Get the SLEPc solver object and set initial guess for one basis vector
176  // this has to be done _after_ the EquationSystems object is initialized
177  EigenSolver<Number> & slepc_eps = eigen_system.get_eigen_solver();
178  NumericVector<Number> & initial_space = eigen_system.add_vector("initial_space");
179  initial_space.add(1.0);
180  slepc_eps.set_initial_space(initial_space);
181 #endif
182 
183  // Solve the system "Eigensystem".
184  eigen_system.solve();
185 
186  // Get the number of converged eigen pairs.
187  unsigned int nconv = eigen_system.get_n_converged();
188 
189  libMesh::out << "Number of converged eigenpairs: "
190  << nconv
191  << "\n"
192  << std::endl;
193 
194  // Get the last converged eigenpair
195  if (nconv != 0)
196  {
197  eigen_system.get_eigenpair(nconv-1);
198 
199 #ifdef LIBMESH_HAVE_EXODUS_API
200  // Write the eigen vector to file.
201  ExodusII_IO (mesh).write_equation_systems ("out.e", equation_systems);
202 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
203  }
204  else
205  {
206  libMesh::out << "WARNING: Solver did not converge!\n" << nconv << std::endl;
207  }
208 
209 #endif // LIBMESH_HAVE_SLEPC
210 
211  // All done.
212  return 0;
213 }
214 
215 
216 
218  const std::string & libmesh_dbg_var(system_name))
219 {
220 
221  // It is a good idea to make sure we are assembling
222  // the proper system.
223  libmesh_assert_equal_to (system_name, "Eigensystem");
224 
225 #ifdef LIBMESH_HAVE_SLEPC
226 
227  // Get a constant reference to the mesh object.
228  const MeshBase & mesh = es.get_mesh();
229 
230  // The dimension that we are running.
231  const unsigned int dim = mesh.mesh_dimension();
232 
233  // Get a reference to our system.
234  EigenSystem & eigen_system = es.get_system<EigenSystem> ("Eigensystem");
235 
236  // Get a constant reference to the Finite Element type
237  // for the first (and only) variable in the system.
238  FEType fe_type = eigen_system.get_dof_map().variable_type(0);
239 
240  // A reference to the two system matrices
241  SparseMatrix<Number> & matrix_A = eigen_system.get_matrix_A();
242  SparseMatrix<Number> & matrix_B = eigen_system.get_matrix_B();
243 
244  // Build a Finite Element object of the specified type. Since the
245  // FEBase::build() member dynamically creates memory we will
246  // store the object as a std::unique_ptr<FEBase>. This can be thought
247  // of as a pointer that will clean up after itself.
248  std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
249 
250  // A Gauss quadrature rule for numerical integration.
251  // Use the default quadrature order.
252  QGauss qrule (dim, fe_type.default_quadrature_order());
253 
254  // Tell the finite element object to use our quadrature rule.
255  fe->attach_quadrature_rule (&qrule);
256 
257  // The element Jacobian * quadrature weight at each integration point.
258  const std::vector<Real> & JxW = fe->get_JxW();
259 
260  // The element shape functions evaluated at the quadrature points.
261  const std::vector<std::vector<Real>> & phi = fe->get_phi();
262 
263  // The element shape function gradients evaluated at the quadrature
264  // points.
265  const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
266 
267  // A reference to the DofMap object for this system. The DofMap
268  // object handles the index translation from node and element numbers
269  // to degree of freedom numbers.
270  const DofMap & dof_map = eigen_system.get_dof_map();
271 
272  // The element mass and stiffness matrices.
275 
276  // This vector will hold the degree of freedom indices for
277  // the element. These define where in the global system
278  // the element degrees of freedom get mapped.
279  std::vector<dof_id_type> dof_indices;
280 
281  // Now we will loop over all the elements in the mesh that
282  // live on the local processor. We will compute the element
283  // matrix and right-hand-side contribution. In case users
284  // later modify this program to include refinement, we will
285  // be safe and will only consider the active elements;
286  // hence we use a variant of the active_elem_iterator.
287  for (const auto & elem : mesh.active_local_element_ptr_range())
288  {
289  // Get the degree of freedom indices for the
290  // current element. These define where in the global
291  // matrix and right-hand-side this element will
292  // contribute to.
293  dof_map.dof_indices (elem, dof_indices);
294 
295  // Compute the element-specific data for the current
296  // element. This involves computing the location of the
297  // quadrature points (q_point) and the shape functions
298  // (phi, dphi) for the current element.
299  fe->reinit (elem);
300 
301  // Zero the element matrices before
302  // summing them. We use the resize member here because
303  // the number of degrees of freedom might have changed from
304  // the last element. Note that this will be the case if the
305  // element type is different (i.e. the last element was a
306  // triangle, now we are on a quadrilateral).
307  const unsigned int n_dofs =
308  cast_int<unsigned int>(dof_indices.size());
309  Ke.resize (n_dofs, n_dofs);
310  Me.resize (n_dofs, n_dofs);
311 
312  // Now loop over the quadrature points. This handles
313  // the numeric integration.
314  //
315  // We will build the element matrix. This involves
316  // a double loop to integrate the test functions (i) against
317  // the trial functions (j).
318  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
319  for (unsigned int i=0; i<n_dofs; i++)
320  for (unsigned int j=0; j<n_dofs; j++)
321  {
322  Me(i,j) += JxW[qp]*phi[i][qp]*phi[j][qp];
323  Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
324  }
325 
326  // On an unrefined mesh, constrain_element_matrix does
327  // nothing. If this assembly function is ever repurposed to
328  // run on a refined mesh, getting the hanging node constraints
329  // right will be important. Note that, even with
330  // asymmetric_constraint_rows = false, the constrained dof
331  // diagonals still exist in the matrix, with diagonal entries
332  // that are there to ensure non-singular matrices for linear
333  // solves but which would generate positive non-physical
334  // eigenvalues for eigensolves.
335  // dof_map.constrain_element_matrix(Ke, dof_indices, false);
336  // dof_map.constrain_element_matrix(Me, dof_indices, false);
337 
338  // Finally, simply add the element contribution to the
339  // overall matrices A and B.
340  matrix_A.add_matrix (Ke, dof_indices);
341  matrix_B.add_matrix (Me, dof_indices);
342  } // end of element loop
343 
344 
345 #else
346  // Avoid compiler warnings
347  libmesh_ignore(es);
348 #endif // LIBMESH_HAVE_SLEPC
349 }
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:196
OStreamProxy err
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.
virtual void set_initial_space(NumericVector< T > &initial_space_in)=0
Provide one basis vector for the initial guess.
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 SparseMatrix< Number > & get_matrix_B() const
Definition: eigen_system.C:351
void assemble_mass(EquationSystems &es, const std::string &system_name)
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
NumericVector< Number > & add_vector(std::string_view vec_name, const bool projections=true, const ParallelType type=PARALLEL)
Adds the additional vector vec_name to this system.
Definition: system.C:768
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 EigenSolver< Number > & get_eigen_solver() const
Definition: eigen_system.C:482
const T_sys & get_system(std::string_view name) const
This is the MeshBase class.
Definition: mesh_base.h:75
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
int main(int argc, char **argv)
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.
void set_eigenproblem_type(EigenProblemType ept)
Sets the type of the current eigen problem.
Definition: eigen_system.C:85
unsigned int get_n_converged() const
Definition: eigen_system.h:130
virtual void init()
Initialize all the systems.
virtual void add(const numeric_index_type i, const T value)=0
Adds value to the vector entry specified by i.
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
This class provides an interface to solvers for eigenvalue problems.
Definition: eigen_solver.h:58