libMesh
Functions
eigenproblems_ex3.C File Reference

Go to the source code of this file.

Functions

void assemble_matrices (EquationSystems &es, const std::string &system_name)
 
void get_dirichlet_dofs (EquationSystems &es, const std::string &system_name, std::set< unsigned int > &global_dirichlet_dofs_set)
 
int main (int argc, char **argv)
 
void assemble_matrices (EquationSystems &es, const std::string &libmesh_dbg_var(system_name))
 

Function Documentation

◆ 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.

293 {
294 
295  // It is a good idea to make sure we are assembling
296  // the proper system.
297  libmesh_assert_equal_to (system_name, "Eigensystem");
298 
299 #ifdef LIBMESH_HAVE_SLEPC
300 
301  // Get a constant reference to the mesh object.
302  const MeshBase & mesh = es.get_mesh();
303 
304  // The dimension that we are running.
305  const unsigned int dim = mesh.mesh_dimension();
306 
307  // Get a reference to our system.
308  EigenSystem & eigen_system = es.get_system<EigenSystem> ("Eigensystem");
309 
310  // Get a constant reference to the Finite Element type
311  // for the first (and only) variable in the system.
312  FEType fe_type = eigen_system.get_dof_map().variable_type(0);
313 
314  // A reference to the two system matrices
315  SparseMatrix<Number> & matrix_A = *eigen_system.matrix_A;
316  SparseMatrix<Number> & matrix_B = *eigen_system.matrix_B;
317 
318  // Build a Finite Element object of the specified type. Since the
319  // FEBase::build() member dynamically creates memory we will
320  // store the object as a std::unique_ptr<FEBase>. This can be thought
321  // of as a pointer that will clean up after itself.
322  std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
323 
324  // A Gauss quadrature rule for numerical integration.
325  // Use the default quadrature order.
326  QGauss qrule (dim, fe_type.default_quadrature_order());
327 
328  // Tell the finite element object to use our quadrature rule.
329  fe->attach_quadrature_rule (&qrule);
330 
331  // The element Jacobian * quadrature weight at each integration point.
332  const std::vector<Real> & JxW = fe->get_JxW();
333 
334  // The element shape functions evaluated at the quadrature points.
335  const std::vector<std::vector<Real>> & phi = fe->get_phi();
336 
337  // The element shape function gradients evaluated at the quadrature
338  // points.
339  const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
340 
341  // A reference to the DofMap object for this system. The DofMap
342  // object handles the index translation from node and element numbers
343  // to degree of freedom numbers.
344  const DofMap & dof_map = eigen_system.get_dof_map();
345 
346  // The element mass and stiffness matrices.
349 
350  // This vector will hold the degree of freedom indices for
351  // the element. These define where in the global system
352  // the element degrees of freedom get mapped.
353  std::vector<dof_id_type> dof_indices;
354 
355 
356  // Now we will loop over all the elements in the mesh that
357  // live on the local processor. We will compute the element
358  // matrix and right-hand-side contribution. In case users
359  // later modify this program to include refinement, we will
360  // be safe and will only consider the active elements;
361  // hence we use a variant of the active_elem_iterator.
362  for (const auto & elem : mesh.active_local_element_ptr_range())
363  {
364  // Get the degree of freedom indices for the
365  // current element. These define where in the global
366  // matrix and right-hand-side this element will
367  // contribute to.
368  dof_map.dof_indices (elem, dof_indices);
369 
370  // Compute the element-specific data for the current
371  // element. This involves computing the location of the
372  // quadrature points (q_point) and the shape functions
373  // (phi, dphi) for the current element.
374  fe->reinit (elem);
375 
376  // Zero the element matrices before
377  // summing them. We use the resize member here because
378  // the number of degrees of freedom might have changed from
379  // the last element. Note that this will be the case if the
380  // element type is different (i.e. the last element was a
381  // triangle, now we are on a quadrilateral).
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);
386 
387  // Now loop over the quadrature points. This handles
388  // the numeric integration.
389  //
390  // We will build the element matrix. This involves
391  // a double loop to integrate the test functions (i) against
392  // the trial functions (j).
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++)
396  {
397  Me(i,j) += JxW[qp]*phi[i][qp]*phi[j][qp];
398  Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
399  }
400 
401  // The calls to constrain_element_matrix below have no effect in
402  // the current example. However, if users modify this example to
403  // include hanging nodes due to mesh refinement, for example,
404  // then it is essential to call constrain_element_matrix.
405  // As a result we include constrain_element_matrix here to
406  // ensure this example is ready to be used with hanging nodes.
407  // (Note that constrained rows/cols will be eliminated from
408  // the eigenproblem by the CondensedEigenSystem.)
409  dof_map.constrain_element_matrix(Ke, dof_indices, false);
410  dof_map.constrain_element_matrix(Me, dof_indices, false);
411 
412  // Finally, simply add the element contribution to the
413  // overall matrices A and B.
414  matrix_A.add_matrix (Ke, dof_indices);
415  matrix_B.add_matrix (Me, dof_indices);
416  } // end of element loop
417 
418 
419 #else
420  // Avoid compiler warnings
421  libmesh_ignore(es);
422 #endif // LIBMESH_HAVE_SLEPC
423 }

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 
)

Referenced by main().

◆ 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.

91 {
92  // Initialize libMesh and the dependent libraries.
93  LibMeshInit init (argc, argv);
94 
95  // This example uses an ExodusII input file
96 #ifndef LIBMESH_HAVE_EXODUS_API
97  libmesh_example_requires(false, "--enable-exodus");
98 #endif
99 
100  // This example is designed for the SLEPc eigen solver interface.
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!"
105  << std::endl;
106 
107  return 0;
108 #else
109 
110 #ifdef LIBMESH_DEFAULT_SINGLE_PRECISION
111  // SLEPc currently gives us a nasty crash with Real==float
112  libmesh_example_requires(false, "--disable-singleprecision");
113 #endif
114 
115 #if defined(LIBMESH_USE_COMPLEX_NUMBERS) && SLEPC_VERSION_LESS_THAN(3,6,2)
116  // SLEPc used to give us an "inner product not well defined" with
117  // Number==complex; but this problem seems to be solved in newer versions.
118  libmesh_example_requires(false, "--disable-complex or use SLEPc>=3.6.2");
119 #endif
120 
121  // We use Dirichlet boundary conditions here
122 #ifndef LIBMESH_ENABLE_DIRICHLET
123  libmesh_example_requires(false, "--enable-dirichlet");
124 #endif
125 
126  // Tell the user what we are doing.
127  {
128  libMesh::out << "Running " << argv[0];
129 
130  for (int i=1; i<argc; i++)
131  libMesh::out << " " << argv[i];
132 
133  libMesh::out << std::endl << std::endl;
134  }
135 
136  // Skip this 2D example if libMesh was compiled as 1D-only.
137  libmesh_example_requires(2 <= LIBMESH_DIM, "2D support");
138 
139  // Use GetPot to parse the command line arguments
140  GetPot command_line (argc, argv);
141 
142  // Read the mesh name from the command line
143  std::string mesh_name = "";
144  if (command_line.search(1, "-mesh_name"))
145  mesh_name = command_line.next(mesh_name);
146 
147  // Also, read in the index of the eigenvector that we should plot
148  // (zero-based indexing, as usual!)
149  unsigned int plotting_index = 0;
150  if (command_line.search(1, "-plotting_index"))
151  plotting_index = command_line.next(plotting_index);
152 
153  // Finally, read in the number of eigenpairs we want to compute!
154  unsigned int n_evals = 0;
155  if (command_line.search(1, "-n_evals"))
156  n_evals = command_line.next(n_evals);
157 
158  // Append the .e to mesh_name
159  std::ostringstream mesh_name_exodus;
160  mesh_name_exodus << mesh_name << "_mesh.e";
161 
162  // Create a mesh, with dimension to be overridden by the file, on
163  // the default MPI communicator.
164  Mesh mesh(init.comm());
165 
166  mesh.read(mesh_name_exodus.str());
167 
168  // Add boundary IDs to this mesh so that we can use DirichletBoundary
169  // Each processor should know about each boundary condition it can
170  // see, so we loop over all elements, not just local elements.
171  for (const auto & elem : mesh.element_ptr_range())
172  for (auto side : elem->side_index_range())
173  if (elem->neighbor_ptr (side) == nullptr)
174  mesh.get_boundary_info().add_side(elem, side, BOUNDARY_ID);
175 
176  // Print information about the mesh to the screen.
177  mesh.print_info();
178 
179  // Create an equation systems object.
180  EquationSystems equation_systems (mesh);
181 
182  // Create a CondensedEigenSystem named "Eigensystem" and (for convenience)
183  // use a reference to the system we create.
184  CondensedEigenSystem & eigen_system =
185  equation_systems.add_system<CondensedEigenSystem> ("Eigensystem");
186 
187  // Declare the system variables.
188  // Adds the variable "p" to "Eigensystem". "p"
189  // will be approximated using second-order approximation.
190  eigen_system.add_variable("p", SECOND);
191 
192  // Give the system a pointer to the matrix assembly
193  // function defined below.
195 
196  // Set the number of requested eigenpairs n_evals and the number
197  // of basis vectors used in the solution algorithm.
198  equation_systems.parameters.set<unsigned int>("eigenpairs") = n_evals;
199  equation_systems.parameters.set<unsigned int>("basis vectors") = n_evals*3;
200 
201  // Set the solver tolerance and the maximum number of iterations.
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;
205 
206  // Set the type of the problem, here we deal with
207  // a generalized Hermitian problem.
208  eigen_system.set_eigenproblem_type(GHEP);
209 
210  // Set the target eigenvalue
211  eigen_system.eigen_solver->set_position_of_spectrum(0., TARGET_REAL);
212 
213  {
214  std::set<boundary_id_type> boundary_ids;
215  boundary_ids.insert(BOUNDARY_ID);
216 
217  std::vector<unsigned int> variables;
218  variables.push_back(0);
219 
220  ZeroFunction<> zf;
221 
222 #ifdef LIBMESH_ENABLE_DIRICHLET
223  // Most DirichletBoundary users will want to supply a "locally
224  // indexed" functor
225  DirichletBoundary dirichlet_bc(boundary_ids, variables, zf,
227 
228  eigen_system.get_dof_map().add_dirichlet_boundary(dirichlet_bc);
229 #endif
230  }
231 
232  // Initialize the data structures for the equation system.
233  equation_systems.init();
234 
235  // Prints information about the system to the screen.
236  equation_systems.print_info();
237 
238  eigen_system.initialize_condensed_dofs();
239 
240  // Solve the system "Eigensystem".
241  eigen_system.solve();
242 
243  // Get the number of converged eigen pairs.
244  unsigned int nconv = eigen_system.get_n_converged();
245 
246  libMesh::out << "Number of converged eigenpairs: "
247  << nconv
248  << "\n"
249  << std::endl;
250 
251  if (plotting_index > n_evals)
252  {
253  libMesh::out << "WARNING: Solver did not converge for the requested eigenvector!" << std::endl;
254  }
255 
256  // write out all of the computed eigenvalues and plot the specified eigenvector
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());
260 
261  for (unsigned int i=0; i<nconv; i++)
262  {
263  std::pair<Real,Real> eval = eigen_system.get_eigenpair(i);
264 
265  // The eigenvalues should be real!
266  libmesh_assert_less (eval.second, TOLERANCE);
267  evals_file << eval.first << std::endl;
268 
269  // plot the specified eigenvector
270  if (i == plotting_index)
271  {
272 #ifdef LIBMESH_HAVE_EXODUS_API
273  // Write the eigen vector to file.
274  std::ostringstream eigenvector_output_name;
275  eigenvector_output_name << mesh_name << "_evec.e";
276  ExodusII_IO (mesh).write_equation_systems (eigenvector_output_name.str(), equation_systems);
277 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
278  }
279  }
280 
281  evals_file.close();
282 
283 #endif // LIBMESH_HAVE_SLEPC
284 
285  // All done.
286  return 0;
287 }

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().

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
assemble_matrices
void assemble_matrices(EquationSystems &es, const std::string &system_name)
libMesh::EquationSystems::get_mesh
const MeshBase & get_mesh() const
Definition: equation_systems.h:637
libMesh::MeshBase::read
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.
libMesh::MeshBase::get_boundary_info
const BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
Definition: mesh_base.h:132
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::TARGET_REAL
Definition: enum_eigen_solver_type.h:80
libMesh::CondensedEigenSystem::solve
virtual void solve() override
Override to solve the condensed eigenproblem with the dofs in local_non_condensed_dofs_vector strippe...
Definition: condensed_eigen_system.C:91
libMesh::DofMap::add_dirichlet_boundary
void add_dirichlet_boundary(const DirichletBoundary &dirichlet_boundary)
Adds a copy of the specified Dirichlet boundary to the system.
Definition: dof_map_constraints.C:4390
libMesh::MeshBase::active_local_element_ptr_range
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
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::CondensedEigenSystem::get_eigenpair
virtual std::pair< Real, Real > get_eigenpair(dof_id_type i) override
Override get_eigenpair() to retrieve the eigenpair for the condensed eigensolve.
Definition: condensed_eigen_system.C:177
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::SECOND
Definition: enum_order.h:43
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::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
libMesh::TriangleWrapper::init
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
libMesh::MeshBase::element_ptr_range
virtual SimpleRange< element_iterator > element_ptr_range()=0
libMesh::LOCAL_VARIABLE_ORDER
Definition: dirichlet_boundaries.h:62
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::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::CondensedEigenSystem
This class extends EigenSystem to allow a simple way of solving (standard or generalized) eigenvalue ...
Definition: condensed_eigen_system.h:46
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::ZeroFunction
ConstFunction that simply returns 0.
Definition: zero_function.h:36
libMesh::CondensedEigenSystem::initialize_condensed_dofs
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.
Definition: condensed_eigen_system.C:47
libMesh::EquationSystems
This is the EquationSystems class.
Definition: equation_systems.h:74
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::DofMap::constrain_element_matrix
void constrain_element_matrix(DenseMatrix< Number > &matrix, std::vector< dof_id_type > &elem_dofs, bool asymmetric_constraint_rows=true) const
Constrains the element matrix.
Definition: dof_map.h:2021
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::DirichletBoundary
This class allows one to associate Dirichlet boundary values with a given set of mesh boundary ids an...
Definition: dirichlet_boundaries.h:88
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
libMesh::BoundaryInfo::add_side
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.
Definition: boundary_info.C:886