libMesh::SlepcEigenSolver< T > Class Template Reference

EigenSolver implementation based on SLEPc. More...

#include <slepc_eigen_solver.h>

Inheritance diagram for libMesh::SlepcEigenSolver< T >:

Public Member Functions

 SlepcEigenSolver (const Parallel::Communicator &comm_in)
 
 ~SlepcEigenSolver ()
 
virtual void clear () override
 
virtual void init () override
 
virtual std::pair< unsigned int, unsigned int > solve_standard (SparseMatrix< T > &matrix_A, int nev, int ncv, const double tol, const unsigned int m_its) override
 
virtual std::pair< unsigned int, unsigned int > solve_standard (ShellMatrix< T > &shell_matrix, int nev, int ncv, const double tol, const unsigned int m_its) override
 
virtual std::pair< unsigned int, unsigned int > solve_generalized (SparseMatrix< T > &matrix_A, SparseMatrix< T > &matrix_B, int nev, int ncv, const double tol, const unsigned int m_its) override
 
virtual std::pair< unsigned int, unsigned int > solve_generalized (ShellMatrix< T > &matrix_A, SparseMatrix< T > &matrix_B, int nev, int ncv, const double tol, const unsigned int m_its) override
 
virtual std::pair< unsigned int, unsigned int > solve_generalized (SparseMatrix< T > &matrix_A, ShellMatrix< T > &matrix_B, int nev, int ncv, const double tol, const unsigned int m_its) override
 
virtual std::pair< unsigned int, unsigned int > solve_generalized (ShellMatrix< T > &matrix_A, ShellMatrix< T > &matrix_B, int nev, int ncv, const double tol, const unsigned int m_its) override
 
virtual std::pair< Real, Realget_eigenpair (dof_id_type i, NumericVector< T > &solution_in) override
 
virtual std::pair< Real, Realget_eigenvalue (dof_id_type i) override
 
Real get_relative_error (unsigned int i)
 
virtual void attach_deflation_space (NumericVector< T > &deflation_vector) override
 
virtual void set_initial_space (NumericVector< T > &initial_space_in) override
 
EPS eps ()
 
bool initialized () const
 
bool get_close_matrix_before_solve () const
 
void set_close_matrix_before_solve (bool val)
 
EigenSolverType eigen_solver_type () const
 
EigenProblemType eigen_problem_type () const
 
PositionOfSpectrum position_of_spectrum () const
 
void set_eigensolver_type (const EigenSolverType est)
 
void set_eigenproblem_type (EigenProblemType ept)
 
void set_position_of_spectrum (PositionOfSpectrum pos)
 
void set_position_of_spectrum (Real pos)
 
void set_position_of_spectrum (Real pos, PositionOfSpectrum target)
 
void set_solver_configuration (SolverConfiguration &solver_configuration)
 
const Parallel::Communicatorcomm () const
 
processor_id_type n_processors () const
 
processor_id_type processor_id () const
 

Static Public Member Functions

static std::unique_ptr< EigenSolver< T > > build (const Parallel::Communicator &comm_in, const SolverPackage solver_package=SLEPC_SOLVERS)
 
static std::string get_info ()
 
static void print_info (std::ostream &out=libMesh::out)
 
static unsigned int n_objects ()
 
static void enable_print_counter_info ()
 
static void disable_print_counter_info ()
 

Protected Types

typedef std::map< std::string, std::pair< unsigned int, unsigned int > > Counts
 

Protected Member Functions

void increment_constructor_count (const std::string &name)
 
void increment_destructor_count (const std::string &name)
 

Protected Attributes

EigenSolverType _eigen_solver_type
 
EigenProblemType _eigen_problem_type
 
PositionOfSpectrum _position_of_spectrum
 
bool _is_initialized
 
SolverConfiguration_solver_configuration
 
Real _target_val
 
bool _close_matrix_before_solve
 
const Parallel::Communicator_communicator
 

Static Protected Attributes

static Counts _counts
 
static Threads::atomic< unsigned int > _n_objects
 
static Threads::spin_mutex _mutex
 
static bool _enable_print_counter = true
 

Private Member Functions

std::pair< unsigned int, unsigned int > _solve_standard_helper (Mat mat, int nev, int ncv, const double tol, const unsigned int m_its)
 
std::pair< unsigned int, unsigned int > _solve_generalized_helper (Mat mat_A, Mat mat_B, int nev, int ncv, const double tol, const unsigned int m_its)
 
void set_slepc_solver_type ()
 
void set_slepc_problem_type ()
 
void set_slepc_position_of_spectrum ()
 

Static Private Member Functions

static PetscErrorCode _petsc_shell_matrix_mult (Mat mat, Vec arg, Vec dest)
 
static PetscErrorCode _petsc_shell_matrix_get_diagonal (Mat mat, Vec dest)
 

Private Attributes

EPS _eps
 

Detailed Description

template<typename T>
class libMesh::SlepcEigenSolver< T >

EigenSolver implementation based on SLEPc.

This class provides an interface to the SLEPc eigenvalue solver library from http://slepc.upv.es/.

Author
Steffen Peterson
Date
2005

Definition at line 50 of file slepc_eigen_solver.h.

Member Typedef Documentation

◆ Counts

typedef std::map<std::string, std::pair<unsigned int, unsigned int> > libMesh::ReferenceCounter::Counts
protectedinherited

Data structure to log the information. The log is identified by the class name.

Definition at line 117 of file reference_counter.h.

Constructor & Destructor Documentation

◆ SlepcEigenSolver()

template<typename T >
libMesh::SlepcEigenSolver< T >::SlepcEigenSolver ( const Parallel::Communicator comm_in)

Constructor. Initializes Petsc data structures

Definition at line 43 of file slepc_eigen_solver.C.

References libMesh::EigenSolver< T >::_eigen_problem_type, libMesh::EigenSolver< T >::_eigen_solver_type, libMesh::ARNOLDI, and libMesh::NHEP.

43  :
44  EigenSolver<T>(comm_in)
45 {
47  this->_eigen_problem_type = NHEP;
48 }
EigenSolverType _eigen_solver_type
Definition: eigen_solver.h:280
EigenProblemType _eigen_problem_type
Definition: eigen_solver.h:285

◆ ~SlepcEigenSolver()

template<typename T >
libMesh::SlepcEigenSolver< T >::~SlepcEigenSolver ( )

Destructor.

Definition at line 53 of file slepc_eigen_solver.C.

54 {
55  this->clear ();
56 }
virtual void clear() override

Member Function Documentation

◆ _petsc_shell_matrix_get_diagonal()

template<typename T >
PetscErrorCode libMesh::SlepcEigenSolver< T >::_petsc_shell_matrix_get_diagonal ( Mat  mat,
Vec  dest 
)
staticprivate

Internal function if shell matrix mode is used, this just calls the shell matrix's get_diagonal function. Required in order to use Jacobi preconditioning.

Definition at line 889 of file slepc_eigen_solver.C.

References libMesh::ParallelObject::comm(), libMesh::ShellMatrix< T >::get_diagonal(), and ierr.

890 {
891  // Get the matrix context.
892  PetscErrorCode ierr=0;
893  void * ctx;
894  ierr = MatShellGetContext(mat,&ctx);
895 
897  PetscObjectGetComm((PetscObject)mat,&comm);
898  CHKERRABORT(comm,ierr);
899 
900  // Get user shell matrix object.
901  const ShellMatrix<T> & shell_matrix = *static_cast<const ShellMatrix<T> *>(ctx);
902 
903  // Make \p NumericVector instances around the vector.
904  PetscVector<T> dest_global(dest, shell_matrix.comm());
905 
906  // Call the user function.
907  shell_matrix.get_diagonal(dest_global);
908 
909  return ierr;
910 }
MPI_Comm communicator
Definition: communicator.h:57
const Parallel::Communicator & comm() const
PetscErrorCode ierr

◆ _petsc_shell_matrix_mult()

template<typename T >
PetscErrorCode libMesh::SlepcEigenSolver< T >::_petsc_shell_matrix_mult ( Mat  mat,
Vec  arg,
Vec  dest 
)
staticprivate

Internal function if shell matrix mode is used, this just calls the shell matrix's matrix multiplication function. See PetscLinearSolver for a similar implementation.

Definition at line 864 of file slepc_eigen_solver.C.

References libMesh::ParallelObject::comm(), ierr, and libMesh::ShellMatrix< T >::vector_mult().

865 {
866  // Get the matrix context.
867  PetscErrorCode ierr=0;
868  void * ctx;
869  ierr = MatShellGetContext(mat,&ctx);
870 
872  PetscObjectGetComm((PetscObject)mat,&comm);
873  CHKERRABORT(comm,ierr);
874 
875  // Get user shell matrix object.
876  const ShellMatrix<T> & shell_matrix = *static_cast<const ShellMatrix<T> *>(ctx);
877 
878  // Make \p NumericVector instances around the vectors.
879  PetscVector<T> arg_global(arg, shell_matrix.comm());
880  PetscVector<T> dest_global(dest, shell_matrix.comm());
881 
882  // Call the user function.
883  shell_matrix.vector_mult(dest_global,arg_global);
884 
885  return ierr;
886 }
MPI_Comm communicator
Definition: communicator.h:57
const Parallel::Communicator & comm() const
PetscErrorCode ierr

◆ _solve_generalized_helper()

template<typename T >
std::pair< unsigned int, unsigned int > libMesh::SlepcEigenSolver< T >::_solve_generalized_helper ( Mat  mat_A,
Mat  mat_B,
int  nev,
int  ncv,
const double  tol,
const unsigned int  m_its 
)
private

Helper function that actually performs the generalized eigensolve.

Definition at line 459 of file slepc_eigen_solver.C.

References ierr.

465 {
466  LOG_SCOPE("solve_generalized()", "SlepcEigenSolver");
467 
468  PetscErrorCode ierr=0;
469 
470  // converged eigen pairs and number of iterations
471  PetscInt nconv=0;
472  PetscInt its=0;
473 
474 #ifdef DEBUG
475  // The relative error.
476  PetscReal error, re, im;
477 
478  // Pointer to vectors of the real parts, imaginary parts.
479  PetscScalar kr, ki;
480 #endif
481 
482  // Set operators.
483  ierr = EPSSetOperators (_eps, mat_A, mat_B);
484  LIBMESH_CHKERR(ierr);
485 
486  //set the problem type and the position of the spectrum
489 
490  // Set eigenvalues to be computed.
491 #if SLEPC_VERSION_LESS_THAN(3,0,0)
492  ierr = EPSSetDimensions (_eps, nev, ncv);
493 #else
494  ierr = EPSSetDimensions (_eps, nev, ncv, PETSC_DECIDE);
495 #endif
496  LIBMESH_CHKERR(ierr);
497 
498 
499  // Set the tolerance and maximum iterations.
500  ierr = EPSSetTolerances (_eps, tol, m_its);
501  LIBMESH_CHKERR(ierr);
502 
503  // Set runtime options, e.g.,
504  // -eps_type <type>, -eps_nev <nev>, -eps_ncv <ncv>
505  // Similar to PETSc, these options will override those specified
506  // above as long as EPSSetFromOptions() is called _after_ any
507  // other customization routines.
508  ierr = EPSSetFromOptions (_eps);
509  LIBMESH_CHKERR(ierr);
510 
511  // If the SolverConfiguration object is provided, use it to override
512  // solver options.
513  if (this->_solver_configuration)
514  {
516  }
517 
518  // Solve the eigenproblem.
519  ierr = EPSSolve (_eps);
520  LIBMESH_CHKERR(ierr);
521 
522  // Get the number of iterations.
523  ierr = EPSGetIterationNumber (_eps, &its);
524  LIBMESH_CHKERR(ierr);
525 
526  // Get number of converged eigenpairs.
527  ierr = EPSGetConverged(_eps,&nconv);
528  LIBMESH_CHKERR(ierr);
529 
530 
531 #ifdef DEBUG
532  // ierr = PetscPrintf(this->comm().get(),
533  // "\n Number of iterations: %d\n"
534  // " Number of converged eigenpairs: %d\n\n", its, nconv);
535 
536  // Display eigenvalues and relative errors.
537  ierr = PetscPrintf(this->comm().get(),
538  " k ||Ax-kx||/|kx|\n"
539  " ----------------- -----------------\n" );
540  LIBMESH_CHKERR(ierr);
541 
542  for (PetscInt i=0; i<nconv; i++ )
543  {
544  ierr = EPSGetEigenpair(_eps, i, &kr, &ki, PETSC_NULL, PETSC_NULL);
545  LIBMESH_CHKERR(ierr);
546 
547 #if SLEPC_VERSION_LESS_THAN(3,6,0)
548  ierr = EPSComputeRelativeError(_eps, i, &error);
549 #else
550  ierr = EPSComputeError(_eps, i, EPS_ERROR_RELATIVE, &error);
551 #endif
552  LIBMESH_CHKERR(ierr);
553 
554 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
555  re = PetscRealPart(kr);
556  im = PetscImaginaryPart(kr);
557 #else
558  re = kr;
559  im = ki;
560 #endif
561 
562  if (im != .0)
563  {
564  ierr = PetscPrintf(this->comm().get()," %9f%+9f i %12f\n", re, im, error);
565  LIBMESH_CHKERR(ierr);
566  }
567  else
568  {
569  ierr = PetscPrintf(this->comm().get()," %12f %12f\n", re, error);
570  LIBMESH_CHKERR(ierr);
571  }
572  }
573 
574  ierr = PetscPrintf(this->comm().get(),"\n" );
575  LIBMESH_CHKERR(ierr);
576 #endif // DEBUG
577 
578  // return the number of converged eigenpairs
579  // and the number of iterations
580  return std::make_pair(nconv, its);
581 }
SolverConfiguration * _solver_configuration
Definition: eigen_solver.h:301
const Parallel::Communicator & comm() const
virtual void configure_solver()=0
PetscErrorCode ierr

◆ _solve_standard_helper()

template<typename T >
std::pair< unsigned int, unsigned int > libMesh::SlepcEigenSolver< T >::_solve_standard_helper ( Mat  mat,
int  nev,
int  ncv,
const double  tol,
const unsigned int  m_its 
)
private

Helper function that actually performs the standard eigensolve.

Definition at line 163 of file slepc_eigen_solver.C.

References ierr.

168 {
169  LOG_SCOPE("solve_standard()", "SlepcEigenSolver");
170 
171  PetscErrorCode ierr=0;
172 
173  // converged eigen pairs and number of iterations
174  PetscInt nconv=0;
175  PetscInt its=0;
176 
177 #ifdef DEBUG
178  // The relative error.
179  PetscReal error, re, im;
180 
181  // Pointer to vectors of the real parts, imaginary parts.
182  PetscScalar kr, ki;
183 #endif
184 
185  // Set operators.
186  ierr = EPSSetOperators (_eps, mat, PETSC_NULL);
187  LIBMESH_CHKERR(ierr);
188 
189  //set the problem type and the position of the spectrum
192 
193  // Set eigenvalues to be computed.
194 #if SLEPC_VERSION_LESS_THAN(3,0,0)
195  ierr = EPSSetDimensions (_eps, nev, ncv);
196 #else
197  ierr = EPSSetDimensions (_eps, nev, ncv, PETSC_DECIDE);
198 #endif
199  LIBMESH_CHKERR(ierr);
200  // Set the tolerance and maximum iterations.
201  ierr = EPSSetTolerances (_eps, tol, m_its);
202  LIBMESH_CHKERR(ierr);
203 
204  // Set runtime options, e.g.,
205  // -eps_type <type>, -eps_nev <nev>, -eps_ncv <ncv>
206  // Similar to PETSc, these options will override those specified
207  // above as long as EPSSetFromOptions() is called _after_ any
208  // other customization routines.
209  ierr = EPSSetFromOptions (_eps);
210  LIBMESH_CHKERR(ierr);
211 
212  // If the SolverConfiguration object is provided, use it to override
213  // solver options.
214  if (this->_solver_configuration)
215  {
217  }
218 
219  // Solve the eigenproblem.
220  ierr = EPSSolve (_eps);
221  LIBMESH_CHKERR(ierr);
222 
223  // Get the number of iterations.
224  ierr = EPSGetIterationNumber (_eps, &its);
225  LIBMESH_CHKERR(ierr);
226 
227  // Get number of converged eigenpairs.
228  ierr = EPSGetConverged(_eps,&nconv);
229  LIBMESH_CHKERR(ierr);
230 
231 
232 #ifdef DEBUG
233  // ierr = PetscPrintf(this->comm().get(),
234  // "\n Number of iterations: %d\n"
235  // " Number of converged eigenpairs: %d\n\n", its, nconv);
236 
237  // Display eigenvalues and relative errors.
238  ierr = PetscPrintf(this->comm().get(),
239  " k ||Ax-kx||/|kx|\n"
240  " ----------------- -----------------\n" );
241  LIBMESH_CHKERR(ierr);
242 
243  for (PetscInt i=0; i<nconv; i++ )
244  {
245  ierr = EPSGetEigenpair(_eps, i, &kr, &ki, PETSC_NULL, PETSC_NULL);
246  LIBMESH_CHKERR(ierr);
247 
248 #if SLEPC_VERSION_LESS_THAN(3,6,0)
249  ierr = EPSComputeRelativeError(_eps, i, &error);
250 #else
251  ierr = EPSComputeError(_eps, i, EPS_ERROR_RELATIVE, &error);
252 #endif
253  LIBMESH_CHKERR(ierr);
254 
255 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
256  re = PetscRealPart(kr);
257  im = PetscImaginaryPart(kr);
258 #else
259  re = kr;
260  im = ki;
261 #endif
262 
263  if (im != .0)
264  {
265  ierr = PetscPrintf(this->comm().get()," %9f%+9f i %12f\n", re, im, error);
266  LIBMESH_CHKERR(ierr);
267  }
268  else
269  {
270  ierr = PetscPrintf(this->comm().get()," %12f %12f\n", re, error);
271  LIBMESH_CHKERR(ierr);
272  }
273  }
274 
275  ierr = PetscPrintf(this->comm().get(),"\n" );
276  LIBMESH_CHKERR(ierr);
277 #endif // DEBUG
278 
279  // return the number of converged eigenpairs
280  // and the number of iterations
281  return std::make_pair(nconv, its);
282 }
SolverConfiguration * _solver_configuration
Definition: eigen_solver.h:301
const Parallel::Communicator & comm() const
virtual void configure_solver()=0
PetscErrorCode ierr

◆ attach_deflation_space()

template<typename T >
void libMesh::SlepcEigenSolver< T >::attach_deflation_space ( NumericVector< T > &  deflation_vector)
overridevirtual

Attach a deflation space defined by a single vector.

Implements libMesh::EigenSolver< T >.

Definition at line 814 of file slepc_eigen_solver.C.

References ierr, libMesh::TriangleWrapper::init(), and libMesh::PetscVector< T >::vec().

815 {
816  this->init();
817 
818  PetscErrorCode ierr = 0;
819 
820  // Make sure the input vector is actually a PetscVector
821  PetscVector<T> * deflation_vector_petsc_vec =
822  dynamic_cast<PetscVector<T> *>(&deflation_vector_in);
823 
824  if (!deflation_vector_petsc_vec)
825  libmesh_error_msg("Error attaching deflation space: input vector must be a PetscVector.");
826 
827  // Get a handle for the underlying Vec.
828  Vec deflation_vector = deflation_vector_petsc_vec->vec();
829 
830 #if SLEPC_VERSION_LESS_THAN(3,1,0)
831  ierr = EPSAttachDeflationSpace(_eps, 1, &deflation_vector, PETSC_FALSE);
832 #else
833  ierr = EPSSetDeflationSpace(_eps, 1, &deflation_vector);
834 #endif
835  LIBMESH_CHKERR(ierr);
836 }
PetscErrorCode ierr
virtual void init() override

◆ build()

template<typename T >
std::unique_ptr< EigenSolver< T > > libMesh::EigenSolver< T >::build ( const Parallel::Communicator comm_in,
const SolverPackage  solver_package = SLEPC_SOLVERS 
)
staticinherited

Builds an EigenSolver using the linear solver package specified by solver_package

Definition at line 59 of file eigen_solver.C.

References libMesh::SLEPC_SOLVERS.

61 {
62  // Build the appropriate solver
63  switch (solver_package)
64  {
65 
66 #ifdef LIBMESH_HAVE_SLEPC
67  case SLEPC_SOLVERS:
68  return libmesh_make_unique<SlepcEigenSolver<T>>(comm);
69 #endif
70 
71  default:
72  libmesh_error_msg("ERROR: Unrecognized eigen solver package: " << solver_package);
73  }
74 
75  return std::unique_ptr<EigenSolver<T>>();
76 }
const Parallel::Communicator & comm() const

◆ clear()

template<typename T >
void libMesh::SlepcEigenSolver< T >::clear ( )
overridevirtual

Release all memory and clear data structures.

Reimplemented from libMesh::EigenSolver< T >.

Definition at line 61 of file slepc_eigen_solver.C.

References libMesh::libMeshPrivateData::_is_initialized, ierr, libMesh::initialized(), and libMesh::KRYLOVSCHUR.

62 {
63  if (this->initialized())
64  {
65  this->_is_initialized = false;
66 
67  PetscErrorCode ierr=0;
68 
69  ierr = LibMeshEPSDestroy(&_eps);
70  LIBMESH_CHKERR(ierr);
71 
72  // SLEPc default eigenproblem solver
74  }
75 }
bool initialized() const
Definition: eigen_solver.h:94
PetscErrorCode ierr
EigenSolverType _eigen_solver_type
Definition: eigen_solver.h:280

◆ comm()

const Parallel::Communicator& libMesh::ParallelObject::comm ( ) const
inlineinherited
Returns
A reference to the Parallel::Communicator object used by this mesh.

Definition at line 89 of file parallel_object.h.

References libMesh::ParallelObject::_communicator.

Referenced by libMesh::__libmesh_petsc_diff_solver_jacobian(), libMesh::__libmesh_petsc_diff_solver_monitor(), libMesh::__libmesh_petsc_diff_solver_residual(), libMesh::__libmesh_tao_equality_constraints(), libMesh::__libmesh_tao_equality_constraints_jacobian(), libMesh::__libmesh_tao_gradient(), libMesh::__libmesh_tao_hessian(), libMesh::__libmesh_tao_inequality_constraints(), libMesh::__libmesh_tao_inequality_constraints_jacobian(), libMesh::__libmesh_tao_objective(), libMesh::MeshRefinement::_coarsen_elements(), libMesh::ExactSolution::_compute_error(), libMesh::UniformRefinementEstimator::_estimate_error(), libMesh::BoundaryInfo::_find_id_maps(), libMesh::SlepcEigenSolver< T >::_petsc_shell_matrix_get_diagonal(), libMesh::PetscLinearSolver< T >::_petsc_shell_matrix_get_diagonal(), libMesh::SlepcEigenSolver< T >::_petsc_shell_matrix_mult(), libMesh::PetscLinearSolver< T >::_petsc_shell_matrix_mult(), libMesh::PetscLinearSolver< T >::_petsc_shell_matrix_mult_add(), libMesh::EquationSystems::_read_impl(), libMesh::MeshRefinement::_refine_elements(), libMesh::MeshRefinement::_smooth_flags(), libMesh::PetscDMWrapper::add_dofs_helper(), libMesh::PetscDMWrapper::add_dofs_to_section(), libMesh::ImplicitSystem::add_matrix(), libMesh::System::add_vector(), libMesh::UnstructuredMesh::all_second_order(), libMesh::MeshTools::Modification::all_tri(), libMesh::LaplaceMeshSmoother::allgather_graph(), libMesh::FEMSystem::assemble_qoi(), libMesh::MeshCommunication::assign_global_indices(), libMesh::DofMap::attach_matrix(), libMesh::MeshTools::Generation::build_extrusion(), libMesh::BoundaryInfo::build_node_list_from_side_list(), libMesh::EquationSystems::build_parallel_elemental_solution_vector(), libMesh::EquationSystems::build_parallel_solution_vector(), libMesh::PetscDMWrapper::build_section(), libMesh::PetscDMWrapper::build_sf(), libMesh::MeshBase::cache_elem_dims(), libMesh::System::calculate_norm(), libMesh::DofMap::check_dirichlet_bcid_consistency(), libMesh::PetscDMWrapper::check_section_n_dofs(), libMesh::Nemesis_IO_Helper::compute_num_global_elem_blocks(), libMesh::Nemesis_IO_Helper::compute_num_global_nodesets(), libMesh::Nemesis_IO_Helper::compute_num_global_sidesets(), libMesh::Problem_Interface::computeF(), libMesh::Problem_Interface::computeJacobian(), libMesh::Problem_Interface::computePreconditioner(), libMesh::ExodusII_IO::copy_elemental_solution(), libMesh::MeshTools::correct_node_proc_ids(), libMesh::MeshTools::create_bounding_box(), libMesh::MeshTools::create_nodal_bounding_box(), libMesh::MeshRefinement::create_parent_error_vector(), libMesh::MeshTools::create_processor_bounding_box(), libMesh::MeshTools::create_subdomain_bounding_box(), libMesh::MeshCommunication::delete_remote_elements(), libMesh::DofMap::distribute_dofs(), DMlibMeshFunction(), DMlibMeshJacobian(), DMlibMeshSetSystem_libMesh(), DMVariableBounds_libMesh(), libMesh::MeshRefinement::eliminate_unrefined_patches(), libMesh::EpetraVector< T >::EpetraVector(), libMesh::WeightedPatchRecoveryErrorEstimator::estimate_error(), libMesh::PatchRecoveryErrorEstimator::estimate_error(), libMesh::JumpErrorEstimator::estimate_error(), libMesh::AdjointRefinementEstimator::estimate_error(), libMesh::ExactErrorEstimator::estimate_error(), libMesh::MeshRefinement::flag_elements_by_elem_fraction(), libMesh::MeshRefinement::flag_elements_by_error_fraction(), libMesh::MeshRefinement::flag_elements_by_nelem_target(), libMesh::CondensedEigenSystem::get_eigenpair(), libMesh::DofMap::get_info(), libMesh::ImplicitSystem::get_linear_solver(), libMesh::LocationMap< T >::init(), libMesh::TimeSolver::init(), libMesh::SystemSubsetBySubdomain::init(), libMesh::PetscDMWrapper::init_and_attach_petscdm(), libMesh::EigenSystem::init_data(), libMesh::EigenSystem::init_matrices(), libMesh::OptimizationSystem::initialize_equality_constraints_storage(), libMesh::OptimizationSystem::initialize_inequality_constraints_storage(), libMesh::MeshTools::libmesh_assert_consistent_distributed(), libMesh::MeshTools::libmesh_assert_consistent_distributed_nodes(), libMesh::MeshTools::libmesh_assert_contiguous_dof_ids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_new_node_procids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Elem >(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_topology_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_valid_boundary_ids(), libMesh::MeshTools::libmesh_assert_valid_dof_ids(), libMesh::MeshTools::libmesh_assert_valid_neighbors(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_flags(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_object_ids(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_p_levels(), libMesh::MeshTools::libmesh_assert_valid_refinement_flags(), libMesh::MeshTools::libmesh_assert_valid_unique_ids(), libMesh::libmesh_petsc_snes_fd_residual(), libMesh::libmesh_petsc_snes_jacobian(), libMesh::libmesh_petsc_snes_mffd_residual(), libMesh::libmesh_petsc_snes_postcheck(), libMesh::libmesh_petsc_snes_residual(), libMesh::libmesh_petsc_snes_residual_helper(), libMesh::MeshRefinement::limit_level_mismatch_at_edge(), libMesh::MeshRefinement::limit_level_mismatch_at_node(), libMesh::MeshRefinement::limit_overrefined_boundary(), libMesh::MeshRefinement::limit_underrefined_boundary(), libMesh::MeshRefinement::make_coarsening_compatible(), libMesh::MeshCommunication::make_elems_parallel_consistent(), libMesh::MeshRefinement::make_flags_parallel_consistent(), libMesh::MeshCommunication::make_new_node_proc_ids_parallel_consistent(), libMesh::MeshCommunication::make_new_nodes_parallel_consistent(), libMesh::MeshCommunication::make_node_ids_parallel_consistent(), libMesh::MeshCommunication::make_node_proc_ids_parallel_consistent(), libMesh::MeshCommunication::make_node_unique_ids_parallel_consistent(), libMesh::MeshCommunication::make_nodes_parallel_consistent(), libMesh::MeshCommunication::make_p_levels_parallel_consistent(), libMesh::MeshRefinement::make_refinement_compatible(), libMesh::FEMSystem::mesh_position_set(), libMesh::DistributedMesh::n_active_elem(), libMesh::MeshTools::n_active_levels(), libMesh::BoundaryInfo::n_boundary_conds(), libMesh::BoundaryInfo::n_edge_conds(), libMesh::CondensedEigenSystem::n_global_non_condensed_dofs(), libMesh::MeshTools::n_levels(), libMesh::BoundaryInfo::n_nodeset_conds(), libMesh::MeshTools::n_p_levels(), libMesh::BoundaryInfo::n_shellface_conds(), libMesh::DistributedMesh::parallel_max_elem_id(), libMesh::DistributedMesh::parallel_max_node_id(), libMesh::ReplicatedMesh::parallel_max_unique_id(), libMesh::DistributedMesh::parallel_max_unique_id(), libMesh::DistributedMesh::parallel_n_elem(), libMesh::DistributedMesh::parallel_n_nodes(), libMesh::SparsityPattern::Build::parallel_sync(), libMesh::MeshTools::paranoid_n_levels(), libMesh::petsc_auto_fieldsplit(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::MeshBase::prepare_for_use(), libMesh::Nemesis_IO::read(), libMesh::XdrIO::read(), libMesh::CheckpointIO::read_header(), libMesh::XdrIO::read_header(), libMesh::System::read_header(), libMesh::System::read_legacy_data(), libMesh::System::read_SCALAR_dofs(), libMesh::XdrIO::read_serialized_bc_names(), libMesh::XdrIO::read_serialized_bcs_helper(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::XdrIO::read_serialized_nodes(), libMesh::XdrIO::read_serialized_nodesets(), libMesh::XdrIO::read_serialized_subdomain_names(), libMesh::System::read_serialized_vector(), libMesh::MeshBase::recalculate_n_partitions(), libMesh::MeshRefinement::refine_and_coarsen_elements(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::CheckpointIO::select_split_config(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::PetscDMWrapper::set_point_range_in_section(), libMesh::PetscDiffSolver::setup_petsc_data(), libMesh::LaplaceMeshSmoother::smooth(), libMesh::split_mesh(), libMesh::MeshBase::subdomain_ids(), libMesh::BoundaryInfo::sync(), libMesh::MeshRefinement::test_level_one(), libMesh::MeshRefinement::test_unflagged(), libMesh::MeshTools::total_weight(), libMesh::MeshRefinement::uniformly_coarsen(), libMesh::NameBasedIO::write(), libMesh::XdrIO::write(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), and libMesh::XdrIO::write_serialized_nodesets().

90  { return _communicator; }
const Parallel::Communicator & _communicator

◆ disable_print_counter_info()

void libMesh::ReferenceCounter::disable_print_counter_info ( )
staticinherited

Definition at line 106 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

Referenced by libMesh::LibMeshInit::LibMeshInit().

107 {
108  _enable_print_counter = false;
109  return;
110 }

◆ eigen_problem_type()

template<typename T >
EigenProblemType libMesh::EigenSolver< T >::eigen_problem_type ( ) const
inlineinherited
Returns
The type of the eigen problem.

Definition at line 135 of file eigen_solver.h.

References libMesh::EigenSolver< T >::_eigen_problem_type.

135 { return _eigen_problem_type;}
EigenProblemType _eigen_problem_type
Definition: eigen_solver.h:285

◆ eigen_solver_type()

template<typename T >
EigenSolverType libMesh::EigenSolver< T >::eigen_solver_type ( ) const
inlineinherited
Returns
The type of eigensolver to use.

Definition at line 130 of file eigen_solver.h.

References libMesh::EigenSolver< T >::_eigen_solver_type.

130 { return _eigen_solver_type; }
EigenSolverType _eigen_solver_type
Definition: eigen_solver.h:280

◆ enable_print_counter_info()

void libMesh::ReferenceCounter::enable_print_counter_info ( )
staticinherited

Methods to enable/disable the reference counter output from print_info()

Definition at line 100 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

101 {
102  _enable_print_counter = true;
103  return;
104 }

◆ eps()

template<typename T >
EPS libMesh::SlepcEigenSolver< T >::eps ( )
inline
Returns
The raw SLEPc EPS pointer.

Definition at line 214 of file slepc_eigen_solver.h.

References libMesh::SlepcEigenSolver< T >::_eps, and libMesh::SlepcEigenSolver< T >::init().

214 { this->init(); return _eps; }
virtual void init() override

◆ get_close_matrix_before_solve()

template<typename T >
bool libMesh::EigenSolver< T >::get_close_matrix_before_solve ( ) const
inlineinherited
Returns
The value of the flag which controls whether libmesh closes the eigenproblem matrices before solving. true by default.

Definition at line 101 of file eigen_solver.h.

References libMesh::EigenSolver< T >::_close_matrix_before_solve.

102  {
103  libmesh_experimental();
105  }

◆ get_eigenpair()

template<typename T >
std::pair< Real, Real > libMesh::SlepcEigenSolver< T >::get_eigenpair ( dof_id_type  i,
NumericVector< T > &  solution_in 
)
overridevirtual
Returns
The real and imaginary part of the ith eigenvalue and copies the respective eigenvector to the solution vector.
Note
The eigenpair may be complex even for real-valued matrices.

Implements libMesh::EigenSolver< T >.

Definition at line 738 of file slepc_eigen_solver.C.

References libMesh::PetscVector< T >::close(), ierr, and libMesh::PetscVector< T >::vec().

740 {
741  PetscErrorCode ierr=0;
742 
743  PetscReal re, im;
744 
745  // Make sure the NumericVector passed in is really a PetscVector
746  PetscVector<T> * solution = dynamic_cast<PetscVector<T> *>(&solution_in);
747 
748  if (!solution)
749  libmesh_error_msg("Error getting eigenvector: input vector must be a PetscVector.");
750 
751  // real and imaginary part of the ith eigenvalue.
752  PetscScalar kr, ki;
753 
754  solution->close();
755 
756  ierr = EPSGetEigenpair(_eps, i, &kr, &ki, solution->vec(), PETSC_NULL);
757  LIBMESH_CHKERR(ierr);
758 
759 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
760  re = PetscRealPart(kr);
761  im = PetscImaginaryPart(kr);
762 #else
763  re = kr;
764  im = ki;
765 #endif
766 
767  return std::make_pair(re, im);
768 }
PetscErrorCode ierr

◆ get_eigenvalue()

template<typename T >
std::pair< Real, Real > libMesh::SlepcEigenSolver< T >::get_eigenvalue ( dof_id_type  i)
overridevirtual

Same as above, but does not copy the eigenvector.

Implements libMesh::EigenSolver< T >.

Definition at line 772 of file slepc_eigen_solver.C.

References ierr.

773 {
774  PetscErrorCode ierr=0;
775 
776  PetscReal re, im;
777 
778  // real and imaginary part of the ith eigenvalue.
779  PetscScalar kr, ki;
780 
781  ierr = EPSGetEigenvalue(_eps, i, &kr, &ki);
782  LIBMESH_CHKERR(ierr);
783 
784 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
785  re = PetscRealPart(kr);
786  im = PetscImaginaryPart(kr);
787 #else
788  re = kr;
789  im = ki;
790 #endif
791 
792  return std::make_pair(re, im);
793 }
PetscErrorCode ierr

◆ get_info()

std::string libMesh::ReferenceCounter::get_info ( )
staticinherited

Gets a string containing the reference information.

Definition at line 47 of file reference_counter.C.

References libMesh::ReferenceCounter::_counts, and libMesh::Quality::name().

Referenced by libMesh::ReferenceCounter::print_info().

48 {
49 #if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)
50 
51  std::ostringstream oss;
52 
53  oss << '\n'
54  << " ---------------------------------------------------------------------------- \n"
55  << "| Reference count information |\n"
56  << " ---------------------------------------------------------------------------- \n";
57 
58  for (const auto & pr : _counts)
59  {
60  const std::string name(pr.first);
61  const unsigned int creations = pr.second.first;
62  const unsigned int destructions = pr.second.second;
63 
64  oss << "| " << name << " reference count information:\n"
65  << "| Creations: " << creations << '\n'
66  << "| Destructions: " << destructions << '\n';
67  }
68 
69  oss << " ---------------------------------------------------------------------------- \n";
70 
71  return oss.str();
72 
73 #else
74 
75  return "";
76 
77 #endif
78 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:42

◆ get_relative_error()

template<typename T >
Real libMesh::SlepcEigenSolver< T >::get_relative_error ( unsigned int  i)
Returns
The relative error $ ||A x - \lambda x|| / |\lambda x| $ of the ith eigenpair (or the equivalent for a general eigenvalue problem).

Definition at line 797 of file slepc_eigen_solver.C.

References ierr.

798 {
799  PetscErrorCode ierr=0;
800  PetscReal error;
801 
802 #if SLEPC_VERSION_LESS_THAN(3,6,0)
803  ierr = EPSComputeRelativeError(_eps, i, &error);
804 #else
805  ierr = EPSComputeError(_eps, i, EPS_ERROR_RELATIVE, &error);
806 #endif
807  LIBMESH_CHKERR(ierr);
808 
809  return error;
810 }
PetscErrorCode ierr

◆ increment_constructor_count()

void libMesh::ReferenceCounter::increment_constructor_count ( const std::string &  name)
inlineprotectedinherited

Increments the construction counter. Should be called in the constructor of any derived class that will be reference counted.

Definition at line 181 of file reference_counter.h.

References libMesh::ReferenceCounter::_counts, libMesh::Quality::name(), and libMesh::Threads::spin_mtx.

Referenced by libMesh::ReferenceCountedObject< RBParametrized >::ReferenceCountedObject().

182 {
183  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
184  std::pair<unsigned int, unsigned int> & p = _counts[name];
185 
186  p.first++;
187 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:42
spin_mutex spin_mtx
Definition: threads.C:29

◆ increment_destructor_count()

void libMesh::ReferenceCounter::increment_destructor_count ( const std::string &  name)
inlineprotectedinherited

Increments the destruction counter. Should be called in the destructor of any derived class that will be reference counted.

Definition at line 194 of file reference_counter.h.

References libMesh::ReferenceCounter::_counts, libMesh::Quality::name(), and libMesh::Threads::spin_mtx.

Referenced by libMesh::ReferenceCountedObject< RBParametrized >::~ReferenceCountedObject().

195 {
196  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
197  std::pair<unsigned int, unsigned int> & p = _counts[name];
198 
199  p.second++;
200 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:42
spin_mutex spin_mtx
Definition: threads.C:29

◆ init()

template<typename T >
void libMesh::SlepcEigenSolver< T >::init ( )
overridevirtual

Initialize data structures if not done so already.

Implements libMesh::EigenSolver< T >.

Definition at line 80 of file slepc_eigen_solver.C.

References libMesh::libMeshPrivateData::_is_initialized, ierr, and libMesh::initialized().

Referenced by libMesh::SlepcEigenSolver< T >::eps().

81 {
82 
83  PetscErrorCode ierr=0;
84 
85  // Initialize the data structures if not done so already.
86  if (!this->initialized())
87  {
88  this->_is_initialized = true;
89 
90  // Create the eigenproblem solver context
91  ierr = EPSCreate (this->comm().get(), &_eps);
92  LIBMESH_CHKERR(ierr);
93 
94  // Set user-specified solver
96  }
97 }
const Parallel::Communicator & comm() const
bool initialized() const
Definition: eigen_solver.h:94
PetscErrorCode ierr

◆ initialized()

template<typename T >
bool libMesh::EigenSolver< T >::initialized ( ) const
inlineinherited
Returns
true if the data structures are initialized, false otherwise.

Definition at line 94 of file eigen_solver.h.

References libMesh::EigenSolver< T >::_is_initialized.

94 { return _is_initialized; }

◆ n_objects()

static unsigned int libMesh::ReferenceCounter::n_objects ( )
inlinestaticinherited

Prints the number of outstanding (created, but not yet destroyed) objects.

Definition at line 83 of file reference_counter.h.

References libMesh::ReferenceCounter::_n_objects.

84  { return _n_objects; }
static Threads::atomic< unsigned int > _n_objects

◆ n_processors()

processor_id_type libMesh::ParallelObject::n_processors ( ) const
inlineinherited
Returns
The number of processors in the group.

Definition at line 95 of file parallel_object.h.

References libMesh::ParallelObject::_communicator, and libMesh::Parallel::Communicator::size().

Referenced by libMesh::BoundaryInfo::_find_id_maps(), libMesh::PetscDMWrapper::add_dofs_to_section(), libMesh::DistributedMesh::add_elem(), libMesh::DistributedMesh::add_node(), libMesh::LaplaceMeshSmoother::allgather_graph(), libMesh::FEMSystem::assembly(), libMesh::AztecLinearSolver< T >::AztecLinearSolver(), libMesh::BoundaryInfo::build_node_list_from_side_list(), libMesh::EquationSystems::build_parallel_elemental_solution_vector(), libMesh::DistributedMesh::clear(), libMesh::Nemesis_IO_Helper::compute_border_node_ids(), libMesh::Nemesis_IO_Helper::construct_nemesis_filename(), libMesh::UnstructuredMesh::create_pid_mesh(), libMesh::MeshTools::create_processor_bounding_box(), libMesh::DofMap::distribute_dofs(), libMesh::DofMap::distribute_local_dofs_node_major(), libMesh::DofMap::distribute_local_dofs_var_major(), libMesh::EnsightIO::EnsightIO(), libMesh::MeshBase::get_info(), libMesh::SystemSubsetBySubdomain::init(), libMesh::PetscDMWrapper::init_and_attach_petscdm(), libMesh::Nemesis_IO_Helper::initialize(), libMesh::DistributedMesh::insert_elem(), libMesh::MeshTools::libmesh_assert_contiguous_dof_ids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_new_node_procids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Elem >(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_topology_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_valid_boundary_ids(), libMesh::MeshTools::libmesh_assert_valid_dof_ids(), libMesh::MeshTools::libmesh_assert_valid_neighbors(), libMesh::MeshTools::libmesh_assert_valid_refinement_flags(), libMesh::DofMap::local_variable_indices(), libMesh::MeshRefinement::make_coarsening_compatible(), libMesh::MeshBase::n_active_elem_on_proc(), libMesh::MeshBase::n_elem_on_proc(), libMesh::MeshBase::n_nodes_on_proc(), libMesh::MeshBase::partition(), libMesh::PetscLinearSolver< T >::PetscLinearSolver(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::NameBasedIO::read(), libMesh::Nemesis_IO::read(), libMesh::CheckpointIO::read(), libMesh::CheckpointIO::read_connectivity(), libMesh::XdrIO::read_header(), libMesh::CheckpointIO::read_nodes(), libMesh::System::read_parallel_data(), libMesh::System::read_SCALAR_dofs(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::System::read_serialized_vector(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::PetscDMWrapper::set_point_range_in_section(), libMesh::MeshRefinement::uniformly_coarsen(), libMesh::DistributedMesh::update_parallel_id_counts(), libMesh::GMVIO::write_binary(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::System::write_parallel_data(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), and libMesh::XdrIO::write_serialized_nodesets().

96  { return cast_int<processor_id_type>(_communicator.size()); }
processor_id_type size() const
Definition: communicator.h:175
const Parallel::Communicator & _communicator

◆ position_of_spectrum()

template<typename T >
PositionOfSpectrum libMesh::EigenSolver< T >::position_of_spectrum ( ) const
inlineinherited
Returns
The position of the spectrum to compute.

Definition at line 140 of file eigen_solver.h.

References libMesh::EigenSolver< T >::_position_of_spectrum.

141  { return _position_of_spectrum;}
PositionOfSpectrum _position_of_spectrum
Definition: eigen_solver.h:290

◆ print_info()

void libMesh::ReferenceCounter::print_info ( std::ostream &  out = libMesh::out)
staticinherited

Prints the reference information, by default to libMesh::out.

Definition at line 87 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter, and libMesh::ReferenceCounter::get_info().

88 {
90  out_stream << ReferenceCounter::get_info();
91 }
static std::string get_info()

◆ processor_id()

processor_id_type libMesh::ParallelObject::processor_id ( ) const
inlineinherited
Returns
The rank of this processor in the group.

Definition at line 101 of file parallel_object.h.

References libMesh::ParallelObject::_communicator, and libMesh::Parallel::Communicator::rank().

Referenced by libMesh::BoundaryInfo::_find_id_maps(), libMesh::EquationSystems::_read_impl(), libMesh::PetscDMWrapper::add_dofs_to_section(), libMesh::DistributedMesh::add_elem(), libMesh::BoundaryInfo::add_elements(), libMesh::DofMap::add_neighbors_to_send_list(), libMesh::DistributedMesh::add_node(), libMesh::UnstructuredMesh::all_second_order(), libMesh::MeshTools::Modification::all_tri(), libMesh::FEMSystem::assembly(), libMesh::EquationSystems::build_discontinuous_solution_vector(), libMesh::Nemesis_IO_Helper::build_element_and_node_maps(), libMesh::InfElemBuilder::build_inf_elem(), libMesh::BoundaryInfo::build_node_list_from_side_list(), libMesh::EquationSystems::build_parallel_elemental_solution_vector(), libMesh::DistributedMesh::clear(), libMesh::ExodusII_IO_Helper::close(), libMesh::Nemesis_IO_Helper::compute_border_node_ids(), libMesh::Nemesis_IO_Helper::compute_communication_map_parameters(), libMesh::Nemesis_IO_Helper::compute_internal_and_border_elems_and_internal_nodes(), libMesh::Nemesis_IO_Helper::compute_node_communication_maps(), libMesh::Nemesis_IO_Helper::compute_num_global_elem_blocks(), libMesh::Nemesis_IO_Helper::compute_num_global_nodesets(), libMesh::Nemesis_IO_Helper::compute_num_global_sidesets(), libMesh::Nemesis_IO_Helper::construct_nemesis_filename(), libMesh::MeshTools::correct_node_proc_ids(), libMesh::ExodusII_IO_Helper::create(), libMesh::DistributedMesh::delete_elem(), libMesh::DistributedMesh::delete_node(), libMesh::MeshCommunication::delete_remote_elements(), libMesh::DofMap::distribute_dofs(), libMesh::DofMap::distribute_local_dofs_node_major(), libMesh::DofMap::distribute_local_dofs_var_major(), libMesh::DistributedMesh::DistributedMesh(), libMesh::DofMap::end_dof(), libMesh::DofMap::end_old_dof(), libMesh::EnsightIO::EnsightIO(), libMesh::MeshFunction::find_element(), libMesh::MeshFunction::find_elements(), libMesh::UnstructuredMesh::find_neighbors(), libMesh::DofMap::first_dof(), libMesh::DofMap::first_old_dof(), libMesh::Nemesis_IO_Helper::get_cmap_params(), libMesh::Nemesis_IO_Helper::get_eb_info_global(), libMesh::Nemesis_IO_Helper::get_elem_cmap(), libMesh::Nemesis_IO_Helper::get_elem_map(), libMesh::MeshBase::get_info(), libMesh::DofMap::get_info(), libMesh::Nemesis_IO_Helper::get_init_global(), libMesh::Nemesis_IO_Helper::get_init_info(), libMesh::Nemesis_IO_Helper::get_loadbal_param(), libMesh::Nemesis_IO_Helper::get_node_cmap(), libMesh::Nemesis_IO_Helper::get_node_map(), libMesh::Nemesis_IO_Helper::get_ns_param_global(), libMesh::Nemesis_IO_Helper::get_ss_param_global(), libMesh::SparsityPattern::Build::handle_vi_vj(), libMesh::SystemSubsetBySubdomain::init(), libMesh::ExodusII_IO_Helper::initialize(), libMesh::ExodusII_IO_Helper::initialize_element_variables(), libMesh::ExodusII_IO_Helper::initialize_global_variables(), libMesh::ExodusII_IO_Helper::initialize_nodal_variables(), libMesh::DistributedMesh::insert_elem(), libMesh::DofMap::is_evaluable(), libMesh::SparsityPattern::Build::join(), libMesh::DofMap::last_dof(), libMesh::MeshTools::libmesh_assert_consistent_distributed(), libMesh::MeshTools::libmesh_assert_consistent_distributed_nodes(), libMesh::MeshTools::libmesh_assert_contiguous_dof_ids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Elem >(), libMesh::MeshTools::libmesh_assert_valid_neighbors(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_object_ids(), libMesh::DofMap::local_variable_indices(), libMesh::MeshRefinement::make_coarsening_compatible(), libMesh::MeshBase::n_active_local_elem(), libMesh::BoundaryInfo::n_boundary_conds(), libMesh::BoundaryInfo::n_edge_conds(), libMesh::DofMap::n_local_dofs(), libMesh::System::n_local_dofs(), libMesh::MeshBase::n_local_elem(), libMesh::MeshBase::n_local_nodes(), libMesh::BoundaryInfo::n_nodeset_conds(), libMesh::BoundaryInfo::n_shellface_conds(), libMesh::SparsityPattern::Build::operator()(), libMesh::DistributedMesh::own_node(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::Nemesis_IO_Helper::put_cmap_params(), libMesh::Nemesis_IO_Helper::put_elem_cmap(), libMesh::Nemesis_IO_Helper::put_elem_map(), libMesh::Nemesis_IO_Helper::put_loadbal_param(), libMesh::Nemesis_IO_Helper::put_node_cmap(), libMesh::Nemesis_IO_Helper::put_node_map(), libMesh::NameBasedIO::read(), libMesh::Nemesis_IO::read(), libMesh::XdrIO::read(), libMesh::CheckpointIO::read(), libMesh::ExodusII_IO_Helper::read_elem_num_map(), libMesh::ExodusII_IO_Helper::read_global_values(), libMesh::CheckpointIO::read_header(), libMesh::XdrIO::read_header(), libMesh::System::read_header(), libMesh::System::read_legacy_data(), libMesh::ExodusII_IO_Helper::read_node_num_map(), libMesh::System::read_parallel_data(), libMesh::System::read_SCALAR_dofs(), libMesh::XdrIO::read_serialized_bc_names(), libMesh::XdrIO::read_serialized_bcs_helper(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::System::read_serialized_data(), libMesh::XdrIO::read_serialized_nodes(), libMesh::XdrIO::read_serialized_nodesets(), libMesh::XdrIO::read_serialized_subdomain_names(), libMesh::System::read_serialized_vector(), libMesh::System::read_serialized_vectors(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::CheckpointIO::select_split_config(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::PetscDMWrapper::set_point_range_in_section(), libMesh::LaplaceMeshSmoother::smooth(), libMesh::MeshTools::total_weight(), libMesh::MeshRefinement::uniformly_coarsen(), libMesh::Parallel::Packing< T >::unpack(), libMesh::DistributedMesh::update_parallel_id_counts(), libMesh::NameBasedIO::write(), libMesh::XdrIO::write(), libMesh::CheckpointIO::write(), libMesh::EquationSystems::write(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::ExodusII_IO::write_element_data(), libMesh::ExodusII_IO_Helper::write_element_values(), libMesh::ExodusII_IO_Helper::write_elements(), libMesh::ExodusII_IO::write_global_data(), libMesh::ExodusII_IO_Helper::write_global_values(), libMesh::System::write_header(), libMesh::ExodusII_IO::write_information_records(), libMesh::ExodusII_IO_Helper::write_information_records(), libMesh::ExodusII_IO_Helper::write_nodal_coordinates(), libMesh::UCDIO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data_discontinuous(), libMesh::ExodusII_IO_Helper::write_nodal_values(), libMesh::Nemesis_IO_Helper::write_nodesets(), libMesh::ExodusII_IO_Helper::write_nodesets(), libMesh::System::write_parallel_data(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bc_names(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::System::write_serialized_data(), libMesh::XdrIO::write_serialized_nodes(), libMesh::XdrIO::write_serialized_nodesets(), libMesh::XdrIO::write_serialized_subdomain_names(), libMesh::System::write_serialized_vector(), libMesh::System::write_serialized_vectors(), libMesh::Nemesis_IO_Helper::write_sidesets(), libMesh::ExodusII_IO_Helper::write_sidesets(), libMesh::ExodusII_IO::write_timestep(), libMesh::ExodusII_IO_Helper::write_timestep(), and libMesh::ExodusII_IO::write_timestep_discontinuous().

102  { return cast_int<processor_id_type>(_communicator.rank()); }
const Parallel::Communicator & _communicator
processor_id_type rank() const
Definition: communicator.h:173

◆ set_close_matrix_before_solve()

template<typename T >
void libMesh::EigenSolver< T >::set_close_matrix_before_solve ( bool  val)
inlineinherited

Set the flag which controls whether libmesh closes the eigenproblem matrices before solving.

Definition at line 111 of file eigen_solver.h.

References libMesh::EigenSolver< T >::_close_matrix_before_solve.

112  {
113  libmesh_experimental();
115  }

◆ set_eigenproblem_type()

template<typename T >
void libMesh::EigenSolver< T >::set_eigenproblem_type ( EigenProblemType  ept)
inlineinherited

Sets the type of the eigenproblem.

Definition at line 152 of file eigen_solver.h.

References libMesh::EigenSolver< T >::_eigen_problem_type.

153  {_eigen_problem_type = ept;}
EigenProblemType _eigen_problem_type
Definition: eigen_solver.h:285

◆ set_eigensolver_type()

template<typename T >
void libMesh::EigenSolver< T >::set_eigensolver_type ( const EigenSolverType  est)
inlineinherited

Sets the type of eigensolver to use.

Definition at line 146 of file eigen_solver.h.

References libMesh::EigenSolver< T >::_eigen_solver_type.

147  { _eigen_solver_type = est; }
EigenSolverType _eigen_solver_type
Definition: eigen_solver.h:280

◆ set_initial_space()

template<typename T >
void libMesh::SlepcEigenSolver< T >::set_initial_space ( NumericVector< T > &  initial_space_in)
overridevirtual

Use initial_space_in as the initial guess.

Implements libMesh::EigenSolver< T >.

Definition at line 839 of file slepc_eigen_solver.C.

References ierr, libMesh::TriangleWrapper::init(), and libMesh::PetscVector< T >::vec().

840 {
841 #if SLEPC_VERSION_LESS_THAN(3,1,0)
842  libmesh_error_msg("SLEPc 3.1 is required to call EigenSolver::set_initial_space()");
843 #else
844  this->init();
845 
846  PetscErrorCode ierr = 0;
847 
848  // Make sure the input vector is actually a PetscVector
849  PetscVector<T> * initial_space_petsc_vec =
850  dynamic_cast<PetscVector<T> *>(&initial_space_in);
851 
852  if (!initial_space_petsc_vec)
853  libmesh_error_msg("Error attaching initial space: input vector must be a PetscVector.");
854 
855  // Get a handle for the underlying Vec.
856  Vec initial_vector = initial_space_petsc_vec->vec();
857 
858  ierr = EPSSetInitialSpace(_eps, 1, &initial_vector);
859  LIBMESH_CHKERR(ierr);
860 #endif
861 }
PetscErrorCode ierr
virtual void init() override

◆ set_position_of_spectrum() [1/3]

template<typename T >
void libMesh::EigenSolver< T >::set_position_of_spectrum ( PositionOfSpectrum  pos)
inlineinherited

Sets the position of the spectrum.

Definition at line 158 of file eigen_solver.h.

References libMesh::EigenSolver< T >::_position_of_spectrum.

159  {_position_of_spectrum= pos;}
PositionOfSpectrum _position_of_spectrum
Definition: eigen_solver.h:290

◆ set_position_of_spectrum() [2/3]

template<typename T >
void libMesh::EigenSolver< T >::set_position_of_spectrum ( Real  pos)
inherited

Definition at line 86 of file eigen_solver.C.

References libMesh::TARGET_MAGNITUDE, and libMesh::TARGET_REAL.

◆ set_position_of_spectrum() [3/3]

template<typename T >
void libMesh::EigenSolver< T >::set_position_of_spectrum ( Real  pos,
PositionOfSpectrum  target 
)
inherited

Definition at line 97 of file eigen_solver.C.

98 {
99  _position_of_spectrum = target;
100  _target_val = pos;
101 }
PositionOfSpectrum _position_of_spectrum
Definition: eigen_solver.h:290

◆ set_slepc_position_of_spectrum()

template<typename T >
void libMesh::SlepcEigenSolver< T >::set_slepc_position_of_spectrum ( )
private

Tells Slepc to compute the spectrum at the position stored in _position_of_spectrum

Definition at line 656 of file slepc_eigen_solver.C.

References ierr, libMesh::LARGEST_IMAGINARY, libMesh::LARGEST_MAGNITUDE, libMesh::LARGEST_REAL, libMesh::SMALLEST_IMAGINARY, libMesh::SMALLEST_MAGNITUDE, libMesh::SMALLEST_REAL, libMesh::TARGET_IMAGINARY, libMesh::TARGET_MAGNITUDE, and libMesh::TARGET_REAL.

657 {
658  PetscErrorCode ierr = 0;
659 
660  switch (this->_position_of_spectrum)
661  {
662  case LARGEST_MAGNITUDE:
663  {
664  ierr = EPSSetWhichEigenpairs (_eps, EPS_LARGEST_MAGNITUDE);
665  LIBMESH_CHKERR(ierr);
666  return;
667  }
668  case SMALLEST_MAGNITUDE:
669  {
670  ierr = EPSSetWhichEigenpairs (_eps, EPS_SMALLEST_MAGNITUDE);
671  LIBMESH_CHKERR(ierr);
672  return;
673  }
674  case LARGEST_REAL:
675  {
676  ierr = EPSSetWhichEigenpairs (_eps, EPS_LARGEST_REAL);
677  LIBMESH_CHKERR(ierr);
678  return;
679  }
680  case SMALLEST_REAL:
681  {
682  ierr = EPSSetWhichEigenpairs (_eps, EPS_SMALLEST_REAL);
683  LIBMESH_CHKERR(ierr);
684  return;
685  }
686  case LARGEST_IMAGINARY:
687  {
688  ierr = EPSSetWhichEigenpairs (_eps, EPS_LARGEST_IMAGINARY);
689  LIBMESH_CHKERR(ierr);
690  return;
691  }
692  case SMALLEST_IMAGINARY:
693  {
694  ierr = EPSSetWhichEigenpairs (_eps, EPS_SMALLEST_IMAGINARY);
695  LIBMESH_CHKERR(ierr);
696  return;
697  }
698 
699  // The EPS_TARGET_XXX enums were added in SLEPc 3.1
700 #if !SLEPC_VERSION_LESS_THAN(3,1,0)
701  case TARGET_MAGNITUDE:
702  {
703  ierr = EPSSetTarget(_eps, this->_target_val);
704  LIBMESH_CHKERR(ierr);
705  ierr = EPSSetWhichEigenpairs (_eps, EPS_TARGET_MAGNITUDE);
706  LIBMESH_CHKERR(ierr);
707  return;
708  }
709  case TARGET_REAL:
710  {
711  ierr = EPSSetTarget(_eps, this->_target_val);
712  LIBMESH_CHKERR(ierr);
713  ierr = EPSSetWhichEigenpairs (_eps, EPS_TARGET_REAL);
714  LIBMESH_CHKERR(ierr);
715  return;
716  }
717  case TARGET_IMAGINARY:
718  {
719  ierr = EPSSetTarget(_eps, this->_target_val);
720  LIBMESH_CHKERR(ierr);
721  ierr = EPSSetWhichEigenpairs (_eps, EPS_TARGET_IMAGINARY);
722  LIBMESH_CHKERR(ierr);
723  return;
724  }
725 #endif
726 
727  default:
728  libmesh_error_msg("ERROR: Unsupported SLEPc position of spectrum: " << this->_position_of_spectrum);
729  }
730 }
PositionOfSpectrum _position_of_spectrum
Definition: eigen_solver.h:290
PetscErrorCode ierr

◆ set_slepc_problem_type()

template<typename T >
void libMesh::SlepcEigenSolver< T >::set_slepc_problem_type ( )
private

Tells Slepc to deal with the type of problem stored in _eigen_problem_type

Definition at line 626 of file slepc_eigen_solver.C.

References libMesh::err, libMesh::GHEP, libMesh::GHIEP, libMesh::GNHEP, libMesh::HEP, ierr, and libMesh::NHEP.

627 {
628  PetscErrorCode ierr = 0;
629 
630  switch (this->_eigen_problem_type)
631  {
632  case NHEP:
633  ierr = EPSSetProblemType (_eps, EPS_NHEP); LIBMESH_CHKERR(ierr); return;
634  case GNHEP:
635  ierr = EPSSetProblemType (_eps, EPS_GNHEP); LIBMESH_CHKERR(ierr); return;
636  case HEP:
637  ierr = EPSSetProblemType (_eps, EPS_HEP); LIBMESH_CHKERR(ierr); return;
638  case GHEP:
639  ierr = EPSSetProblemType (_eps, EPS_GHEP); LIBMESH_CHKERR(ierr); return;
640 #if !SLEPC_VERSION_LESS_THAN(3,3,0)
641  // EPS_GHIEP added in 3.3.0
642  case GHIEP:
643  ierr = EPSSetProblemType (_eps, EPS_GHIEP); LIBMESH_CHKERR(ierr); return;
644 #endif
645 
646  default:
647  libMesh::err << "ERROR: Unsupported SLEPc Eigen Problem: "
648  << this->_eigen_problem_type << std::endl
649  << "Continuing with SLEPc defaults" << std::endl;
650  }
651 }
OStreamProxy err(std::cerr)
PetscErrorCode ierr
EigenProblemType _eigen_problem_type
Definition: eigen_solver.h:285

◆ set_slepc_solver_type()

template<typename T >
void libMesh::SlepcEigenSolver< T >::set_slepc_solver_type ( )
private

Tells Slepc to use the user-specified solver stored in _eigen_solver_type

Definition at line 594 of file slepc_eigen_solver.C.

References libMesh::ARNOLDI, libMesh::Utility::enum_to_string(), libMesh::err, ierr, libMesh::KRYLOVSCHUR, libMesh::LANCZOS, libMesh::LAPACK, libMesh::POWER, and libMesh::SUBSPACE.

595 {
596  PetscErrorCode ierr = 0;
597 
598  switch (this->_eigen_solver_type)
599  {
600  case POWER:
601  ierr = EPSSetType (_eps, EPSPOWER); LIBMESH_CHKERR(ierr); return;
602  case SUBSPACE:
603  ierr = EPSSetType (_eps, EPSSUBSPACE); LIBMESH_CHKERR(ierr); return;
604  case LAPACK:
605  ierr = EPSSetType (_eps, EPSLAPACK); LIBMESH_CHKERR(ierr); return;
606  case ARNOLDI:
607  ierr = EPSSetType (_eps, EPSARNOLDI); LIBMESH_CHKERR(ierr); return;
608  case LANCZOS:
609  ierr = EPSSetType (_eps, EPSLANCZOS); LIBMESH_CHKERR(ierr); return;
610  case KRYLOVSCHUR:
611  ierr = EPSSetType (_eps, EPSKRYLOVSCHUR); LIBMESH_CHKERR(ierr); return;
612  // case ARPACK:
613  // ierr = EPSSetType (_eps, (char *) EPSARPACK); LIBMESH_CHKERR(ierr); return;
614 
615  default:
616  libMesh::err << "ERROR: Unsupported SLEPc Eigen Solver: "
617  << Utility::enum_to_string(this->_eigen_solver_type) << std::endl
618  << "Continuing with SLEPc defaults" << std::endl;
619  }
620 }
OStreamProxy err(std::cerr)
std::string enum_to_string(const T e)
PetscErrorCode ierr
EigenSolverType _eigen_solver_type
Definition: eigen_solver.h:280

◆ set_solver_configuration()

template<typename T >
void libMesh::EigenSolver< T >::set_solver_configuration ( SolverConfiguration solver_configuration)
inherited

Set the solver configuration object.

Definition at line 80 of file eigen_solver.C.

81 {
82  _solver_configuration = &solver_configuration;
83 }
SolverConfiguration * _solver_configuration
Definition: eigen_solver.h:301

◆ solve_generalized() [1/4]

template<typename T >
std::pair< unsigned int, unsigned int > libMesh::SlepcEigenSolver< T >::solve_generalized ( SparseMatrix< T > &  matrix_A,
SparseMatrix< T > &  matrix_B,
int  nev,
int  ncv,
const double  tol,
const unsigned int  m_its 
)
overridevirtual

This function calls the SLEPc solver to compute the eigenpairs for the generalized eigenproblem defined by the matrix_A and matrix_B, which are of type SparseMatrix. The argument nev is the number of eigenpairs to be computed and ncv is the number of basis vectors to be used in the solution procedure. Return values are the number of converged eigen values and the number of the iterations carried out by the eigen solver.

Implements libMesh::EigenSolver< T >.

Definition at line 290 of file slepc_eigen_solver.C.

References libMesh::PetscMatrix< T >::close(), libMesh::TriangleWrapper::init(), and libMesh::PetscMatrix< T >::mat().

296 {
297  this->init ();
298 
299  // Make sure the data passed in are really of Petsc types
300  PetscMatrix<T> * matrix_A = dynamic_cast<PetscMatrix<T> *>(&matrix_A_in);
301  PetscMatrix<T> * matrix_B = dynamic_cast<PetscMatrix<T> *>(&matrix_B_in);
302 
303  if (!matrix_A || !matrix_B)
304  libmesh_error_msg("Error: inputs to solve_generalized() must be of type PetscMatrix.");
305 
306  // Close the matrix and vectors in case this wasn't already done.
307  if (this->_close_matrix_before_solve)
308  {
309  matrix_A->close ();
310  matrix_B->close ();
311  }
312 
313  return _solve_generalized_helper (matrix_A->mat(), matrix_B->mat(), nev, ncv, tol, m_its);
314 }
std::pair< unsigned int, unsigned int > _solve_generalized_helper(Mat mat_A, Mat mat_B, int nev, int ncv, const double tol, const unsigned int m_its)
virtual void init() override

◆ solve_generalized() [2/4]

template<typename T >
std::pair< unsigned int, unsigned int > libMesh::SlepcEigenSolver< T >::solve_generalized ( ShellMatrix< T > &  matrix_A,
SparseMatrix< T > &  matrix_B,
int  nev,
int  ncv,
const double  tol,
const unsigned int  m_its 
)
overridevirtual

Solve generalized eigenproblem when matrix_A is of type ShellMatrix, matrix_B is of type SparseMatrix.

Implements libMesh::EigenSolver< T >.

Definition at line 318 of file slepc_eigen_solver.C.

References libMesh::PetscMatrix< T >::close(), ierr, libMesh::TriangleWrapper::init(), libMesh::ShellMatrix< T >::m(), libMesh::PetscMatrix< T >::mat(), and libMesh::ShellMatrix< T >::n().

324 {
325  this->init ();
326 
327  PetscErrorCode ierr=0;
328 
329  // Prepare the matrix. Note that the const_cast is only necessary
330  // because PETSc does not accept a const void *. Inside the member
331  // function _petsc_shell_matrix() below, the pointer is casted back
332  // to a const ShellMatrix<T> *.
333  Mat mat_A;
334  ierr = MatCreateShell(this->comm().get(),
335  shell_matrix_A.m(), // Specify the number of local rows
336  shell_matrix_A.n(), // Specify the number of local columns
337  PETSC_DETERMINE,
338  PETSC_DETERMINE,
339  const_cast<void *>(static_cast<const void *>(&shell_matrix_A)),
340  &mat_A);
341  LIBMESH_CHKERR(ierr);
342 
343  PetscMatrix<T> * matrix_B = dynamic_cast<PetscMatrix<T> *>(&matrix_B_in);
344 
345  if (!matrix_B)
346  libmesh_error_msg("Error: inputs to solve_generalized() must be of type PetscMatrix.");
347 
348  // Close the matrix and vectors in case this wasn't already done.
349  if (this->_close_matrix_before_solve)
350  matrix_B->close ();
351 
352  ierr = MatShellSetOperation(mat_A,MATOP_MULT,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_mult));
353  LIBMESH_CHKERR(ierr);
354  ierr = MatShellSetOperation(mat_A,MATOP_GET_DIAGONAL,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_get_diagonal));
355  LIBMESH_CHKERR(ierr);
356 
357  return _solve_generalized_helper (mat_A, matrix_B->mat(), nev, ncv, tol, m_its);
358 }
std::pair< unsigned int, unsigned int > _solve_generalized_helper(Mat mat_A, Mat mat_B, int nev, int ncv, const double tol, const unsigned int m_its)
const Parallel::Communicator & comm() const
static PetscErrorCode _petsc_shell_matrix_get_diagonal(Mat mat, Vec dest)
PetscErrorCode ierr
static PetscErrorCode _petsc_shell_matrix_mult(Mat mat, Vec arg, Vec dest)
virtual void init() override

◆ solve_generalized() [3/4]

template<typename T >
std::pair< unsigned int, unsigned int > libMesh::SlepcEigenSolver< T >::solve_generalized ( SparseMatrix< T > &  matrix_A,
ShellMatrix< T > &  matrix_B,
int  nev,
int  ncv,
const double  tol,
const unsigned int  m_its 
)
overridevirtual

Solve generalized eigenproblem when matrix_A is of type SparseMatrix, matrix_B is of type ShellMatrix. When using this function, one should use the command line options: -st_ksp_type gmres -st_pc_type none or -st_ksp_type gmres -st_pc_type jacobi or similar.

Implements libMesh::EigenSolver< T >.

Definition at line 362 of file slepc_eigen_solver.C.

References libMesh::PetscMatrix< T >::close(), ierr, libMesh::TriangleWrapper::init(), libMesh::ShellMatrix< T >::m(), libMesh::PetscMatrix< T >::mat(), and libMesh::ShellMatrix< T >::n().

368 {
369  this->init ();
370 
371  PetscErrorCode ierr=0;
372 
373  PetscMatrix<T> * matrix_A = dynamic_cast<PetscMatrix<T> *>(&matrix_A_in);
374 
375  if (!matrix_A)
376  libmesh_error_msg("Error: inputs to solve_generalized() must be of type PetscMatrix.");
377 
378  // Close the matrix and vectors in case this wasn't already done.
379  if (this->_close_matrix_before_solve)
380  matrix_A->close ();
381 
382  // Prepare the matrix. Note that the const_cast is only necessary
383  // because PETSc does not accept a const void *. Inside the member
384  // function _petsc_shell_matrix() below, the pointer is casted back
385  // to a const ShellMatrix<T> *.
386  Mat mat_B;
387  ierr = MatCreateShell(this->comm().get(),
388  shell_matrix_B.m(), // Specify the number of local rows
389  shell_matrix_B.n(), // Specify the number of local columns
390  PETSC_DETERMINE,
391  PETSC_DETERMINE,
392  const_cast<void *>(static_cast<const void *>(&shell_matrix_B)),
393  &mat_B);
394  LIBMESH_CHKERR(ierr);
395 
396 
397  ierr = MatShellSetOperation(mat_B,MATOP_MULT,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_mult));
398  LIBMESH_CHKERR(ierr);
399  ierr = MatShellSetOperation(mat_B,MATOP_GET_DIAGONAL,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_get_diagonal));
400  LIBMESH_CHKERR(ierr);
401 
402  return _solve_generalized_helper (matrix_A->mat(), mat_B, nev, ncv, tol, m_its);
403 }
std::pair< unsigned int, unsigned int > _solve_generalized_helper(Mat mat_A, Mat mat_B, int nev, int ncv, const double tol, const unsigned int m_its)
const Parallel::Communicator & comm() const
static PetscErrorCode _petsc_shell_matrix_get_diagonal(Mat mat, Vec dest)
PetscErrorCode ierr
static PetscErrorCode _petsc_shell_matrix_mult(Mat mat, Vec arg, Vec dest)
virtual void init() override

◆ solve_generalized() [4/4]

template<typename T >
std::pair< unsigned int, unsigned int > libMesh::SlepcEigenSolver< T >::solve_generalized ( ShellMatrix< T > &  matrix_A,
ShellMatrix< T > &  matrix_B,
int  nev,
int  ncv,
const double  tol,
const unsigned int  m_its 
)
overridevirtual

Solve generalized eigenproblem when both matrix_A and matrix_B are of type ShellMatrix. When using this function, one should use the command line options: -st_ksp_type gmres -st_pc_type none or -st_ksp_type gmres -st_pc_type jacobi or similar.

Implements libMesh::EigenSolver< T >.

Definition at line 407 of file slepc_eigen_solver.C.

References ierr, libMesh::TriangleWrapper::init(), libMesh::ShellMatrix< T >::m(), and libMesh::ShellMatrix< T >::n().

413 {
414  this->init ();
415 
416  PetscErrorCode ierr=0;
417 
418  // Prepare the matrices. Note that the const_casts are only
419  // necessary because PETSc does not accept a const void *. Inside
420  // the member function _petsc_shell_matrix() below, the pointer is
421  // casted back to a const ShellMatrix<T> *.
422  Mat mat_A;
423  ierr = MatCreateShell(this->comm().get(),
424  shell_matrix_A.m(), // Specify the number of local rows
425  shell_matrix_A.n(), // Specify the number of local columns
426  PETSC_DETERMINE,
427  PETSC_DETERMINE,
428  const_cast<void *>(static_cast<const void *>(&shell_matrix_A)),
429  &mat_A);
430  LIBMESH_CHKERR(ierr);
431 
432  Mat mat_B;
433  ierr = MatCreateShell(this->comm().get(),
434  shell_matrix_B.m(), // Specify the number of local rows
435  shell_matrix_B.n(), // Specify the number of local columns
436  PETSC_DETERMINE,
437  PETSC_DETERMINE,
438  const_cast<void *>(static_cast<const void *>(&shell_matrix_B)),
439  &mat_B);
440  LIBMESH_CHKERR(ierr);
441 
442  ierr = MatShellSetOperation(mat_A,MATOP_MULT,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_mult));
443  LIBMESH_CHKERR(ierr);
444  ierr = MatShellSetOperation(mat_A,MATOP_GET_DIAGONAL,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_get_diagonal));
445  LIBMESH_CHKERR(ierr);
446 
447  ierr = MatShellSetOperation(mat_B,MATOP_MULT,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_mult));
448  LIBMESH_CHKERR(ierr);
449  ierr = MatShellSetOperation(mat_B,MATOP_GET_DIAGONAL,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_get_diagonal));
450  LIBMESH_CHKERR(ierr);
451 
452  return _solve_generalized_helper (mat_A, mat_B, nev, ncv, tol, m_its);
453 }
std::pair< unsigned int, unsigned int > _solve_generalized_helper(Mat mat_A, Mat mat_B, int nev, int ncv, const double tol, const unsigned int m_its)
const Parallel::Communicator & comm() const
static PetscErrorCode _petsc_shell_matrix_get_diagonal(Mat mat, Vec dest)
PetscErrorCode ierr
static PetscErrorCode _petsc_shell_matrix_mult(Mat mat, Vec arg, Vec dest)
virtual void init() override

◆ solve_standard() [1/2]

template<typename T >
std::pair< unsigned int, unsigned int > libMesh::SlepcEigenSolver< T >::solve_standard ( SparseMatrix< T > &  matrix_A,
int  nev,
int  ncv,
const double  tol,
const unsigned int  m_its 
)
overridevirtual

This function calls the SLEPc solver to compute the eigenpairs of the SparseMatrix matrix_A. nev is the number of eigenpairs to be computed and ncv is the number of basis vectors to be used in the solution procedure. Return values are the number of converged eigen values and the number of the iterations carried out by the eigen solver.

Implements libMesh::EigenSolver< T >.

Definition at line 103 of file slepc_eigen_solver.C.

References libMesh::PetscMatrix< T >::close(), libMesh::TriangleWrapper::init(), and libMesh::PetscMatrix< T >::mat().

108 {
109  LOG_SCOPE("solve_standard()", "SlepcEigenSolver");
110 
111  this->init ();
112 
113  // Make sure the SparseMatrix passed in is really a PetscMatrix
114  PetscMatrix<T> * matrix_A = dynamic_cast<PetscMatrix<T> *>(&matrix_A_in);
115 
116  if (!matrix_A)
117  libmesh_error_msg("Error: input matrix to solve_standard() must be a PetscMatrix.");
118 
119  // Close the matrix and vectors in case this wasn't already done.
120  if (this->_close_matrix_before_solve)
121  matrix_A->close ();
122 
123  return _solve_standard_helper(matrix_A->mat(), nev, ncv, tol, m_its);
124 }
std::pair< unsigned int, unsigned int > _solve_standard_helper(Mat mat, int nev, int ncv, const double tol, const unsigned int m_its)
virtual void init() override

◆ solve_standard() [2/2]

template<typename T >
std::pair< unsigned int, unsigned int > libMesh::SlepcEigenSolver< T >::solve_standard ( ShellMatrix< T > &  shell_matrix,
int  nev,
int  ncv,
const double  tol,
const unsigned int  m_its 
)
overridevirtual

Same as above except that matrix_A is a ShellMatrix in this case.

Implements libMesh::EigenSolver< T >.

Definition at line 129 of file slepc_eigen_solver.C.

References ierr, libMesh::TriangleWrapper::init(), libMesh::ShellMatrix< T >::m(), and libMesh::ShellMatrix< T >::n().

134 {
135  this->init ();
136 
137  PetscErrorCode ierr=0;
138 
139  // Prepare the matrix. Note that the const_cast is only necessary
140  // because PETSc does not accept a const void *. Inside the member
141  // function _petsc_shell_matrix() below, the pointer is casted back
142  // to a const ShellMatrix<T> *.
143  Mat mat;
144  ierr = MatCreateShell(this->comm().get(),
145  shell_matrix.m(), // Specify the number of local rows
146  shell_matrix.n(), // Specify the number of local columns
147  PETSC_DETERMINE,
148  PETSC_DETERMINE,
149  const_cast<void *>(static_cast<const void *>(&shell_matrix)),
150  &mat);
151  LIBMESH_CHKERR(ierr);
152 
153  ierr = MatShellSetOperation(mat,MATOP_MULT,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_mult));
154  LIBMESH_CHKERR(ierr);
155  ierr = MatShellSetOperation(mat,MATOP_GET_DIAGONAL,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_get_diagonal));
156  LIBMESH_CHKERR(ierr);
157 
158  return _solve_standard_helper(mat, nev, ncv, tol, m_its);
159 }
const Parallel::Communicator & comm() const
static PetscErrorCode _petsc_shell_matrix_get_diagonal(Mat mat, Vec dest)
PetscErrorCode ierr
static PetscErrorCode _petsc_shell_matrix_mult(Mat mat, Vec arg, Vec dest)
std::pair< unsigned int, unsigned int > _solve_standard_helper(Mat mat, int nev, int ncv, const double tol, const unsigned int m_its)
virtual void init() override

Member Data Documentation

◆ _close_matrix_before_solve

template<typename T >
bool libMesh::EigenSolver< T >::_close_matrix_before_solve
protectedinherited

◆ _communicator

◆ _counts

ReferenceCounter::Counts libMesh::ReferenceCounter::_counts
staticprotectedinherited

◆ _eigen_problem_type

template<typename T >
EigenProblemType libMesh::EigenSolver< T >::_eigen_problem_type
protectedinherited

Enum stating which type of eigen problem we deal with.

Definition at line 285 of file eigen_solver.h.

Referenced by libMesh::EigenSolver< T >::eigen_problem_type(), libMesh::EigenSolver< T >::set_eigenproblem_type(), and libMesh::SlepcEigenSolver< T >::SlepcEigenSolver().

◆ _eigen_solver_type

template<typename T >
EigenSolverType libMesh::EigenSolver< T >::_eigen_solver_type
protectedinherited

◆ _enable_print_counter

bool libMesh::ReferenceCounter::_enable_print_counter = true
staticprotectedinherited

Flag to control whether reference count information is printed when print_info is called.

Definition at line 141 of file reference_counter.h.

Referenced by libMesh::ReferenceCounter::disable_print_counter_info(), libMesh::ReferenceCounter::enable_print_counter_info(), and libMesh::ReferenceCounter::print_info().

◆ _eps

template<typename T >
EPS libMesh::SlepcEigenSolver< T >::_eps
private

Eigenproblem solver context

Definition at line 272 of file slepc_eigen_solver.h.

Referenced by libMesh::SlepcEigenSolver< T >::eps().

◆ _is_initialized

template<typename T >
bool libMesh::EigenSolver< T >::_is_initialized
protectedinherited

Flag indicating if the data structures have been initialized.

Definition at line 295 of file eigen_solver.h.

Referenced by libMesh::EigenSolver< T >::initialized().

◆ _mutex

Threads::spin_mutex libMesh::ReferenceCounter::_mutex
staticprotectedinherited

Mutual exclusion object to enable thread-safe reference counting.

Definition at line 135 of file reference_counter.h.

◆ _n_objects

Threads::atomic< unsigned int > libMesh::ReferenceCounter::_n_objects
staticprotectedinherited

The number of objects. Print the reference count information when the number returns to 0.

Definition at line 130 of file reference_counter.h.

Referenced by libMesh::ReferenceCounter::n_objects(), libMesh::ReferenceCounter::ReferenceCounter(), and libMesh::ReferenceCounter::~ReferenceCounter().

◆ _position_of_spectrum

template<typename T >
PositionOfSpectrum libMesh::EigenSolver< T >::_position_of_spectrum
protectedinherited

Enum stating where to evaluate the spectrum.

Definition at line 290 of file eigen_solver.h.

Referenced by libMesh::EigenSolver< T >::position_of_spectrum(), and libMesh::EigenSolver< T >::set_position_of_spectrum().

◆ _solver_configuration

template<typename T >
SolverConfiguration* libMesh::EigenSolver< T >::_solver_configuration
protectedinherited

Optionally store a SolverOptions object that can be used to set parameters like solver type, tolerances and iteration limits.

Definition at line 301 of file eigen_solver.h.

◆ _target_val

template<typename T >
Real libMesh::EigenSolver< T >::_target_val
protectedinherited

Definition at line 303 of file eigen_solver.h.


The documentation for this class was generated from the following files: