Go to the source code of this file.
◆ assemble_matrices() [1/2]
void assemble_matrices |
( |
EquationSystems & |
es, |
|
|
const std::string & |
libmesh_dbg_varsystem_name |
|
) |
| |
Definition at line 291 of file eigenproblems_ex3.C.
297 libmesh_assert_equal_to (system_name,
"Eigensystem");
299 #ifdef LIBMESH_HAVE_SLEPC
322 std::unique_ptr<FEBase> fe (FEBase::build(
dim, fe_type));
326 QGauss qrule (
dim, fe_type.default_quadrature_order());
329 fe->attach_quadrature_rule (&qrule);
332 const std::vector<Real> & JxW = fe->get_JxW();
335 const std::vector<std::vector<Real>> & phi = fe->get_phi();
339 const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
353 std::vector<dof_id_type> dof_indices;
382 const unsigned int n_dofs =
383 cast_int<unsigned int>(dof_indices.size());
384 Ke.
resize (n_dofs, n_dofs);
385 Me.
resize (n_dofs, n_dofs);
393 for (
unsigned int qp=0; qp<qrule.n_points(); qp++)
394 for (
unsigned int i=0; i<n_dofs; i++)
395 for (
unsigned int j=0; j<n_dofs; j++)
397 Me(i,j) += JxW[qp]*phi[i][qp]*phi[j][qp];
398 Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
422 #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, libMesh::EigenSystem::matrix_B, mesh, libMesh::MeshBase::mesh_dimension(), libMesh::DenseMatrix< T >::resize(), and libMesh::DofMap::variable_type().
◆ assemble_matrices() [2/2]
void assemble_matrices |
( |
EquationSystems & |
es, |
|
|
const std::string & |
system_name |
|
) |
| |
◆ get_dirichlet_dofs()
void get_dirichlet_dofs |
( |
EquationSystems & |
es, |
|
|
const std::string & |
system_name, |
|
|
std::set< unsigned int > & |
global_dirichlet_dofs_set |
|
) |
| |
◆ main()
int main |
( |
int |
argc, |
|
|
char ** |
argv |
|
) |
| |
Definition at line 90 of file eigenproblems_ex3.C.
96 #ifndef LIBMESH_HAVE_EXODUS_API
97 libmesh_example_requires(
false,
"--enable-exodus");
101 #ifndef LIBMESH_HAVE_SLEPC
102 if (
init.comm().rank() == 0)
103 libMesh::err <<
"ERROR: This example requires libMesh to be\n"
104 <<
"compiled with SLEPc eigen solvers support!"
110 #ifdef LIBMESH_DEFAULT_SINGLE_PRECISION
112 libmesh_example_requires(
false,
"--disable-singleprecision");
115 #if defined(LIBMESH_USE_COMPLEX_NUMBERS) && SLEPC_VERSION_LESS_THAN(3,6,2)
118 libmesh_example_requires(
false,
"--disable-complex or use SLEPc>=3.6.2");
122 #ifndef LIBMESH_ENABLE_DIRICHLET
123 libmesh_example_requires(
false,
"--enable-dirichlet");
130 for (
int i=1; i<argc; i++)
137 libmesh_example_requires(2 <= LIBMESH_DIM,
"2D support");
140 GetPot command_line (argc, argv);
143 std::string mesh_name =
"";
144 if (command_line.search(1,
"-mesh_name"))
145 mesh_name = command_line.next(mesh_name);
149 unsigned int plotting_index = 0;
150 if (command_line.search(1,
"-plotting_index"))
151 plotting_index = command_line.next(plotting_index);
154 unsigned int n_evals = 0;
155 if (command_line.search(1,
"-n_evals"))
156 n_evals = command_line.next(n_evals);
159 std::ostringstream mesh_name_exodus;
160 mesh_name_exodus << mesh_name <<
"_mesh.e";
172 for (
auto side : elem->side_index_range())
173 if (elem->neighbor_ptr (side) ==
nullptr)
198 equation_systems.parameters.set<
unsigned int>(
"eigenpairs") = n_evals;
199 equation_systems.parameters.set<
unsigned int>(
"basis vectors") = n_evals*3;
202 equation_systems.parameters.set<
Real>(
"linear solver tolerance") =
pow(
TOLERANCE, 5./3.);
203 equation_systems.parameters.set<
unsigned int>
204 (
"linear solver maximum iterations") = 1000;
214 std::set<boundary_id_type> boundary_ids;
215 boundary_ids.insert(BOUNDARY_ID);
217 std::vector<unsigned int> variables;
218 variables.push_back(0);
222 #ifdef LIBMESH_ENABLE_DIRICHLET
233 equation_systems.init();
236 equation_systems.print_info();
241 eigen_system.
solve();
251 if (plotting_index > n_evals)
253 libMesh::out <<
"WARNING: Solver did not converge for the requested eigenvector!" << std::endl;
257 std::ostringstream eigenvalue_output_name;
258 eigenvalue_output_name << mesh_name <<
"_evals.txt";
259 std::ofstream evals_file(eigenvalue_output_name.str().c_str());
261 for (
unsigned int i=0; i<nconv; i++)
266 libmesh_assert_less (eval.second,
TOLERANCE);
267 evals_file << eval.first << std::endl;
270 if (i == plotting_index)
272 #ifdef LIBMESH_HAVE_EXODUS_API
274 std::ostringstream eigenvector_output_name;
275 eigenvector_output_name << mesh_name <<
"_evec.e";
277 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
283 #endif // LIBMESH_HAVE_SLEPC
References libMesh::DofMap::add_dirichlet_boundary(), libMesh::BoundaryInfo::add_side(), libMesh::EquationSystems::add_system(), libMesh::System::add_variable(), assemble_matrices(), libMesh::System::attach_assemble_function(), libMesh::EigenSystem::eigen_solver, libMesh::MeshBase::element_ptr_range(), libMesh::err, libMesh::MeshBase::get_boundary_info(), libMesh::System::get_dof_map(), libMesh::CondensedEigenSystem::get_eigenpair(), libMesh::EigenSystem::get_n_converged(), libMesh::GHEP, libMesh::TriangleWrapper::init(), libMesh::EquationSystems::init(), libMesh::CondensedEigenSystem::initialize_condensed_dofs(), libMesh::LOCAL_VARIABLE_ORDER, mesh, libMesh::out, libMesh::EquationSystems::parameters, std::pow(), libMesh::EquationSystems::print_info(), libMesh::MeshBase::print_info(), libMesh::MeshBase::read(), libMesh::Real, libMesh::SECOND, libMesh::Parameters::set(), libMesh::EigenSystem::set_eigenproblem_type(), libMesh::CondensedEigenSystem::solve(), libMesh::TARGET_REAL, 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.
void assemble_matrices(EquationSystems &es, const std::string &system_name)
const MeshBase & get_mesh() const
virtual void read(const std::string &name, void *mesh_data=nullptr, bool skip_renumber_nodes_and_elements=false, bool skip_find_neighbors=false)=0
Interfaces for reading/writing a mesh to/from a file.
const BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
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 void solve() override
Override to solve the condensed eigenproblem with the dofs in local_non_condensed_dofs_vector strippe...
void add_dirichlet_boundary(const DirichletBoundary &dirichlet_boundary)
Adds a copy of the specified Dirichlet boundary to the system.
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
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.
virtual std::pair< Real, Real > get_eigenpair(dof_id_type i) override
Override get_eigenpair() to retrieve the eigenpair for the condensed eigensolve.
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.
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
virtual SimpleRange< element_iterator > element_ptr_range()=0
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.
This class extends EigenSystem to allow a simple way of solving (standard or generalized) eigenvalue ...
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.
ConstFunction that simply returns 0.
void initialize_condensed_dofs(const std::set< dof_id_type > &global_condensed_dofs_set=std::set< dof_id_type >())
Loop over the dofs on each processor to initialize the list of non-condensed dofs.
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.
This class allows one to associate Dirichlet boundary values with a given set of mesh boundary ids an...
const DofMap & get_dof_map() const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void add_side(const dof_id_type elem, const unsigned short int side, const boundary_id_type id)
Add side side of element number elem with boundary id id to the boundary information data structure.