Go to the source code of this file.
◆ assemble_mass() [1/2]
void assemble_mass |
( |
EquationSystems & |
es, |
|
|
const std::string & |
libmesh_dbg_varsystem_name |
|
) |
| |
Definition at line 183 of file eigenproblems_ex1.C.
188 libmesh_assert_equal_to (system_name,
"Eigensystem");
190 #ifdef LIBMESH_HAVE_SLEPC
212 std::unique_ptr<FEBase> fe (FEBase::build(
dim, fe_type));
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
References libMesh::MeshBase::active_local_element_ptr_range(), libMesh::SparseMatrix< T >::add_matrix(), libMesh::FEGenericBase< OutputType >::build(), libMesh::DofMap::constrain_element_matrix(), dim, libMesh::DofMap::dof_indices(), libMesh::System::get_dof_map(), libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), libMesh::libmesh_ignore(), libMesh::EigenSystem::matrix_A, mesh, libMesh::MeshBase::mesh_dimension(), libMesh::DenseMatrix< T >::resize(), and libMesh::DofMap::variable_type().
◆ assemble_mass() [2/2]
void assemble_mass |
( |
EquationSystems & |
es, |
|
|
const std::string & |
system_name |
|
) |
| |
◆ main()
int main |
( |
int |
argc, |
|
|
char ** |
argv |
|
) |
| |
Definition at line 59 of file eigenproblems_ex1.C.
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");
114 equation_systems.add_system<
EigenSystem> (
"Eigensystem");
129 equation_systems.parameters.set<
unsigned int>(
"eigenpairs") = nev;
130 equation_systems.parameters.set<
unsigned int>(
"basis vectors") = nev*3;
141 equation_systems.parameters.set<
Real>
143 equation_systems.parameters.set<
unsigned int>
144 (
"linear solver maximum iterations") = 1000;
147 equation_systems.init();
150 equation_systems.print_info();
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
References libMesh::EquationSystems::add_system(), libMesh::System::add_variable(), assemble_mass(), libMesh::System::attach_assemble_function(), libMesh::MeshTools::Generation::build_square(), libMesh::FIRST, libMesh::EigenSystem::get_eigenpair(), libMesh::EigenSystem::get_n_converged(), libMesh::TriangleWrapper::init(), libMesh::EquationSystems::init(), mesh, libMesh::out, libMesh::EquationSystems::parameters, std::pow(), libMesh::EquationSystems::print_info(), libMesh::MeshBase::print_info(), libMesh::QUAD4, libMesh::Real, libMesh::Parameters::set(), libMesh::EigenSystem::solve(), libMesh::TOLERANCE, and libMesh::MeshOutput< MT >::write_equation_systems().
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.
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.
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 &...)
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 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.
This is the EquationSystems class.
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.
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