libMesh::EigenSolver< T > Class Template Referenceabstract

Base class which defines the interface for solving eigenproblems. More...

#include <eigen_solver.h>

Inheritance diagram for libMesh::EigenSolver< T >:

Public Member Functions

 EigenSolver (const Parallel::Communicator &comm_in)
 
virtual ~EigenSolver ()
 
bool initialized () const
 
bool get_close_matrix_before_solve () const
 
void set_close_matrix_before_solve (bool val)
 
virtual void clear ()
 
virtual void init ()=0
 
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)
 
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)=0
 
virtual std::pair< unsigned int, unsigned int > solve_standard (ShellMatrix< T > &matrix_A, int nev, int ncv, const double tol, const unsigned int m_its)=0
 
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)=0
 
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)=0
 
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)=0
 
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)=0
 
virtual std::pair< Real, Realget_eigenpair (dof_id_type i, NumericVector< T > &solution)=0
 
virtual std::pair< Real, Realget_eigenvalue (dof_id_type i)=0
 
virtual void attach_deflation_space (NumericVector< T > &deflation_vector)=0
 
virtual void set_initial_space (NumericVector< T > &initial_space_in)=0
 
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
 

Detailed Description

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

Base class which defines the interface for solving eigenproblems.

This class provides an interface to solvers for eigenvalue problems.

Author
Steffen Peterson
Date
2005

Definition at line 68 of file 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

◆ EigenSolver()

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

Constructor. Initializes Solver data structures

Definition at line 36 of file eigen_solver.C.

36  :
37  ParallelObject(comm_in),
41  _is_initialized (false),
42  _solver_configuration(nullptr),
44 {
45 }
ParallelObject(const Parallel::Communicator &comm_in)
SolverConfiguration * _solver_configuration
Definition: eigen_solver.h:301
PositionOfSpectrum _position_of_spectrum
Definition: eigen_solver.h:290
EigenSolverType _eigen_solver_type
Definition: eigen_solver.h:280
EigenProblemType _eigen_problem_type
Definition: eigen_solver.h:285

◆ ~EigenSolver()

template<typename T >
libMesh::EigenSolver< T >::~EigenSolver ( )
virtual

Destructor.

Definition at line 50 of file eigen_solver.C.

51 {
52  this->clear ();
53 }
virtual void clear()
Definition: eigen_solver.h:120

Member Function Documentation

◆ attach_deflation_space()

template<typename T >
virtual void libMesh::EigenSolver< T >::attach_deflation_space ( NumericVector< T > &  deflation_vector)
pure virtual

Attach a deflation space defined by a single vector.

Implemented in libMesh::SlepcEigenSolver< T >.

◆ build()

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

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 >
virtual void libMesh::EigenSolver< T >::clear ( )
inlinevirtual

Release all memory and clear data structures.

Reimplemented in libMesh::SlepcEigenSolver< T >.

Definition at line 120 of file eigen_solver.h.

120 {}

◆ 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
inline
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
inline
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 }

◆ get_close_matrix_before_solve()

template<typename T >
bool libMesh::EigenSolver< T >::get_close_matrix_before_solve ( ) const
inline
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 >
virtual std::pair<Real, Real> libMesh::EigenSolver< T >::get_eigenpair ( dof_id_type  i,
NumericVector< T > &  solution 
)
pure virtual
Returns
The ith eigenvalue (real and imaginary part), and copies the ith eigenvector into the solution vector.

Implemented in libMesh::SlepcEigenSolver< T >.

◆ get_eigenvalue()

template<typename T >
virtual std::pair<Real, Real> libMesh::EigenSolver< T >::get_eigenvalue ( dof_id_type  i)
pure virtual
Returns
The ith eigenvalue (real and imaginary part).

Same as above function, except it does not copy the eigenvector.

Implemented in libMesh::SlepcEigenSolver< T >.

◆ 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

◆ 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 >
virtual void libMesh::EigenSolver< T >::init ( )
pure virtual

Initialize data structures if not done so already.

Implemented in libMesh::SlepcEigenSolver< T >.

◆ initialized()

template<typename T >
bool libMesh::EigenSolver< T >::initialized ( ) const
inline
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
inline
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)
inline

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)
inline

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)
inline

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 >
virtual void libMesh::EigenSolver< T >::set_initial_space ( NumericVector< T > &  initial_space_in)
pure virtual

Provide one basis vector for the initial guess

Implemented in libMesh::SlepcEigenSolver< T >.

◆ set_position_of_spectrum() [1/3]

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

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)

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 
)

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

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

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 >
virtual std::pair<unsigned int, unsigned int> libMesh::EigenSolver< T >::solve_generalized ( SparseMatrix< T > &  matrix_A,
SparseMatrix< T > &  matrix_B,
int  nev,
int  ncv,
const double  tol,
const unsigned int  m_its 
)
pure virtual

Solves the generalized eigenproblem involving SparseMatrices matrix_A and matrix_B.

Returns
The number of converged eigenpairs and the number of iterations.

Implemented in libMesh::SlepcEigenSolver< T >.

◆ solve_generalized() [2/4]

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

Solves the generalized eigenproblem with ShellMatrix matrix_A and SparseMatrix matrix_B.

Returns
The number of converged eigenpairs and the number of iterations.

Implemented in libMesh::SlepcEigenSolver< T >.

◆ solve_generalized() [3/4]

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

Solves the generalized eigenproblem with SparseMatrix matrix_A and ShellMatrix matrix_B.

Returns
The number of converged eigenpairs and the number of iterations.

Implemented in libMesh::SlepcEigenSolver< T >.

◆ solve_generalized() [4/4]

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

Solves the generalized eigenproblem involving ShellMatrices matrix_A and matrix_B.

Returns
The number of converged eigenpairs and the number of iterations.

Implemented in libMesh::SlepcEigenSolver< T >.

◆ solve_standard() [1/2]

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

Solves the standard eigenproblem involving the SparseMatrix matrix_A.

Returns
The number of converged eigenpairs and the number of iterations.

Implemented in libMesh::SlepcEigenSolver< T >.

◆ solve_standard() [2/2]

template<typename T >
virtual std::pair<unsigned int, unsigned int> libMesh::EigenSolver< T >::solve_standard ( ShellMatrix< T > &  matrix_A,
int  nev,
int  ncv,
const double  tol,
const unsigned int  m_its 
)
pure virtual

Solves the standard eigenproblem involving the ShellMatrix matrix_A.

Returns
The number of converged eigenpairs and the number of iterations.

Implemented in libMesh::SlepcEigenSolver< T >.

Member Data Documentation

◆ _close_matrix_before_solve

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

◆ _communicator

◆ _counts

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

◆ _eigen_problem_type

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

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
protected

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

◆ _is_initialized

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

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
protected

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
protected

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
protected

Definition at line 303 of file eigen_solver.h.


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