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/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"
57 const std::string & system_name);
59 int main (
int argc,
char ** argv)
64 #ifdef LIBMESH_DEFAULT_SINGLE_PRECISION
66 libmesh_example_requires(
false,
"--disable-singleprecision");
69 #ifndef LIBMESH_HAVE_SLEPC
70 libmesh_example_requires(
false,
"--enable-slepc");
74 libmesh_error_msg(
"\nUsage: " << argv[0] <<
" -n <number of eigen values>");
81 for (
int i=1; i<argc; i++)
88 const unsigned int nev = std::atoi(argv[2]);
91 libmesh_example_requires(2 <= LIBMESH_DIM,
"2D support");
129 equation_systems.
parameters.
set<
unsigned int>(
"eigenpairs") = nev;
130 equation_systems.
parameters.
set<
unsigned int>(
"basis vectors") = nev*3;
144 (
"linear solver maximum iterations") = 1000;
147 equation_systems.
init();
153 eigen_system.
solve();
158 libMesh::out <<
"Number of converged eigenpairs: " << nconv
159 <<
"\n" << std::endl;
166 #ifdef LIBMESH_HAVE_EXODUS_API
169 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
173 libMesh::out <<
"WARNING: Solver did not converge!\n" << nconv << std::endl;
178 #endif // LIBMESH_HAVE_SLEPC
184 const std::string & libmesh_dbg_var(system_name))
188 libmesh_assert_equal_to (system_name,
"Eigensystem");
190 #ifdef LIBMESH_HAVE_SLEPC
216 QGauss qrule (
dim, fe_type.default_quadrature_order());
219 fe->attach_quadrature_rule (&qrule);
222 const std::vector<Real> & JxW = fe->get_JxW();
225 const std::vector<std::vector<Real>> & phi = fe->get_phi();
238 std::vector<dof_id_type> dof_indices;
266 const unsigned int n_dofs =
267 cast_int<unsigned int>(dof_indices.size());
268 Me.
resize (n_dofs, n_dofs);
276 for (
unsigned int qp=0; qp<qrule.n_points(); qp++)
277 for (
unsigned int i=0; i != n_dofs; i++)
278 for (
unsigned int j=0; j != n_dofs; j++)
279 Me(i,j) += JxW[qp]*phi[i][qp]*phi[j][qp];
300 #endif // LIBMESH_HAVE_SLEPC
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.
void assemble_mass(EquationSystems &es, const std::string &system_name)
const T_sys & get_system(const std::string &name) const
static const Real TOLERANCE
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
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 ...
void constrain_element_matrix(DenseMatrix< Number > &matrix, std::vector< dof_id_type > &elem_dofs, bool asymmetric_constraint_rows=true) const
Constrains the element 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.
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.
int main(int argc, char **argv)
void print_info(std::ostream &os=libMesh::out) const
Prints relevant information about the mesh.
const DofMap & get_dof_map() const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Parameters parameters
Data structure holding arbitrary parameters.