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 &  system_name 
)

Referenced by main().

◆ assemble_matrices() [2/2]

void assemble_matrices ( EquationSystems es,
const std::string &  libmesh_dbg_varsystem_name 
)

Definition at line 282 of file eigenproblems_ex3.C.

References 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::EigenSystem::get_matrix_A(), libMesh::EigenSystem::get_matrix_B(), libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), libMesh::libmesh_ignore(), mesh, libMesh::MeshBase::mesh_dimension(), libMesh::DenseMatrix< T >::resize(), and libMesh::DofMap::variable_type().

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

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

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::command_line_next(), libMesh::err, libMesh::MeshBase::get_boundary_info(), libMesh::System::get_dof_map(), libMesh::EigenSystem::get_eigen_solver(), 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, libMesh::Utility::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::ExodusII_IO::write_equation_systems().

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  const std::string mesh_name =
144  libMesh::command_line_next("-mesh_name", std::string());
145 
146  // Also, read in the index of the eigenvector that we should plot
147  // (zero-based indexing, as usual!)
148  const unsigned int plotting_index =
149  libMesh::command_line_next("-plotting_index", 0u);
150 
151  // Finally, read in the number of eigenpairs we want to compute!
152  const unsigned int n_evals =
153  libMesh::command_line_next("-n_evals", 0u);
154 
155  // Append the .e to mesh_name
156  std::ostringstream mesh_name_exodus;
157  mesh_name_exodus << mesh_name << "_mesh.e";
158 
159  // Create a mesh, with dimension to be overridden by the file, on
160  // the default MPI communicator.
161  Mesh mesh(init.comm());
162 
163  mesh.read(mesh_name_exodus.str());
164 
165  // Add boundary IDs to this mesh so that we can use DirichletBoundary
166  // Each processor should know about each boundary condition it can
167  // see, so we loop over all elements, not just local elements.
168  for (const auto & elem : mesh.element_ptr_range())
169  for (auto side : elem->side_index_range())
170  if (elem->neighbor_ptr (side) == nullptr)
171  mesh.get_boundary_info().add_side(elem, side, BOUNDARY_ID);
172 
173  // Print information about the mesh to the screen.
174  mesh.print_info();
175 
176  // Create an equation systems object.
177  EquationSystems equation_systems (mesh);
178 
179  // Create a CondensedEigenSystem named "Eigensystem" and (for convenience)
180  // use a reference to the system we create.
181  CondensedEigenSystem & eigen_system =
182  equation_systems.add_system<CondensedEigenSystem> ("Eigensystem");
183 
184  // Declare the system variables.
185  // Adds the variable "p" to "Eigensystem". "p"
186  // will be approximated using second-order approximation.
187  eigen_system.add_variable("p", SECOND);
188 
189  // Give the system a pointer to the matrix assembly
190  // function defined below.
192 
193  // Set the number of requested eigenpairs n_evals and the number
194  // of basis vectors used in the solution algorithm.
195  equation_systems.parameters.set<unsigned int>("eigenpairs") = n_evals;
196  equation_systems.parameters.set<unsigned int>("basis vectors") = n_evals*3;
197 
198  // Set the solver tolerance and the maximum number of iterations.
199  equation_systems.parameters.set<Real>("linear solver tolerance") = pow(TOLERANCE, 5./3.);
200  equation_systems.parameters.set<unsigned int>
201  ("linear solver maximum iterations") = 1000;
202 
203  // Set the type of the problem, here we deal with
204  // a generalized Hermitian problem.
205  eigen_system.set_eigenproblem_type(GHEP);
206 
207  // Set the target eigenvalue
208  eigen_system.get_eigen_solver().set_position_of_spectrum(0., TARGET_REAL);
209 
210  {
211  ZeroFunction<> zf;
212 
213 #ifdef LIBMESH_ENABLE_DIRICHLET
214  // Most DirichletBoundary users will want to supply a "locally
215  // indexed" functor
216  DirichletBoundary dirichlet_bc({BOUNDARY_ID}, {0}, zf,
218 
219  eigen_system.get_dof_map().add_dirichlet_boundary(dirichlet_bc);
220 #endif
221  }
222 
223  // Initialize the data structures for the equation system.
224  equation_systems.init();
225 
226  // Prints information about the system to the screen.
227  equation_systems.print_info();
228 
229  eigen_system.initialize_condensed_dofs();
230 
231  // Solve the system "Eigensystem".
232  eigen_system.solve();
233 
234  // Get the number of converged eigen pairs.
235  unsigned int nconv = eigen_system.get_n_converged();
236 
237  libMesh::out << "Number of converged eigenpairs: "
238  << nconv
239  << "\n"
240  << std::endl;
241 
242  if (plotting_index > n_evals)
243  {
244  libMesh::out << "WARNING: Solver did not converge for the requested eigenvector!" << std::endl;
245  }
246 
247  // write out all of the computed eigenvalues and plot the specified eigenvector
248  std::ostringstream eigenvalue_output_name;
249  eigenvalue_output_name << mesh_name << "_evals.txt";
250  std::ofstream evals_file(eigenvalue_output_name.str().c_str());
251 
252  for (unsigned int i=0; i<nconv; i++)
253  {
254  std::pair<Real,Real> eval = eigen_system.get_eigenpair(i);
255 
256  // The eigenvalues should be real!
257  libmesh_assert_less (eval.second, TOLERANCE);
258  evals_file << eval.first << std::endl;
259 
260  // plot the specified eigenvector
261  if (i == plotting_index)
262  {
263 #ifdef LIBMESH_HAVE_EXODUS_API
264  // Write the eigen vector to file.
265  std::ostringstream eigenvector_output_name;
266  eigenvector_output_name << mesh_name << "_evec.e";
267  ExodusII_IO (mesh).write_equation_systems (eigenvector_output_name.str(), equation_systems);
268 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
269  }
270  }
271 
272  evals_file.close();
273 
274 #endif // LIBMESH_HAVE_SLEPC
275 
276  // All done.
277  return 0;
278 }
OStreamProxy err
T command_line_next(std::string name, T default_value)
Use GetPot&#39;s search()/next() functions to get following arguments from the command line...
Definition: libmesh.C:1078
This is the EquationSystems class.
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.
ConstFunction that simply returns 0.
Definition: zero_function.h:38
static constexpr Real TOLERANCE
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Definition: exodusII_io.h:52
MeshBase & mesh
This class allows one to associate Dirichlet boundary values with a given set of mesh boundary ids an...
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:90
const EigenSolver< Number > & get_eigen_solver() const
Definition: eigen_system.C:482
const BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
Definition: mesh_base.h:165
virtual void solve() override
Override to solve the condensed eigenproblem with the dofs in local_non_condensed_dofs_vector strippe...
T pow(const T &x)
Definition: utility.h:328
virtual void write_equation_systems(const std::string &fname, const EquationSystems &es, const std::set< std::string > *system_names=nullptr) override
Writes out the solution for no specific time or timestep.
Definition: exodusII_io.C:2033
void print_info(std::ostream &os=libMesh::out, const unsigned int verbosity=0, const bool global=true) const
Prints relevant information about the mesh.
Definition: mesh_base.C:1562
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
unsigned int add_variable(std::string_view 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:1357
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.
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:2161
This class extends EigenSystem to allow a simple way of solving (standard or generalized) eigenvalue ...
void assemble_matrices(EquationSystems &es, const std::string &system_name)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
OStreamProxy out
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...
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.
Definition: eigen_system.C:85
void add_dirichlet_boundary(const DirichletBoundary &dirichlet_boundary)
Adds a copy of the specified Dirichlet boundary to the system.
unsigned int get_n_converged() const
Definition: eigen_system.h:130
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50
const DofMap & get_dof_map() const
Definition: system.h:2374