Go to the documentation of this file.
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"
60 const std::string & system_name);
64 int main (
int argc,
char ** argv)
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!"
79 #ifdef LIBMESH_DEFAULT_SINGLE_PRECISION
81 libmesh_example_requires(
false,
"--disable-singleprecision");
86 libmesh_error_msg(
"\nUsage: " << argv[0] <<
" -n <number of eigen values>");
93 for (
int i=1; i<argc; i++)
100 const unsigned int nev = std::atoi(argv[2]);
103 libmesh_example_requires(2 <= LIBMESH_DIM,
"2D support");
141 equation_systems.
parameters.
set<
unsigned int>(
"eigenpairs") = nev;
142 equation_systems.
parameters.
set<
unsigned int>(
"basis vectors") = nev*3;
154 equation_systems.
parameters.
set<
unsigned int>(
"linear solver maximum iterations") = 1000;
162 eigen_system.
eigen_solver->set_position_of_spectrum(2.3);
165 equation_systems.
init();
170 #if SLEPC_VERSION_LESS_THAN(3,1,0)
171 libmesh_error_msg(
"SLEPc 3.1 is required to call EigenSolver::set_initial_space()");
175 std::unique_ptr<EigenSolver<Number>> & slepc_eps = eigen_system.
eigen_solver;
177 initial_space.
add(1.0);
178 slepc_eps->set_initial_space(initial_space);
182 eigen_system.
solve();
197 #ifdef LIBMESH_HAVE_EXODUS_API
200 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
204 libMesh::out <<
"WARNING: Solver did not converge!\n" << nconv << std::endl;
207 #endif // LIBMESH_HAVE_SLEPC
216 const std::string & libmesh_dbg_var(system_name))
221 libmesh_assert_equal_to (system_name,
"Eigensystem");
223 #ifdef LIBMESH_HAVE_SLEPC
250 QGauss qrule (
dim, fe_type.default_quadrature_order());
253 fe->attach_quadrature_rule (&qrule);
256 const std::vector<Real> & JxW = fe->get_JxW();
259 const std::vector<std::vector<Real>> & phi = fe->get_phi();
263 const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
277 std::vector<dof_id_type> dof_indices;
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);
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++)
320 Me(i,j) += JxW[qp]*phi[i][qp]*phi[j][qp];
321 Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
346 #endif // LIBMESH_HAVE_SLEPC
virtual void add(const numeric_index_type i, const T value)=0
Adds value to each entry of the vector.
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
std::unique_ptr< SparseMatrix< Number > > matrix_A
The system matrix for standard eigenvalue problems.
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.
const MeshBase & get_mesh() const
This class implements specific orders of Gauss quadrature.
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.
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
virtual void solve() override
Assembles & solves the eigen system.
The libMesh namespace provides an interface to certain functionality in the library.
const T_sys & get_system(const std::string &name) const
static const Real TOLERANCE
std::unique_ptr< EigenSolver< Number > > eigen_solver
The EigenSolver, defining which interface, i.e solver package to use.
void set_eigenproblem_type(EigenProblemType ept)
Sets the type of the current eigen problem.
unsigned int mesh_dimension() const
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
void resize(const unsigned int new_m, const unsigned int new_n)
Resize the matrix.
unsigned int get_n_converged() const
std::unique_ptr< SparseMatrix< Number > > matrix_B
A second system matrix for generalized eigenvalue problems.
int main(int argc, char **argv)
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
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.
This is the MeshBase class.
void libmesh_ignore(const Args &...)
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
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.
virtual void init()
Initialize all the systems.
virtual std::pair< Real, Real > get_eigenpair(dof_id_type i)
double pow(double a, int b)
const FEType & variable_type(const unsigned int c) const
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
This is the EquationSystems class.
T & set(const std::string &)
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 ...
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.
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Manages consistently variables, degrees of freedom, and coefficient vectors for eigenvalue problems.
This class handles the numbering of degrees of freedom on a mesh.
void print_info(std::ostream &os=libMesh::out) const
Prints relevant information about the mesh.
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.
const DofMap & get_dof_map() const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void assemble_mass(EquationSystems &es, const std::string &system_name)
Parameters parameters
Data structure holding arbitrary parameters.