libMesh::EigenSparseMatrix< T > Class Template Referencefinal

#include <eigen_sparse_matrix.h>

Inheritance diagram for libMesh::EigenSparseMatrix< T >:

Public Types

typedef EigenSM DataType
 
typedef Eigen::Triplet< T, eigen_idx_typeTripletType
 

Public Member Functions

 EigenSparseMatrix (const Parallel::Communicator &comm)
 
 EigenSparseMatrix (EigenSparseMatrix &&)=default
 
 EigenSparseMatrix (const EigenSparseMatrix &)=default
 
EigenSparseMatrixoperator= (const EigenSparseMatrix &)=default
 
EigenSparseMatrixoperator= (EigenSparseMatrix &&)=default
 
virtual ~EigenSparseMatrix ()=default
 
virtual void init (const numeric_index_type m, const numeric_index_type n, const numeric_index_type m_l, const numeric_index_type n_l, const numeric_index_type nnz=30, const numeric_index_type noz=10, const numeric_index_type blocksize=1) override
 
virtual void init () override
 
virtual void clear () override
 
virtual void zero () override
 
virtual void close () override
 
virtual numeric_index_type m () const override
 
virtual numeric_index_type n () const override
 
virtual numeric_index_type row_start () const override
 
virtual numeric_index_type row_stop () const override
 
virtual void set (const numeric_index_type i, const numeric_index_type j, const T value) override
 
virtual void add (const numeric_index_type i, const numeric_index_type j, const T value) override
 
virtual void add_matrix (const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols) override
 
virtual void add_matrix (const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &dof_indices) override
 
virtual void add (const T a, const SparseMatrix< T > &X) override
 
virtual T operator() (const numeric_index_type i, const numeric_index_type j) const override
 
virtual Real l1_norm () const override
 
virtual Real linfty_norm () const override
 
virtual bool closed () const override
 
virtual void print_personal (std::ostream &os=libMesh::out) const override
 
virtual void get_diagonal (NumericVector< T > &dest) const override
 
virtual void get_transpose (SparseMatrix< T > &dest) const override
 
virtual bool initialized () const
 
void attach_dof_map (const DofMap &dof_map)
 
virtual bool need_full_sparsity_pattern () const
 
virtual void update_sparsity_pattern (const SparsityPattern::Graph &)
 
virtual void zero_rows (std::vector< numeric_index_type > &rows, T diag_value=0.0)
 
virtual void flush ()
 
virtual void add_block_matrix (const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &brows, const std::vector< numeric_index_type > &bcols)
 
virtual void add_block_matrix (const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &dof_indices)
 
void print (std::ostream &os=libMesh::out, const bool sparse=false) const
 
template<>
void print (std::ostream &os, const bool sparse) const
 
virtual void print_matlab (const std::string &="") const
 
virtual void create_submatrix (SparseMatrix< T > &submatrix, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols) const
 
virtual void reinit_submatrix (SparseMatrix< T > &submatrix, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols) const
 
void vector_mult (NumericVector< T > &dest, const NumericVector< T > &arg) const
 
void vector_mult_add (NumericVector< T > &dest, const NumericVector< T > &arg) const
 
const Parallel::Communicatorcomm () const
 
processor_id_type n_processors () const
 
processor_id_type processor_id () const
 

Static Public Member Functions

static std::unique_ptr< SparseMatrix< T > > build (const Parallel::Communicator &comm, const SolverPackage solver_package=libMesh::default_solver_package())
 
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

virtual void _get_submatrix (SparseMatrix< T > &, const std::vector< numeric_index_type > &, const std::vector< numeric_index_type > &, const bool) const
 
void increment_constructor_count (const std::string &name)
 
void increment_destructor_count (const std::string &name)
 

Protected Attributes

DofMap const * _dof_map
 
bool _is_initialized
 
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 Attributes

DataType _mat
 
bool _closed
 

Friends

class EigenSparseVector< T >
 
class EigenSparseLinearSolver< T >
 

Detailed Description

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

The EigenSparseMatrix class wraps a sparse matrix object from the Eigen library. All overridden virtual functions are documented in sparse_matrix.h.

Author
Benjamin S. Kirk
Date
2013

Definition at line 54 of file eigen_sparse_matrix.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.

◆ DataType

template<typename T>
typedef EigenSM libMesh::EigenSparseMatrix< T >::DataType

Convenient typedefs

Definition at line 83 of file eigen_sparse_matrix.h.

◆ TripletType

template<typename T>
typedef Eigen::Triplet<T,eigen_idx_type> libMesh::EigenSparseMatrix< T >::TripletType

Definition at line 84 of file eigen_sparse_matrix.h.

Constructor & Destructor Documentation

◆ EigenSparseMatrix() [1/3]

template<typename T >
libMesh::EigenSparseMatrix< T >::EigenSparseMatrix ( const Parallel::Communicator comm)

Constructor; initializes the matrix to be empty, without any structure, i.e. the matrix is not usable at all. This constructor is therefore only useful for matrices which are members of a class. All other matrices should be created at a point in the data flow where all necessary information is available.

You have to initialize the matrix before usage with init(...).

Definition at line 161 of file eigen_sparse_matrix.C.

161  :
162  SparseMatrix<T>(comm_in),
163  _closed (false)
164 {
165 }

◆ EigenSparseMatrix() [2/3]

template<typename T>
libMesh::EigenSparseMatrix< T >::EigenSparseMatrix ( EigenSparseMatrix< T > &&  )
default

The 5 special functions can be defaulted for this class, as it does not manage any memory itself.

◆ EigenSparseMatrix() [3/3]

template<typename T>
libMesh::EigenSparseMatrix< T >::EigenSparseMatrix ( const EigenSparseMatrix< T > &  )
default

◆ ~EigenSparseMatrix()

template<typename T>
virtual libMesh::EigenSparseMatrix< T >::~EigenSparseMatrix ( )
virtualdefault

Member Function Documentation

◆ _get_submatrix()

template<typename T>
virtual void libMesh::SparseMatrix< T >::_get_submatrix ( SparseMatrix< T > &  ,
const std::vector< numeric_index_type > &  ,
const std::vector< numeric_index_type > &  ,
const bool   
) const
inlineprotectedvirtualinherited

Protected implementation of the create_submatrix and reinit_submatrix routines.

Note
This function must be overridden in derived classes for it to work properly!

Reimplemented in libMesh::PetscMatrix< T >.

Definition at line 414 of file sparse_matrix.h.

Referenced by libMesh::SparseMatrix< ValOut >::create_submatrix(), and libMesh::SparseMatrix< ValOut >::reinit_submatrix().

418  {
419  libmesh_not_implemented();
420  }

◆ add() [1/2]

template<typename T >
void libMesh::EigenSparseMatrix< T >::add ( const numeric_index_type  i,
const numeric_index_type  j,
const T  value 
)
overridevirtual

Add value to the element (i,j). Throws an error if the entry does not exist. Zero values can be "added" to non-existent entries.

Implements libMesh::SparseMatrix< T >.

Definition at line 239 of file eigen_sparse_matrix.C.

References libMesh::initialized(), and value.

242 {
243  libmesh_assert (this->initialized());
244  libmesh_assert_less (i, this->m());
245  libmesh_assert_less (j, this->n());
246 
247  _mat.coeffRef(i,j) += value;
248 }
virtual bool initialized() const
virtual numeric_index_type m() const override
virtual numeric_index_type n() const override
static const bool value
Definition: xdr_io.C:109

◆ add() [2/2]

template<typename T >
void libMesh::EigenSparseMatrix< T >::add ( const T  a,
const SparseMatrix< T > &  X 
)
overridevirtual

Compute $ A \leftarrow A + a*X $ for scalar a, matrix X.

Implements libMesh::SparseMatrix< T >.

Definition at line 262 of file eigen_sparse_matrix.C.

References libMesh::EigenSparseMatrix< T >::_mat, libMesh::initialized(), libMesh::SparseMatrix< T >::m(), and libMesh::SparseMatrix< T >::n().

263 {
264  libmesh_assert (this->initialized());
265  libmesh_assert_equal_to (this->m(), X_in.m());
266  libmesh_assert_equal_to (this->n(), X_in.n());
267 
268  const EigenSparseMatrix<T> & X =
269  cast_ref<const EigenSparseMatrix<T> &> (X_in);
270 
271  _mat += X._mat*a_in;
272 }
virtual bool initialized() const
virtual numeric_index_type m() const override
virtual numeric_index_type n() const override

◆ add_block_matrix() [1/2]

template<typename T>
void libMesh::SparseMatrix< T >::add_block_matrix ( const DenseMatrix< T > &  dm,
const std::vector< numeric_index_type > &  brows,
const std::vector< numeric_index_type > &  bcols 
)
virtualinherited

Add the full matrix dm to the SparseMatrix. This is useful for adding an element matrix at assembly time. The matrix is assumed blocked, and brow, bcol correspond to the block row and column indices.

Reimplemented in libMesh::PetscMatrix< T >.

Definition at line 56 of file sparse_matrix.C.

Referenced by libMesh::SparseMatrix< ValOut >::add_block_matrix().

59 {
60  libmesh_assert_equal_to (dm.m() / brows.size(), dm.n() / bcols.size());
61 
62  const numeric_index_type blocksize = cast_int<numeric_index_type>
63  (dm.m() / brows.size());
64 
65  libmesh_assert_equal_to (dm.m()%blocksize, 0);
66  libmesh_assert_equal_to (dm.n()%blocksize, 0);
67 
68  std::vector<numeric_index_type> rows, cols;
69 
70  rows.reserve(blocksize*brows.size());
71  cols.reserve(blocksize*bcols.size());
72 
73  for (std::size_t ib=0; ib<brows.size(); ib++)
74  {
75  numeric_index_type i=brows[ib]*blocksize;
76 
77  for (unsigned int v=0; v<blocksize; v++)
78  rows.push_back(i++);
79  }
80 
81  for (std::size_t jb=0; jb<bcols.size(); jb++)
82  {
83  numeric_index_type j=bcols[jb]*blocksize;
84 
85  for (unsigned int v=0; v<blocksize; v++)
86  cols.push_back(j++);
87  }
88 
89  this->add_matrix (dm, rows, cols);
90 }
virtual void add_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols)=0
dof_id_type numeric_index_type
Definition: id_types.h:92

◆ add_block_matrix() [2/2]

template<typename T>
virtual void libMesh::SparseMatrix< T >::add_block_matrix ( const DenseMatrix< T > &  dm,
const std::vector< numeric_index_type > &  dof_indices 
)
inlinevirtualinherited

Same as add_block_matrix(), but assumes the row and column maps are the same. Thus the matrix dm must be square.

Reimplemented in libMesh::PetscMatrix< T >.

Definition at line 253 of file sparse_matrix.h.

255  { this->add_block_matrix (dm, dof_indices, dof_indices); }
virtual void add_block_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &brows, const std::vector< numeric_index_type > &bcols)
Definition: sparse_matrix.C:56

◆ add_matrix() [1/2]

template<typename T >
void libMesh::EigenSparseMatrix< T >::add_matrix ( const DenseMatrix< T > &  dm,
const std::vector< numeric_index_type > &  rows,
const std::vector< numeric_index_type > &  cols 
)
overridevirtual

Add the full matrix dm to the SparseMatrix. This is useful for adding an element matrix at assembly time.

Implements libMesh::SparseMatrix< T >.

Definition at line 121 of file eigen_sparse_matrix.C.

References libMesh::initialized(), libMesh::DenseMatrixBase< T >::m(), and libMesh::DenseMatrixBase< T >::n().

125 {
126  libmesh_assert (this->initialized());
127  unsigned int n_rows = cast_int<unsigned int>(rows.size());
128  unsigned int n_cols = cast_int<unsigned int>(cols.size());
129  libmesh_assert_equal_to (dm.m(), n_rows);
130  libmesh_assert_equal_to (dm.n(), n_cols);
131 
132 
133  for (unsigned int i=0; i<n_rows; i++)
134  for (unsigned int j=0; j<n_cols; j++)
135  this->add(rows[i],cols[j],dm(i,j));
136 }
virtual bool initialized() const
virtual void add(const numeric_index_type i, const numeric_index_type j, const T value) override

◆ add_matrix() [2/2]

template<typename T >
void libMesh::EigenSparseMatrix< T >::add_matrix ( const DenseMatrix< T > &  dm,
const std::vector< numeric_index_type > &  dof_indices 
)
overridevirtual

Same as add_matrix, but assumes the row and column maps are the same. Thus the matrix dm must be square.

Implements libMesh::SparseMatrix< T >.

Definition at line 253 of file eigen_sparse_matrix.C.

255 {
256  this->add_matrix (dm, dof_indices, dof_indices);
257 }
virtual void add_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols) override

◆ attach_dof_map()

template<typename T>
void libMesh::SparseMatrix< T >::attach_dof_map ( const DofMap dof_map)
inlineinherited

Get a pointer to the DofMap to use.

Definition at line 109 of file sparse_matrix.h.

Referenced by libMesh::__libmesh_tao_hessian(), and libMesh::DofMap::attach_matrix().

110  { _dof_map = &dof_map; }
DofMap const * _dof_map

◆ build()

template<typename T >
std::unique_ptr< SparseMatrix< T > > libMesh::SparseMatrix< T >::build ( const Parallel::Communicator comm,
const SolverPackage  solver_package = libMesh::default_solver_package() 
)
staticinherited

Builds a SparseMatrix<T> using the linear solver package specified by solver_package

Definition at line 130 of file sparse_matrix.C.

Referenced by libMesh::ImplicitSystem::add_matrix(), libMesh::EigenSystem::init_data(), and libMesh::EigenSystem::init_matrices().

132 {
133  // Avoid unused parameter warnings when no solver packages are enabled.
135 
136  // Build the appropriate vector
137  switch (solver_package)
138  {
139 
140 #ifdef LIBMESH_HAVE_LASPACK
141  case LASPACK_SOLVERS:
142  return libmesh_make_unique<LaspackMatrix<T>>(comm);
143 #endif
144 
145 
146 #ifdef LIBMESH_HAVE_PETSC
147  case PETSC_SOLVERS:
148  return libmesh_make_unique<PetscMatrix<T>>(comm);
149 #endif
150 
151 
152 #ifdef LIBMESH_TRILINOS_HAVE_EPETRA
153  case TRILINOS_SOLVERS:
154  return libmesh_make_unique<EpetraMatrix<T>>(comm);
155 #endif
156 
157 
158 #ifdef LIBMESH_HAVE_EIGEN
159  case EIGEN_SOLVERS:
160  return libmesh_make_unique<EigenSparseMatrix<T>>(comm);
161 #endif
162 
163  default:
164  libmesh_error_msg("ERROR: Unrecognized solver package: " << solver_package);
165  }
166 }
EIGEN_SOLVERS
Definition: libmesh.C:246
TRILINOS_SOLVERS
Definition: libmesh.C:244
const Parallel::Communicator & comm() const
void libmesh_ignore(const Args &...)
LASPACK_SOLVERS
Definition: libmesh.C:248

◆ clear()

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

Restores the SparseMatrix<T> to a pristine state.

Implements libMesh::SparseMatrix< T >.

Definition at line 170 of file eigen_sparse_matrix.C.

References libMesh::libMeshPrivateData::_is_initialized.

171 {
172  _mat.resize(0,0);
173 
174  _closed = false;
175  this->_is_initialized = false;
176 }

◆ close()

template<typename T>
virtual void libMesh::EigenSparseMatrix< T >::close ( )
inlineoverridevirtual

Calls the SparseMatrix's internal assembly routines, ensuring that the values are consistent across processors.

Implements libMesh::SparseMatrix< T >.

Definition at line 100 of file eigen_sparse_matrix.h.

References libMesh::EigenSparseMatrix< T >::_closed.

Referenced by libMesh::EigenSparseLinearSolver< T >::solve().

100 { this->_closed = true; }

◆ closed()

template<typename T>
virtual bool libMesh::EigenSparseMatrix< T >::closed ( ) const
inlineoverridevirtual
Returns
true if the matrix has been assembled.

Implements libMesh::SparseMatrix< T >.

Definition at line 134 of file eigen_sparse_matrix.h.

References libMesh::EigenSparseMatrix< T >::_closed.

◆ 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

◆ create_submatrix()

template<typename T>
virtual void libMesh::SparseMatrix< T >::create_submatrix ( SparseMatrix< T > &  submatrix,
const std::vector< numeric_index_type > &  rows,
const std::vector< numeric_index_type > &  cols 
) const
inlinevirtualinherited

This function creates a matrix called "submatrix" which is defined by the row and column indices given in the "rows" and "cols" entries. Currently this operation is only defined for the PetscMatrix type.

Definition at line 354 of file sparse_matrix.h.

357  {
358  this->_get_submatrix(submatrix,
359  rows,
360  cols,
361  false); // false means DO NOT REUSE submatrix
362  }
virtual void _get_submatrix(SparseMatrix< T > &, const std::vector< numeric_index_type > &, const std::vector< numeric_index_type > &, const bool) const

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

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

◆ flush()

template<typename T>
virtual void libMesh::SparseMatrix< T >::flush ( )
inlinevirtualinherited

For PETSc matrix , this function is similar to close but without shrinking memory. This is useful when we want to switch between ADD_VALUES and INSERT_VALUES. close should be called before using the matrix.

Reimplemented in libMesh::PetscMatrix< T >.

Definition at line 182 of file sparse_matrix.h.

182 { close(); }
virtual void close()=0

◆ get_diagonal()

template<typename T >
void libMesh::EigenSparseMatrix< T >::get_diagonal ( NumericVector< T > &  dest) const
overridevirtual

Copies the diagonal part of the matrix into dest.

Implements libMesh::SparseMatrix< T >.

Definition at line 141 of file eigen_sparse_matrix.C.

References libMesh::EigenSparseVector< T >::_vec.

142 {
143  EigenSparseVector<T> & dest = cast_ref<EigenSparseVector<T> &>(dest_in);
144 
145  dest._vec = _mat.diagonal();
146 }

◆ 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_transpose()

template<typename T >
void libMesh::EigenSparseMatrix< T >::get_transpose ( SparseMatrix< T > &  dest) const
overridevirtual

Copies the transpose of the matrix into dest, which may be *this.

Implements libMesh::SparseMatrix< T >.

Definition at line 151 of file eigen_sparse_matrix.C.

References libMesh::EigenSparseMatrix< T >::_mat.

152 {
153  EigenSparseMatrix<T> & dest = cast_ref<EigenSparseMatrix<T> &>(dest_in);
154 
155  dest._mat = _mat.transpose();
156 }

◆ 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() [1/2]

template<typename T>
virtual void libMesh::EigenSparseMatrix< T >::init ( const numeric_index_type  m,
const numeric_index_type  n,
const numeric_index_type  m_l,
const numeric_index_type  n_l,
const numeric_index_type  nnz = 30,
const numeric_index_type  noz = 10,
const numeric_index_type  blocksize = 1 
)
overridevirtual

Initialize SparseMatrix with the specified sizes.

Parameters
mThe global number of rows.
nThe global number of columns.
m_lThe local number of rows.
n_lThe local number of columns.
nnzThe number of on-diagonal nonzeros per row (defaults to 30).
nozThe number of off-diagonal nonzeros per row (defaults to 10).
blocksizeOptional value indicating dense coupled blocks for systems with multiple variables all of the same type.

Implements libMesh::SparseMatrix< T >.

◆ init() [2/2]

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

Initialize this matrix using the sparsity structure computed by dof_map.

Implements libMesh::SparseMatrix< T >.

Definition at line 63 of file eigen_sparse_matrix.C.

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

64 {
65  // Ignore calls on initialized objects
66  if (this->initialized())
67  return;
68 
69  // We need the DofMap for this!
70  libmesh_assert(this->_dof_map);
71 
72  // Clear initialized matrices
73  if (this->initialized())
74  this->clear();
75 
76  const numeric_index_type n_rows = this->_dof_map->n_dofs();
77  const numeric_index_type n_cols = n_rows;
78 
79 #ifndef NDEBUG
80  // The following variables are only used for assertions,
81  // so avoid declaring them when asserts are inactive.
82  const numeric_index_type n_l = this->_dof_map->n_dofs_on_processor(0);
83  const numeric_index_type m_l = n_l;
84 #endif
85 
86  // Laspack Matrices only work for uniprocessor cases
87  libmesh_assert_equal_to (n_rows, n_cols);
88  libmesh_assert_equal_to (m_l, n_rows);
89  libmesh_assert_equal_to (n_l, n_cols);
90 
91  const std::vector<numeric_index_type> & n_nz = this->_dof_map->get_n_nz();
92 
93 #ifndef NDEBUG
94  // The following variables are only used for assertions,
95  // so avoid declaring them when asserts are inactive.
96  const std::vector<numeric_index_type> & n_oz = this->_dof_map->get_n_oz();
97 #endif
98 
99  // Make sure the sparsity pattern isn't empty
100  libmesh_assert_equal_to (n_nz.size(), n_l);
101  libmesh_assert_equal_to (n_oz.size(), n_l);
102 
103  if (n_rows==0)
104  {
105  _mat.resize(0,0);
106  return;
107  }
108 
109  _mat.resize(n_rows,n_cols);
110  _mat.reserve(n_nz);
111 
112  this->_is_initialized = true;
113 
114  libmesh_assert_equal_to (n_rows, this->m());
115  libmesh_assert_equal_to (n_cols, this->n());
116 }
virtual bool initialized() const
const std::vector< dof_id_type > & get_n_oz() const
Definition: dof_map.h:472
virtual numeric_index_type m() const override
dof_id_type n_dofs_on_processor(const processor_id_type proc) const
Definition: dof_map.h:590
virtual numeric_index_type n() const override
dof_id_type n_dofs() const
Definition: dof_map.h:574
dof_id_type numeric_index_type
Definition: id_types.h:92
DofMap const * _dof_map
const std::vector< dof_id_type > & get_n_nz() const
Definition: dof_map.h:459
virtual void clear() override

◆ initialized()

template<typename T>
virtual bool libMesh::SparseMatrix< T >::initialized ( ) const
inlinevirtualinherited
Returns
true if the matrix has been initialized, false otherwise.

Definition at line 104 of file sparse_matrix.h.

Referenced by libMesh::PetscMatrix< T >::_get_submatrix(), libMesh::ImplicitSystem::assemble(), and libMesh::ImplicitSystem::init_matrices().

104 { return _is_initialized; }

◆ l1_norm()

template<typename T >
Real libMesh::EigenSparseMatrix< T >::l1_norm ( ) const
overridevirtual
Returns
The $ \ell_1 $-norm of the matrix, that is the max column sum: $ |M|_1 = \max_{j} \sum_{i} |M_{ij}| $

This is the natural matrix norm that is compatible with the $ \ell_1 $-norm for vectors, i.e. $ |Mv|_1 \leq |M|_1 |v|_1 $. (cf. Haemmerlin-Hoffmann : Numerische Mathematik)

Implements libMesh::SparseMatrix< T >.

Definition at line 291 of file eigen_sparse_matrix.C.

References std::abs().

292 {
293  // There does not seem to be a straightforward way to iterate over
294  // the columns of an EigenSparseMatrix. So we use some extra
295  // storage and keep track of the column sums while going over the
296  // row entries...
297  std::vector<Real> abs_col_sums(this->n());
298 
299  // For a row-major Eigen SparseMatrix like we're using, the
300  // InnerIterator iterates over the non-zero entries of rows.
301  for (unsigned row=0; row<this->m(); ++row)
302  {
303  EigenSM::InnerIterator it(_mat, row);
304  for (; it; ++it)
305  abs_col_sums[it.col()] += std::abs(it.value());
306  }
307 
308  return *(std::max_element(abs_col_sums.begin(), abs_col_sums.end()));
309 }
double abs(double a)
virtual numeric_index_type m() const override
virtual numeric_index_type n() const override

◆ linfty_norm()

template<typename T >
Real libMesh::EigenSparseMatrix< T >::linfty_norm ( ) const
overridevirtual
Returns
The $ \ell_{\infty} $-norm of the matrix, that is the max row sum:

$ |M|_{\infty} = \max_{i} \sum_{j} |M_{ij}| $

This is the natural matrix norm that is compatible to the $ \ell_{\infty} $-norm of vectors, i.e. $ |Mv|_{\infty} \leq |M|_{\infty} |v|_{\infty} $. (cf. Haemmerlin-Hoffmann : Numerische Mathematik)

Implements libMesh::SparseMatrix< T >.

Definition at line 314 of file eigen_sparse_matrix.C.

References std::abs(), std::max(), and libMesh::Real.

315 {
316  Real max_abs_row_sum = 0.;
317 
318  // For a row-major Eigen SparseMatrix like we're using, the
319  // InnerIterator iterates over the non-zero entries of rows.
320  for (unsigned row=0; row<this->m(); ++row)
321  {
322  Real current_abs_row_sum = 0.;
323  EigenSM::InnerIterator it(_mat, row);
324  for (; it; ++it)
325  current_abs_row_sum += std::abs(it.value());
326 
327  max_abs_row_sum = std::max(max_abs_row_sum, current_abs_row_sum);
328  }
329 
330  return max_abs_row_sum;
331 }
double abs(double a)
virtual numeric_index_type m() const override
long double max(long double a, double b)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ m()

template<typename T >
numeric_index_type libMesh::EigenSparseMatrix< T >::m ( ) const
overridevirtual
Returns
The row-dimension of the matrix.

Implements libMesh::SparseMatrix< T >.

Definition at line 189 of file eigen_sparse_matrix.C.

References libMesh::initialized().

190 {
191  libmesh_assert (this->initialized());
192 
193  return cast_int<numeric_index_type>(_mat.rows());
194 }
virtual bool initialized() const

◆ n()

template<typename T >
numeric_index_type libMesh::EigenSparseMatrix< T >::n ( ) const
overridevirtual
Returns
The column-dimension of the matrix.

Implements libMesh::SparseMatrix< T >.

Definition at line 199 of file eigen_sparse_matrix.C.

References libMesh::initialized().

200 {
201  libmesh_assert (this->initialized());
202 
203  return cast_int<numeric_index_type>(_mat.cols());
204 }
virtual bool initialized() const

◆ 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

◆ need_full_sparsity_pattern()

template<typename T>
virtual bool libMesh::SparseMatrix< T >::need_full_sparsity_pattern ( ) const
inlinevirtualinherited
Returns
true if this sparse matrix format needs to be fed the graph of the sparse matrix.

This is true for LaspackMatrix, but not PetscMatrix. In the case where the full graph is not required, we can efficiently approximate it to provide a good estimate of the required size of the sparse matrix.

Reimplemented in libMesh::EpetraMatrix< T >, and libMesh::LaspackMatrix< T >.

Definition at line 121 of file sparse_matrix.h.

Referenced by libMesh::DofMap::attach_matrix().

122  { return false; }

◆ operator()()

template<typename T >
T libMesh::EigenSparseMatrix< T >::operator() ( const numeric_index_type  i,
const numeric_index_type  j 
) const
overridevirtual
Returns
A copy of matrix entry (i,j).
Note
This may be an expensive operation, and you should always be careful where you call this function.

Implements libMesh::SparseMatrix< T >.

Definition at line 278 of file eigen_sparse_matrix.C.

References libMesh::initialized().

280 {
281  libmesh_assert (this->initialized());
282  libmesh_assert_less (i, this->m());
283  libmesh_assert_less (j, this->n());
284 
285  return _mat.coeff(i,j);
286 }
virtual bool initialized() const
virtual numeric_index_type m() const override
virtual numeric_index_type n() const override

◆ operator=() [1/2]

template<typename T>
EigenSparseMatrix& libMesh::EigenSparseMatrix< T >::operator= ( const EigenSparseMatrix< T > &  )
default

◆ operator=() [2/2]

template<typename T>
EigenSparseMatrix& libMesh::EigenSparseMatrix< T >::operator= ( EigenSparseMatrix< T > &&  )
default

◆ print() [1/2]

template<>
void libMesh::SparseMatrix< Complex >::print ( std::ostream &  os,
const bool  sparse 
) const
inherited

Definition at line 96 of file sparse_matrix.C.

97 {
98  // std::complex<>::operator<<() is defined, but use this form
99 
100  if (sparse)
101  {
102  libmesh_not_implemented();
103  }
104 
105  os << "Real part:" << std::endl;
106  for (numeric_index_type i=0; i<this->m(); i++)
107  {
108  for (numeric_index_type j=0; j<this->n(); j++)
109  os << std::setw(8) << (*this)(i,j).real() << " ";
110  os << std::endl;
111  }
112 
113  os << std::endl << "Imaginary part:" << std::endl;
114  for (numeric_index_type i=0; i<this->m(); i++)
115  {
116  for (numeric_index_type j=0; j<this->n(); j++)
117  os << std::setw(8) << (*this)(i,j).imag() << " ";
118  os << std::endl;
119  }
120 }
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual numeric_index_type m() const =0
virtual numeric_index_type n() const =0

◆ print() [2/2]

template<typename T >
void libMesh::SparseMatrix< T >::print ( std::ostream &  os = libMesh::out,
const bool  sparse = false 
) const
inherited

Print the contents of the matrix to the screen in a uniform style, regardless of matrix/solver package being used.

Definition at line 200 of file sparse_matrix.C.

Referenced by libMesh::EigenSparseMatrix< T >::print_personal(), and libMesh::LaspackMatrix< T >::print_personal().

201 {
202  parallel_object_only();
203 
204  libmesh_assert (this->initialized());
205 
206  if (!this->_dof_map)
207  libmesh_error_msg("Error! Trying to print a matrix with no dof_map set!");
208 
209  // We'll print the matrix from processor 0 to make sure
210  // it's serialized properly
211  if (this->processor_id() == 0)
212  {
213  libmesh_assert_equal_to (this->_dof_map->first_dof(), 0);
214  for (numeric_index_type i=this->_dof_map->first_dof();
215  i!=this->_dof_map->end_dof(); ++i)
216  {
217  if (sparse)
218  {
219  for (numeric_index_type j=0; j<this->n(); j++)
220  {
221  T c = (*this)(i,j);
222  if (c != static_cast<T>(0.0))
223  {
224  os << i << " " << j << " " << c << std::endl;
225  }
226  }
227  }
228  else
229  {
230  for (numeric_index_type j=0; j<this->n(); j++)
231  os << (*this)(i,j) << " ";
232  os << std::endl;
233  }
234  }
235 
236  std::vector<numeric_index_type> ibuf, jbuf;
237  std::vector<T> cbuf;
238  numeric_index_type currenti = this->_dof_map->end_dof();
239  for (processor_id_type p=1; p < this->n_processors(); ++p)
240  {
241  this->comm().receive(p, ibuf);
242  this->comm().receive(p, jbuf);
243  this->comm().receive(p, cbuf);
244  libmesh_assert_equal_to (ibuf.size(), jbuf.size());
245  libmesh_assert_equal_to (ibuf.size(), cbuf.size());
246 
247  if (ibuf.empty())
248  continue;
249  libmesh_assert_greater_equal (ibuf.front(), currenti);
250  libmesh_assert_greater_equal (ibuf.back(), ibuf.front());
251 
252  std::size_t currentb = 0;
253  for (;currenti <= ibuf.back(); ++currenti)
254  {
255  if (sparse)
256  {
257  for (numeric_index_type j=0; j<this->n(); j++)
258  {
259  if (currentb < ibuf.size() &&
260  ibuf[currentb] == currenti &&
261  jbuf[currentb] == j)
262  {
263  os << currenti << " " << j << " " << cbuf[currentb] << std::endl;
264  currentb++;
265  }
266  }
267  }
268  else
269  {
270  for (numeric_index_type j=0; j<this->n(); j++)
271  {
272  if (currentb < ibuf.size() &&
273  ibuf[currentb] == currenti &&
274  jbuf[currentb] == j)
275  {
276  os << cbuf[currentb] << " ";
277  currentb++;
278  }
279  else
280  os << static_cast<T>(0.0) << " ";
281  }
282  os << std::endl;
283  }
284  }
285  }
286  if (!sparse)
287  {
288  for (; currenti != this->m(); ++currenti)
289  {
290  for (numeric_index_type j=0; j<this->n(); j++)
291  os << static_cast<T>(0.0) << " ";
292  os << std::endl;
293  }
294  }
295  }
296  else
297  {
298  std::vector<numeric_index_type> ibuf, jbuf;
299  std::vector<T> cbuf;
300 
301  // We'll assume each processor has access to entire
302  // matrix rows, so (*this)(i,j) is valid if i is a local index.
303  for (numeric_index_type i=this->_dof_map->first_dof();
304  i!=this->_dof_map->end_dof(); ++i)
305  {
306  for (numeric_index_type j=0; j<this->n(); j++)
307  {
308  T c = (*this)(i,j);
309  if (c != static_cast<T>(0.0))
310  {
311  ibuf.push_back(i);
312  jbuf.push_back(j);
313  cbuf.push_back(c);
314  }
315  }
316  }
317  this->comm().send(0,ibuf);
318  this->comm().send(0,jbuf);
319  this->comm().send(0,cbuf);
320  }
321 }
virtual bool initialized() const
void send(const unsigned int dest_processor_id, const T &buf, const MessageTag &tag=no_tag) const
uint8_t processor_id_type
Definition: id_types.h:99
const Parallel::Communicator & comm() const
processor_id_type n_processors() const
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual numeric_index_type m() const =0
DofMap const * _dof_map
Status receive(const unsigned int dest_processor_id, T &buf, const MessageTag &tag=any_tag) const
dof_id_type first_dof(const processor_id_type proc) const
Definition: dof_map.h:599
dof_id_type end_dof(const processor_id_type proc) const
Definition: dof_map.h:641
processor_id_type processor_id() const
virtual numeric_index_type n() const =0

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

◆ print_matlab()

template<typename T>
virtual void libMesh::SparseMatrix< T >::print_matlab ( const std::string &  = "") const
inlinevirtualinherited

Print the contents of the matrix in Matlab's sparse matrix format. Optionally prints the matrix to the file named name. If name is not specified it is dumped to the screen.

Reimplemented in libMesh::PetscMatrix< T >.

Definition at line 344 of file sparse_matrix.h.

345  {
346  libmesh_not_implemented();
347  }

◆ print_personal()

template<typename T>
virtual void libMesh::EigenSparseMatrix< T >::print_personal ( std::ostream &  os = libMesh::out) const
inlineoverridevirtual

Print the contents of the matrix to the screen in a package-personalized style, if available.

Implements libMesh::SparseMatrix< T >.

Definition at line 136 of file eigen_sparse_matrix.h.

References libMesh::SparseMatrix< T >::print().

136 { this->print(os); }
void print(std::ostream &os=libMesh::out, const bool sparse=false) const

◆ 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

◆ reinit_submatrix()

template<typename T>
virtual void libMesh::SparseMatrix< T >::reinit_submatrix ( SparseMatrix< T > &  submatrix,
const std::vector< numeric_index_type > &  rows,
const std::vector< numeric_index_type > &  cols 
) const
inlinevirtualinherited

This function is similar to the one above, but it allows you to reuse the existing sparsity pattern of "submatrix" instead of reallocating it again. This should hopefully be more efficient if you are frequently extracting submatrices of the same size.

Definition at line 370 of file sparse_matrix.h.

373  {
374  this->_get_submatrix(submatrix,
375  rows,
376  cols,
377  true); // true means REUSE submatrix
378  }
virtual void _get_submatrix(SparseMatrix< T > &, const std::vector< numeric_index_type > &, const std::vector< numeric_index_type > &, const bool) const

◆ row_start()

template<typename T >
numeric_index_type libMesh::EigenSparseMatrix< T >::row_start ( ) const
overridevirtual
Returns
The index of the first matrix row stored on this processor.

Implements libMesh::SparseMatrix< T >.

Definition at line 209 of file eigen_sparse_matrix.C.

210 {
211  return 0;
212 }

◆ row_stop()

template<typename T >
numeric_index_type libMesh::EigenSparseMatrix< T >::row_stop ( ) const
overridevirtual
Returns
The index of the last matrix row (+1) stored on this processor.

Implements libMesh::SparseMatrix< T >.

Definition at line 217 of file eigen_sparse_matrix.C.

218 {
219  return this->m();
220 }
virtual numeric_index_type m() const override

◆ set()

template<typename T >
void libMesh::EigenSparseMatrix< T >::set ( const numeric_index_type  i,
const numeric_index_type  j,
const T  value 
)
overridevirtual

Set the element (i,j) to value. Throws an error if the entry does not exist. Zero values can be "stored" in non-existent fields.

Implements libMesh::SparseMatrix< T >.

Definition at line 225 of file eigen_sparse_matrix.C.

References libMesh::initialized(), and value.

228 {
229  libmesh_assert (this->initialized());
230  libmesh_assert_less (i, this->m());
231  libmesh_assert_less (j, this->n());
232 
233  _mat.coeffRef(i,j) = value;
234 }
virtual bool initialized() const
virtual numeric_index_type m() const override
virtual numeric_index_type n() const override
static const bool value
Definition: xdr_io.C:109

◆ update_sparsity_pattern()

template<typename T>
virtual void libMesh::SparseMatrix< T >::update_sparsity_pattern ( const SparsityPattern::Graph )
inlinevirtualinherited

Updates the matrix sparsity pattern. When your SparseMatrix<T> implementation does not need this data, simply do not override this method.

Reimplemented in libMesh::EpetraMatrix< T >, and libMesh::LaspackMatrix< T >.

Definition at line 129 of file sparse_matrix.h.

Referenced by libMesh::DofMap::attach_matrix().

129 {}

◆ vector_mult()

template<typename T>
void libMesh::SparseMatrix< T >::vector_mult ( NumericVector< T > &  dest,
const NumericVector< T > &  arg 
) const
inherited

Multiplies the matrix by the NumericVector arg and stores the result in NumericVector dest.

Definition at line 170 of file sparse_matrix.C.

Referenced by libMesh::ImplicitSystem::qoi_parameter_hessian(), and libMesh::ImplicitSystem::qoi_parameter_hessian_vector_product().

172 {
173  dest.zero();
174  this->vector_mult_add(dest,arg);
175 }
void vector_mult_add(NumericVector< T > &dest, const NumericVector< T > &arg) const

◆ vector_mult_add()

template<typename T>
void libMesh::SparseMatrix< T >::vector_mult_add ( NumericVector< T > &  dest,
const NumericVector< T > &  arg 
) const
inherited

Multiplies the matrix by the NumericVector arg and adds the result to the NumericVector dest.

Definition at line 180 of file sparse_matrix.C.

Referenced by libMesh::ImplicitSystem::weighted_sensitivity_adjoint_solve().

182 {
183  /* This functionality is actually implemented in the \p
184  NumericVector class. */
185  dest.add_vector(arg,*this);
186 }

◆ zero()

template<typename T >
void libMesh::EigenSparseMatrix< T >::zero ( )
overridevirtual

Set all entries to 0.

Implements libMesh::SparseMatrix< T >.

Definition at line 181 of file eigen_sparse_matrix.C.

182 {
183  _mat.setZero();
184 }

◆ zero_rows()

template<typename T>
void libMesh::SparseMatrix< T >::zero_rows ( std::vector< numeric_index_type > &  rows,
diag_value = 0.0 
)
virtualinherited

Sets all row entries to 0 then puts diag_value in the diagonal entry.

Reimplemented in libMesh::PetscMatrix< T >.

Definition at line 191 of file sparse_matrix.C.

192 {
193  /* This functionality isn't implemented or stubbed in every subclass yet */
194  libmesh_not_implemented();
195 }

Friends And Related Function Documentation

◆ EigenSparseLinearSolver< T >

template<typename T>
friend class EigenSparseLinearSolver< T >
friend

Definition at line 158 of file eigen_sparse_matrix.h.

◆ EigenSparseVector< T >

template<typename T>
friend class EigenSparseVector< T >
friend

Make other Eigen datatypes friends

Definition at line 157 of file eigen_sparse_matrix.h.

Member Data Documentation

◆ _closed

template<typename T>
bool libMesh::EigenSparseMatrix< T >::_closed
private

Flag indicating if the matrix has been closed yet.

Definition at line 152 of file eigen_sparse_matrix.h.

Referenced by libMesh::EigenSparseMatrix< T >::close(), and libMesh::EigenSparseMatrix< T >::closed().

◆ _communicator

◆ _counts

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

◆ _dof_map

template<typename T>
DofMap const* libMesh::SparseMatrix< T >::_dof_map
protectedinherited

The DofMap object associated with this object.

Definition at line 425 of file sparse_matrix.h.

Referenced by libMesh::SparseMatrix< ValOut >::attach_dof_map().

◆ _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::SparseMatrix< T >::_is_initialized
protectedinherited

◆ _mat

template<typename T>
DataType libMesh::EigenSparseMatrix< T >::_mat
private

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


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