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