libMesh::LaspackMatrix< T > Class Template Reference

#include <laspack_matrix.h>

Inheritance diagram for libMesh::LaspackMatrix< T >:

Public Member Functions

 LaspackMatrix (const Parallel::Communicator &comm LIBMESH_CAN_DEFAULT_TO_COMMWORLD)
 
 ~LaspackMatrix ()
 
virtual bool need_full_sparsity_pattern () const libmesh_override
 
virtual void update_sparsity_pattern (const SparsityPattern::Graph &) libmesh_override
 
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) libmesh_override
 
virtual void init () libmesh_override
 
virtual void clear () libmesh_override
 
virtual void zero () libmesh_override
 
virtual void close () const libmesh_override
 
virtual numeric_index_type m () const libmesh_override
 
virtual numeric_index_type n () const libmesh_override
 
virtual numeric_index_type row_start () const libmesh_override
 
virtual numeric_index_type row_stop () const libmesh_override
 
virtual void set (const numeric_index_type i, const numeric_index_type j, const T value) libmesh_override
 
virtual void add (const numeric_index_type i, const numeric_index_type j, const T value) libmesh_override
 
virtual void add_matrix (const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols) libmesh_override
 
virtual void add_matrix (const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &dof_indices) libmesh_override
 
virtual void add (const T a, SparseMatrix< T > &X) libmesh_override
 
virtual T operator() (const numeric_index_type i, const numeric_index_type j) const libmesh_override
 
virtual Real l1_norm () const libmesh_override
 
virtual Real linfty_norm () const libmesh_override
 
virtual bool closed () const libmesh_override
 
virtual void print_personal (std::ostream &os=libMesh::out) const libmesh_override
 
virtual void get_diagonal (NumericVector< T > &dest) const libmesh_override
 
virtual void get_transpose (SparseMatrix< T > &dest) const libmesh_override
 
virtual bool initialized () const
 
void attach_dof_map (const DofMap &dof_map)
 
virtual void zero_rows (std::vector< numeric_index_type > &rows, T diag_value=0.0)
 
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 UniquePtr< 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 Member Functions

numeric_index_type pos (const numeric_index_type i, const numeric_index_type j) const
 

Private Attributes

QMatrix _QMat
 
std::vector< numeric_index_type_csr
 
std::vector< std::vector< numeric_index_type >::const_iterator > _row_start
 
bool _closed
 

Friends

class LaspackVector< T >
 
class LaspackLinearSolver< T >
 

Detailed Description

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

The LaspackMatrix class wraps a QMatrix object from the Laspack library. Currently Laspack only supports real datatypes, so this class is a full specialization of SparseMatrix<T> with T = Real.

Author
Benjamin S. Kirk
Date
2003

Definition at line 55 of file laspack_matrix.h.

Member Typedef Documentation

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 110 of file reference_counter.h.

Constructor & Destructor Documentation

template<typename T>
libMesh::LaspackMatrix< T >::LaspackMatrix ( const Parallel::Communicator &comm  LIBMESH_CAN_DEFAULT_TO_COMMWORLD)

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 236 of file laspack_matrix.C.

236  :
237  SparseMatrix<T>(comm),
238  _closed (false)
239 {
240 }
const Parallel::Communicator & comm() const
template<typename T >
libMesh::LaspackMatrix< T >::~LaspackMatrix ( )

Destructor. Free all memory, but do not release the memory of the sparsity structure.

Definition at line 245 of file laspack_matrix.C.

References libMesh::LaspackMatrix< T >::clear().

246 {
247  this->clear ();
248 }
virtual void clear() libmesh_override

Member Function Documentation

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 that this function must be redefined in derived classes for it to work properly!

Reimplemented in libMesh::PetscMatrix< T >.

Definition at line 425 of file sparse_matrix.h.

429  {
430  libmesh_not_implemented();
431  }
template<typename T >
void libMesh::LaspackMatrix< T >::add ( const numeric_index_type  i,
const numeric_index_type  j,
const T  value 
)
virtual

Add value to the element (i,j). Throws an error if the entry does not exist. Still, it is allowed to store zero values in non-existent fields.

Implements libMesh::SparseMatrix< T >.

Definition at line 354 of file laspack_matrix.C.

References libMesh::LaspackMatrix< T >::_QMat, libMesh::LaspackMatrix< T >::_row_start, libMesh::SparseMatrix< T >::initialized(), libMesh::libmesh_assert(), libMesh::LaspackMatrix< T >::m(), libMesh::LaspackMatrix< T >::n(), and libMesh::LaspackMatrix< T >::pos().

Referenced by libMesh::LaspackMatrix< T >::close().

357 {
358  libmesh_assert (this->initialized());
359  libmesh_assert_less (i, this->m());
360  libmesh_assert_less (j, this->n());
361 
362  const numeric_index_type position = this->pos(i,j);
363 
364  // Sanity check
365  libmesh_assert_equal_to (*(_row_start[i]+position), j);
366 
367  Q_AddVal (&_QMat, i+1, position, value);
368 }
virtual numeric_index_type m() const libmesh_override
libmesh_assert(j)
virtual bool initialized() const
dof_id_type numeric_index_type
Definition: id_types.h:92
std::vector< std::vector< numeric_index_type >::const_iterator > _row_start
virtual numeric_index_type n() const libmesh_override
numeric_index_type pos(const numeric_index_type i, const numeric_index_type j) const
template<typename T >
void libMesh::LaspackMatrix< T >::add ( const T  a,
SparseMatrix< T > &  X 
)
virtual

Add a Sparse matrix X, scaled with a, to this, stores the result in this: $\texttt{this} += a*X $. LASPACK does not provide a true axpy for matrices, so a hand-coded version with hopefully acceptable performance is provided.

Implements libMesh::SparseMatrix< T >.

Definition at line 382 of file laspack_matrix.C.

References libMesh::LaspackMatrix< T >::_QMat, libMesh::LaspackMatrix< T >::_row_start, libMesh::SparseMatrix< T >::initialized(), libMesh::libmesh_assert(), libMesh::LaspackMatrix< T >::m(), libMesh::SparseMatrix< T >::m(), libMesh::LaspackMatrix< T >::n(), and libMesh::SparseMatrix< T >::n().

383 {
384  libmesh_assert (this->initialized());
385  libmesh_assert_equal_to (this->m(), X_in.m());
386  libmesh_assert_equal_to (this->n(), X_in.n());
387 
388  LaspackMatrix<T> * X = cast_ptr<LaspackMatrix<T> *> (&X_in);
389  _LPNumber a = static_cast<_LPNumber> (a_in);
390 
391  libmesh_assert(X);
392 
393  // loops taken from LaspackMatrix<T>::zero ()
394 
395  const numeric_index_type n_rows = this->m();
396 
397  for (numeric_index_type row=0; row<n_rows; row++)
398  {
399  const std::vector<numeric_index_type>::const_iterator
400  r_start = _row_start[row];
401 
402  const numeric_index_type len = (_row_start[row+1] - _row_start[row]);
403 
404  // Make sure we agree on the row length
405  libmesh_assert_equal_to (len, Q_GetLen(&_QMat, row+1));
406  // compare matrix sparsity structures
407  libmesh_assert_equal_to (len, Q_GetLen(&(X->_QMat), row+1));
408 
409 
410  for (numeric_index_type l=0; l<len; l++)
411  {
412  const numeric_index_type j = *(r_start + l);
413 
414  // Make sure the data structures are working
415  libmesh_assert_equal_to ((j+1), Q_GetPos (&_QMat, row+1, l));
416 
417  const _LPNumber value = a * Q_GetEl(const_cast<QMatrix*>(&(X->_QMat)), row+1, j+1);
418  Q_AddVal (&_QMat, row+1, l, value);
419  }
420  }
421 }
virtual numeric_index_type m() const libmesh_override
libmesh_assert(j)
virtual bool initialized() const
dof_id_type numeric_index_type
Definition: id_types.h:92
std::vector< std::vector< numeric_index_type >::const_iterator > _row_start
virtual numeric_index_type n() const libmesh_override
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 Sparse matrix. This is useful for adding an element matrix at assembly time. The matrix is assumed blocked, and brow, bcol correspond to the block row, columm indices.

Reimplemented in libMesh::PetscMatrix< T >.

Definition at line 61 of file sparse_matrix.C.

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

64 {
65  libmesh_assert_equal_to (dm.m() / brows.size(), dm.n() / bcols.size());
66 
67  const numeric_index_type blocksize = cast_int<numeric_index_type>
68  (dm.m() / brows.size());
69 
70  libmesh_assert_equal_to (dm.m()%blocksize, 0);
71  libmesh_assert_equal_to (dm.n()%blocksize, 0);
72 
73  std::vector<numeric_index_type> rows, cols;
74 
75  rows.reserve(blocksize*brows.size());
76  cols.reserve(blocksize*bcols.size());
77 
78  for (std::size_t ib=0; ib<brows.size(); ib++)
79  {
80  numeric_index_type i=brows[ib]*blocksize;
81 
82  for (unsigned int v=0; v<blocksize; v++)
83  rows.push_back(i++);
84  }
85 
86  for (std::size_t jb=0; jb<bcols.size(); jb++)
87  {
88  numeric_index_type j=bcols[jb]*blocksize;
89 
90  for (unsigned int v=0; v<blocksize; v++)
91  cols.push_back(j++);
92  }
93 
94  this->add_matrix (dm, rows, cols);
95 }
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
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 258 of file sparse_matrix.h.

260  { 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:61
template<typename T >
void libMesh::LaspackMatrix< T >::add_matrix ( const DenseMatrix< T > &  dm,
const std::vector< numeric_index_type > &  rows,
const std::vector< numeric_index_type > &  cols 
)
virtual

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

Implements libMesh::SparseMatrix< T >.

Definition at line 200 of file laspack_matrix.C.

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

Referenced by libMesh::LaspackMatrix< T >::add_matrix(), and libMesh::LaspackMatrix< T >::close().

204 {
205  libmesh_assert (this->initialized());
206  unsigned int n_rows = cast_int<unsigned int>(rows.size());
207  unsigned int n_cols = cast_int<unsigned int>(cols.size());
208  libmesh_assert_equal_to (dm.m(), n_rows);
209  libmesh_assert_equal_to (dm.n(), n_cols);
210 
211 
212  for (unsigned int i=0; i<n_rows; i++)
213  for (unsigned int j=0; j<n_cols; j++)
214  this->add(rows[i],cols[j],dm(i,j));
215 }
libmesh_assert(j)
virtual bool initialized() const
virtual void add(const numeric_index_type i, const numeric_index_type j, const T value) libmesh_override
template<typename T >
void libMesh::LaspackMatrix< T >::add_matrix ( const DenseMatrix< T > &  dm,
const std::vector< numeric_index_type > &  dof_indices 
)
virtual

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

Implements libMesh::SparseMatrix< T >.

Definition at line 373 of file laspack_matrix.C.

References libMesh::LaspackMatrix< T >::add_matrix().

375 {
376  this->add_matrix (dm, dof_indices, dof_indices);
377 }
virtual void add_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols) libmesh_override
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 112 of file sparse_matrix.h.

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

113  { _dof_map = &dof_map; }
DofMap const * _dof_map
template<typename T >
UniquePtr< 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 135 of file sparse_matrix.C.

References libMesh::ParallelObject::comm(), libMesh::EIGEN_SOLVERS, libMesh::LASPACK_SOLVERS, libMesh::PETSC_SOLVERS, and libMesh::TRILINOS_SOLVERS.

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

137 {
138  // Build the appropriate vector
139  switch (solver_package)
140  {
141 
142 #ifdef LIBMESH_HAVE_LASPACK
143  case LASPACK_SOLVERS:
144  return UniquePtr<SparseMatrix<T> >(new LaspackMatrix<T>(comm));
145 #endif
146 
147 
148 #ifdef LIBMESH_HAVE_PETSC
149  case PETSC_SOLVERS:
150  return UniquePtr<SparseMatrix<T> >(new PetscMatrix<T>(comm));
151 #endif
152 
153 
154 #ifdef LIBMESH_TRILINOS_HAVE_EPETRA
155  case TRILINOS_SOLVERS:
156  return UniquePtr<SparseMatrix<T> >(new EpetraMatrix<T>(comm));
157 #endif
158 
159 
160 #ifdef LIBMESH_HAVE_EIGEN
161  case EIGEN_SOLVERS:
162  return UniquePtr<SparseMatrix<T> >(new EigenSparseMatrix<T>(comm));
163 #endif
164 
165  default:
166  libmesh_error_msg("ERROR: Unrecognized solver package: " << solver_package);
167  }
168 
169  libmesh_error_msg("We'll never get here!");
170  return UniquePtr<SparseMatrix<T> >();
171 }
EIGEN_SOLVERS
Definition: libmesh.C:260
TRILINOS_SOLVERS
Definition: libmesh.C:258
LASPACK_SOLVERS
Definition: libmesh.C:262
const Parallel::Communicator & comm() const
template<typename T >
void libMesh::LaspackMatrix< T >::clear ( )
virtual

Release all memory and return to a state just like after having called the default constructor.

Implements libMesh::SparseMatrix< T >.

Definition at line 253 of file laspack_matrix.C.

References libMesh::LaspackMatrix< T >::_closed, libMesh::LaspackMatrix< T >::_csr, libMesh::SparseMatrix< T >::_is_initialized, libMesh::LaspackMatrix< T >::_QMat, libMesh::LaspackMatrix< T >::_row_start, and libMesh::SparseMatrix< T >::initialized().

Referenced by libMesh::LaspackMatrix< T >::need_full_sparsity_pattern(), and libMesh::LaspackMatrix< T >::~LaspackMatrix().

254 {
255  if (this->initialized())
256  {
257  Q_Destr(&_QMat);
258  }
259 
260  _csr.clear();
261  _row_start.clear();
262  _closed = false;
263  this->_is_initialized = false;
264 }
std::vector< numeric_index_type > _csr
virtual bool initialized() const
std::vector< std::vector< numeric_index_type >::const_iterator > _row_start
template<typename T>
virtual void libMesh::LaspackMatrix< T >::close ( ) const
inlinevirtual
template<typename T>
virtual bool libMesh::LaspackMatrix< T >::closed ( ) const
inlinevirtual

see if Laspack matrix has been closed and fully assembled yet

Implements libMesh::SparseMatrix< T >.

Definition at line 251 of file laspack_matrix.h.

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

251 { return _closed; }
const Parallel::Communicator& libMesh::ParallelObject::comm ( ) const
inlineinherited
Returns
a reference to the Parallel::Communicator object used by this mesh.

Definition at line 87 of file parallel_object.h.

References libMesh::ParallelObject::_communicator.

Referenced by libMesh::__libmesh_petsc_diff_solver_monitor(), libMesh::__libmesh_petsc_diff_solver_residual(), libMesh::__libmesh_petsc_snes_jacobian(), libMesh::__libmesh_petsc_snes_postcheck(), libMesh::__libmesh_petsc_snes_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::ParmetisPartitioner::_do_repartition(), libMesh::UniformRefinementEstimator::_estimate_error(), libMesh::BoundaryInfo::_find_id_maps(), libMesh::PetscLinearSolver< T >::_petsc_shell_matrix_get_diagonal(), libMesh::SlepcEigenSolver< T >::_petsc_shell_matrix_get_diagonal(), libMesh::PetscLinearSolver< T >::_petsc_shell_matrix_mult(), libMesh::SlepcEigenSolver< 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::ImplicitSystem::add_matrix(), libMesh::System::add_vector(), libMesh::EigenSparseLinearSolver< T >::adjoint_solve(), libMesh::UnstructuredMesh::all_second_order(), libMesh::MeshTools::Modification::all_tri(), libMesh::LaplaceMeshSmoother::allgather_graph(), libMesh::FEMSystem::assemble_qoi(), libMesh::MeshCommunication::assign_global_indices(), libMesh::ParmetisPartitioner::assign_partitioning(), libMesh::DofMap::attach_matrix(), libMesh::Parallel::BinSorter< KeyType, IdxType >::binsort(), libMesh::Parallel::Sort< KeyType, IdxType >::binsort(), libMesh::MeshTools::bounding_box(), libMesh::MeshCommunication::broadcast(), libMesh::SparseMatrix< T >::build(), libMesh::MeshTools::Generation::build_extrusion(), libMesh::Parallel::Histogram< KeyType, IdxType >::build_histogram(), libMesh::PetscNonlinearSolver< T >::build_mat_null_space(), libMesh::BoundaryInfo::build_node_list_from_side_list(), libMesh::EquationSystems::build_parallel_solution_vector(), libMesh::MeshBase::cache_elem_dims(), libMesh::System::calculate_norm(), libMesh::DofMap::check_dirichlet_bcid_consistency(), libMesh::DistributedVector< T >::clone(), libMesh::EigenSparseVector< T >::clone(), libMesh::LaspackVector< T >::clone(), libMesh::EpetraVector< T >::clone(), libMesh::PetscVector< T >::clone(), libMesh::EpetraVector< T >::close(), libMesh::Parallel::Sort< KeyType, IdxType >::communicate_bins(), 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::MeshTools::correct_node_proc_ids(), libMesh::MeshRefinement::create_parent_error_vector(), 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::MeshCommunication::gather(), libMesh::MeshCommunication::gather_neighboring_elements(), libMesh::CondensedEigenSystem::get_eigenpair(), libMesh::DofMap::get_info(), libMesh::ImplicitSystem::get_linear_solver(), libMesh::EquationSystems::get_solution(), libMesh::LocationMap< T >::init(), libMesh::PetscDiffSolver::init(), libMesh::TimeSolver::init(), libMesh::TopologyMap::init(), libMesh::TaoOptimizationSolver< T >::init(), libMesh::PetscNonlinearSolver< T >::init(), libMesh::DistributedVector< T >::init(), libMesh::EpetraVector< T >::init(), libMesh::PetscVector< T >::init(), libMesh::SystemSubsetBySubdomain::init(), libMesh::EigenSystem::init_data(), libMesh::EigenSystem::init_matrices(), libMesh::ParmetisPartitioner::initialize(), libMesh::OptimizationSystem::initialize_equality_constraints_storage(), libMesh::OptimizationSystem::initialize_inequality_constraints_storage(), 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::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::DistributedVector< T >::max(), libMesh::FEMSystem::mesh_position_set(), libMesh::MeshSerializer::MeshSerializer(), libMesh::DistributedVector< T >::min(), 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::Partitioner::partition(), libMesh::MetisPartitioner::partition_range(), libMesh::Partitioner::partition_unpartitioned_elements(), libMesh::petsc_auto_fieldsplit(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::MeshBase::prepare_for_use(), libMesh::SparseMatrix< T >::print(), libMesh::MeshTools::processor_bounding_box(), libMesh::Nemesis_IO::read(), libMesh::XdrIO::read(), libMesh::XdrIO::read_serialized_bc_names(), libMesh::XdrIO::read_serialized_bcs_helper(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::XdrIO::read_serialized_nodes(), libMesh::XdrIO::read_serialized_nodesets(), libMesh::XdrIO::read_serialized_subdomain_names(), libMesh::MeshBase::recalculate_n_partitions(), libMesh::MeshCommunication::redistribute(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::MeshCommunication::send_coarse_ghosts(), libMesh::Partitioner::set_node_processor_ids(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::Partitioner::set_parent_processor_ids(), libMesh::LaplaceMeshSmoother::smooth(), libMesh::Parallel::Sort< KeyType, IdxType >::sort(), libMesh::MeshTools::subdomain_bounding_box(), libMesh::MeshBase::subdomain_ids(), libMesh::BoundaryInfo::sync(), libMesh::Parallel::sync_element_data_by_parent_id(), libMesh::Parallel::sync_node_data_by_element_id(), libMesh::MeshRefinement::test_level_one(), libMesh::MeshRefinement::test_unflagged(), libMesh::MeshTools::total_weight(), libMesh::NameBasedIO::write(), libMesh::XdrIO::write(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), libMesh::XdrIO::write_serialized_nodesets(), libMesh::DistributedVector< T >::zero_clone(), libMesh::LaspackVector< T >::zero_clone(), libMesh::EigenSparseVector< T >::zero_clone(), libMesh::EpetraVector< T >::zero_clone(), and libMesh::PetscVector< T >::zero_clone().

88  { return _communicator; }
const Parallel::Communicator & _communicator
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 368 of file sparse_matrix.h.

371  {
372  this->_get_submatrix(submatrix,
373  rows,
374  cols,
375  false); // false means DO NOT REUSE submatrix
376  }
virtual void _get_submatrix(SparseMatrix< T > &, const std::vector< numeric_index_type > &, const std::vector< numeric_index_type > &, const bool) const
void libMesh::ReferenceCounter::disable_print_counter_info ( )
staticinherited
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.

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

101 {
102  _enable_print_counter = true;
103  return;
104 }
template<typename T >
void libMesh::LaspackMatrix< T >::get_diagonal ( NumericVector< T > &  dest) const
virtual

Copies the diagonal part of the matrix into dest.

Implements libMesh::SparseMatrix< T >.

Definition at line 220 of file laspack_matrix.C.

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

221 {
222  libmesh_not_implemented();
223 }
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 (Counts::iterator it = _counts.begin();
59  it != _counts.end(); ++it)
60  {
61  const std::string name(it->first);
62  const unsigned int creations = it->second.first;
63  const unsigned int destructions = it->second.second;
64 
65  oss << "| " << name << " reference count information:\n"
66  << "| Creations: " << creations << '\n'
67  << "| Destructions: " << destructions << '\n';
68  }
69 
70  oss << " ---------------------------------------------------------------------------- \n";
71 
72  return oss.str();
73 
74 #else
75 
76  return "";
77 
78 #endif
79 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:39
template<typename T >
void libMesh::LaspackMatrix< T >::get_transpose ( SparseMatrix< T > &  dest) const
virtual

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

Implements libMesh::SparseMatrix< T >.

Definition at line 228 of file laspack_matrix.C.

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

229 {
230  libmesh_not_implemented();
231 }
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 160 of file reference_counter.h.

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

Referenced by libMesh::ReferenceCounter::n_objects(), and libMesh::ReferenceCountedObject< RBParametrized >::ReferenceCountedObject().

161 {
162  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
163  std::pair<unsigned int, unsigned int> & p = _counts[name];
164 
165  p.first++;
166 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:39
spin_mutex spin_mtx
Definition: threads.C:29
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 173 of file reference_counter.h.

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

Referenced by libMesh::ReferenceCounter::n_objects(), and libMesh::ReferenceCountedObject< RBParametrized >::~ReferenceCountedObject().

174 {
175  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
176  std::pair<unsigned int, unsigned int> & p = _counts[name];
177 
178  p.second++;
179 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:39
spin_mutex spin_mtx
Definition: threads.C:29
template<typename T>
virtual void libMesh::LaspackMatrix< 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 
)
virtual

Initialize a Laspack matrix that is of global dimension $ m \times n $ with local dimensions $ m_l \times n_l $. nnz is the number of on-processor nonzeros per row (defaults to 30). noz is the number of on-processor nonzeros per row (defaults to 30). Optionally supports a block size, which indicates dense coupled blocks for systems with multiple variables all of the same type.

Implements libMesh::SparseMatrix< T >.

template<typename T >
void libMesh::LaspackMatrix< T >::init ( )
virtual

Initialize using sparsity structure computed by dof_map.

Implements libMesh::SparseMatrix< T >.

Definition at line 149 of file laspack_matrix.C.

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

Referenced by libMesh::LaspackMatrix< T >::need_full_sparsity_pattern(), and libMesh::LaspackMatrix< T >::update_sparsity_pattern().

150 {
151  // Ignore calls on initialized objects
152  if (this->initialized())
153  return;
154 
155  // We need the DofMap for this!
156  libmesh_assert(this->_dof_map);
157 
158  // Clear intialized matrices
159  if (this->initialized())
160  this->clear();
161 
162  const numeric_index_type n_rows = this->_dof_map->n_dofs();
163 #ifndef NDEBUG
164  // The following variables are only used for assertions,
165  // so avoid declaring them when asserts are inactive.
166  const numeric_index_type n_cols = n_rows;
167  const numeric_index_type n_l = this->_dof_map->n_dofs_on_processor(0);
168  const numeric_index_type m_l = n_l;
169 #endif
170 
171  // Laspack Matrices only work for uniprocessor cases
172  libmesh_assert_equal_to (n_rows, n_cols);
173  libmesh_assert_equal_to (m_l, n_rows);
174  libmesh_assert_equal_to (n_l, n_cols);
175 
176 #ifndef NDEBUG
177  // The following variables are only used for assertions,
178  // so avoid declaring them when asserts are inactive.
179  const std::vector<numeric_index_type> & n_nz = this->_dof_map->get_n_nz();
180  const std::vector<numeric_index_type> & n_oz = this->_dof_map->get_n_oz();
181 #endif
182 
183  // Make sure the sparsity pattern isn't empty
184  libmesh_assert_equal_to (n_nz.size(), n_l);
185  libmesh_assert_equal_to (n_oz.size(), n_l);
186 
187  if (n_rows==0)
188  return;
189 
190  Q_Constr(&_QMat, const_cast<char *>("Mat"), n_rows, _LPFalse, Rowws, Normal, _LPTrue);
191 
192  this->_is_initialized = true;
193 
194  libmesh_assert_equal_to (n_rows, this->m());
195 }
dof_id_type n_dofs() const
Definition: dof_map.h:506
virtual numeric_index_type m() const libmesh_override
dof_id_type n_dofs_on_processor(const processor_id_type proc) const
Definition: dof_map.h:522
libmesh_assert(j)
virtual bool initialized() const
dof_id_type numeric_index_type
Definition: id_types.h:92
const std::vector< dof_id_type > & get_n_oz() const
Definition: dof_map.h:412
DofMap const * _dof_map
const std::vector< dof_id_type > & get_n_nz() const
Definition: dof_map.h:400
virtual void clear() libmesh_override
template<typename T>
virtual Real libMesh::LaspackMatrix< T >::l1_norm ( ) const
inlinevirtual

Return the l1-norm of the matrix, that is $|M|_1=max_{all columns j}\sum_{all rows i} |M_ij|$, (max. sum of columns). This is the natural matrix norm that is compatible to the l1-norm for vectors, i.e. $|Mv|_1\leq |M|_1 |v|_1$.

Implements libMesh::SparseMatrix< T >.

Definition at line 232 of file laspack_matrix.h.

232 { libmesh_not_implemented(); return 0.; }
template<typename T>
virtual Real libMesh::LaspackMatrix< T >::linfty_norm ( ) const
inlinevirtual

Return the linfty-norm of the matrix, that is $|M|_\infty=max_{all rows i}\sum_{all columns j} |M_ij|$, (max. sum of rows). This is the natural matrix norm that is compatible to the linfty-norm of vectors, i.e. $|Mv|_\infty \leq |M|_\infty |v|_\infty$.

Implements libMesh::SparseMatrix< T >.

Definition at line 245 of file laspack_matrix.h.

245 { libmesh_not_implemented(); return 0.; }
template<typename T >
numeric_index_type libMesh::LaspackMatrix< T >::n ( ) const
virtual
static unsigned int libMesh::ReferenceCounter::n_objects ( )
inlinestaticinherited
processor_id_type libMesh::ParallelObject::n_processors ( ) const
inlineinherited
Returns
the number of processors in the group.

Definition at line 93 of file parallel_object.h.

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

Referenced by libMesh::ParmetisPartitioner::_do_repartition(), libMesh::BoundaryInfo::_find_id_maps(), libMesh::DistributedMesh::add_elem(), libMesh::DistributedMesh::add_node(), libMesh::LaplaceMeshSmoother::allgather_graph(), libMesh::FEMSystem::assembly(), libMesh::ParmetisPartitioner::assign_partitioning(), libMesh::AztecLinearSolver< T >::AztecLinearSolver(), libMesh::MeshCommunication::broadcast(), libMesh::BoundaryInfo::build_node_list_from_side_list(), libMesh::DistributedMesh::clear(), libMesh::Nemesis_IO_Helper::compute_border_node_ids(), libMesh::Nemesis_IO_Helper::construct_nemesis_filename(), libMesh::MeshTools::correct_node_proc_ids(), libMesh::UnstructuredMesh::create_pid_mesh(), libMesh::DofMap::distribute_dofs(), libMesh::DofMap::distribute_local_dofs_node_major(), libMesh::DofMap::distribute_local_dofs_var_major(), libMesh::DistributedMesh::DistributedMesh(), libMesh::EnsightIO::EnsightIO(), libMesh::MeshCommunication::gather(), libMesh::MeshCommunication::gather_neighboring_elements(), libMesh::MeshBase::get_info(), libMesh::EquationSystems::get_solution(), libMesh::DistributedVector< T >::init(), libMesh::SystemSubsetBySubdomain::init(), libMesh::ParmetisPartitioner::initialize(), libMesh::Nemesis_IO_Helper::initialize(), libMesh::DistributedMesh::insert_elem(), 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::MeshBase::n_active_elem_on_proc(), libMesh::MeshBase::n_elem_on_proc(), libMesh::MeshBase::n_nodes_on_proc(), libMesh::SparsityPattern::Build::parallel_sync(), libMesh::Partitioner::partition(), libMesh::MeshBase::partition(), libMesh::Partitioner::partition_unpartitioned_elements(), libMesh::PetscLinearSolver< T >::PetscLinearSolver(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::SparseMatrix< T >::print(), libMesh::MeshTools::processor_bounding_box(), libMesh::NameBasedIO::read(), libMesh::Nemesis_IO::read(), libMesh::XdrIO::read(), libMesh::MeshCommunication::redistribute(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::Partitioner::repartition(), libMesh::MeshCommunication::send_coarse_ghosts(), libMesh::Partitioner::set_node_processor_ids(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::Parallel::Sort< KeyType, IdxType >::sort(), libMesh::DistributedMesh::update_parallel_id_counts(), libMesh::GMVIO::write_binary(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), and libMesh::XdrIO::write_serialized_nodesets().

94  { return cast_int<processor_id_type>(_communicator.size()); }
unsigned int size() const
Definition: parallel.h:679
const Parallel::Communicator & _communicator
template<typename T>
virtual bool libMesh::LaspackMatrix< T >::need_full_sparsity_pattern ( ) const
inlinevirtual
template<typename T >
T libMesh::LaspackMatrix< T >::operator() ( const numeric_index_type  i,
const numeric_index_type  j 
) const
virtual

Return the value of the entry (i,j). This may be an expensive operation, and you should always be careful where you call this function.

Implements libMesh::SparseMatrix< T >.

Definition at line 427 of file laspack_matrix.C.

References libMesh::LaspackMatrix< T >::_QMat, libMesh::SparseMatrix< T >::initialized(), libMesh::libmesh_assert(), libMesh::LaspackMatrix< T >::m(), and libMesh::LaspackMatrix< T >::n().

Referenced by libMesh::LaspackMatrix< T >::close().

429 {
430  libmesh_assert (this->initialized());
431  libmesh_assert_less (i, this->m());
432  libmesh_assert_less (j, this->n());
433 
434  return Q_GetEl (const_cast<QMatrix*>(&_QMat), i+1, j+1);
435 }
virtual numeric_index_type m() const libmesh_override
libmesh_assert(j)
virtual bool initialized() const
virtual numeric_index_type n() const libmesh_override
template<typename T >
numeric_index_type libMesh::LaspackMatrix< T >::pos ( const numeric_index_type  i,
const numeric_index_type  j 
) const
private
Returns
the position in the compressed row storage scheme of the $ (i,j) $ element.

Definition at line 440 of file laspack_matrix.C.

References libMesh::LaspackMatrix< T >::_csr, libMesh::LaspackMatrix< T >::_row_start, libMesh::libmesh_assert(), libMesh::LaspackMatrix< T >::m(), and libMesh::LaspackMatrix< T >::n().

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

442 {
443  libmesh_assert_less (i, this->m());
444  libmesh_assert_less (j, this->n());
445  libmesh_assert_less (i+1, _row_start.size());
446  libmesh_assert (_row_start.back() == _csr.end());
447 
448  // note this requires the _csr to be
449  std::pair<std::vector<numeric_index_type>::const_iterator,
450  std::vector<numeric_index_type>::const_iterator> p =
451  std::equal_range (_row_start[i],
452  _row_start[i+1],
453  j);
454 
455  // Make sure the row contains the element j
456  libmesh_assert (p.first != p.second);
457 
458  // Make sure the values match
459  libmesh_assert (*p.first == j);
460 
461  // Return the position in the compressed row
462  return std::distance (_row_start[i], p.first);
463 }
std::vector< numeric_index_type > _csr
virtual numeric_index_type m() const libmesh_override
libmesh_assert(j)
std::vector< std::vector< numeric_index_type >::const_iterator > _row_start
virtual numeric_index_type n() const libmesh_override
template<>
void libMesh::SparseMatrix< Complex >::print ( std::ostream &  os,
const bool  sparse 
) const
inherited

Definition at line 101 of file sparse_matrix.C.

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

102 {
103  // std::complex<>::operator<<() is defined, but use this form
104 
105  if(sparse)
106  {
107  libmesh_not_implemented();
108  }
109 
110  os << "Real part:" << std::endl;
111  for (numeric_index_type i=0; i<this->m(); i++)
112  {
113  for (numeric_index_type j=0; j<this->n(); j++)
114  os << std::setw(8) << (*this)(i,j).real() << " ";
115  os << std::endl;
116  }
117 
118  os << std::endl << "Imaginary part:" << std::endl;
119  for (numeric_index_type i=0; i<this->m(); i++)
120  {
121  for (numeric_index_type j=0; j<this->n(); j++)
122  os << std::setw(8) << (*this)(i,j).imag() << " ";
123  os << std::endl;
124  }
125 }
virtual numeric_index_type n() const =0
virtual numeric_index_type m() const =0
dof_id_type numeric_index_type
Definition: id_types.h:92
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 205 of file sparse_matrix.C.

References libMesh::SparseMatrix< T >::_dof_map, libMesh::ParallelObject::comm(), libMesh::DofMap::end_dof(), libMesh::DofMap::first_dof(), libMesh::SparseMatrix< T >::initialized(), libMesh::libmesh_assert(), libMesh::SparseMatrix< T >::m(), libMesh::SparseMatrix< T >::n(), libMesh::ParallelObject::n_processors(), libMesh::ParallelObject::processor_id(), libMesh::Parallel::Communicator::receive(), and libMesh::Parallel::Communicator::send().

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

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

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

Definition at line 88 of file reference_counter.C.

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

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

89 {
91 }
static std::string get_info()
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 358 of file sparse_matrix.h.

359  {
360  libmesh_not_implemented();
361  }
template<typename T>
virtual void libMesh::LaspackMatrix< T >::print_personal ( std::ostream &  os = libMesh::out) const
inlinevirtual

Print the contents of the matrix, by default to libMesh::out. Currently identical to print().

Implements libMesh::SparseMatrix< T >.

Definition at line 257 of file laspack_matrix.h.

References libMesh::LaspackMatrix< T >::get_diagonal(), libMesh::LaspackMatrix< T >::get_transpose(), libMesh::LaspackMatrix< T >::pos(), and libMesh::SparseMatrix< T >::print().

257 { this->print(os); }
void print(std::ostream &os=libMesh::out, const bool sparse=false) const
processor_id_type libMesh::ParallelObject::processor_id ( ) const
inlineinherited
Returns
the rank of this processor in the group.

Definition at line 99 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::DistributedMesh::add_elem(), libMesh::BoundaryInfo::add_elements(), libMesh::DofMap::add_neighbors_to_send_list(), libMesh::DistributedMesh::add_node(), libMesh::MeshRefinement::add_node(), libMesh::MeshTools::Modification::all_tri(), libMesh::FEMSystem::assembly(), libMesh::ParmetisPartitioner::assign_partitioning(), libMesh::MeshCommunication::broadcast(), libMesh::EquationSystems::build_discontinuous_solution_vector(), libMesh::Nemesis_IO_Helper::build_element_and_node_maps(), libMesh::ParmetisPartitioner::build_graph(), libMesh::InfElemBuilder::build_inf_elem(), libMesh::BoundaryInfo::build_node_list_from_side_list(), libMesh::DofMap::build_sparsity(), 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::EnsightIO::EnsightIO(), libMesh::MeshFunction::find_element(), libMesh::MeshFunction::find_elements(), libMesh::UnstructuredMesh::find_neighbors(), libMesh::MeshCommunication::gather(), libMesh::MeshCommunication::gather_neighboring_elements(), 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::EquationSystems::get_solution(), libMesh::Nemesis_IO_Helper::get_ss_param_global(), libMesh::SparsityPattern::Build::handle_vi_vj(), libMesh::DistributedVector< T >::init(), libMesh::SystemSubsetBySubdomain::init(), libMesh::ParmetisPartitioner::initialize(), 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::MeshTools::libmesh_assert_parallel_consistent_procids< Elem >(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Node >(), 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::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::WeightedPatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::SparsityPattern::Build::operator()(), libMesh::PatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::SparsityPattern::Build::parallel_sync(), libMesh::MetisPartitioner::partition_range(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::SparseMatrix< T >::print(), libMesh::NumericVector< T >::print_global(), 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::ExodusII_IO_Helper::read_elem_num_map(), libMesh::ExodusII_IO_Helper::read_node_num_map(), libMesh::XdrIO::read_serialized_bc_names(), libMesh::XdrIO::read_serialized_bcs_helper(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::XdrIO::read_serialized_nodes(), libMesh::XdrIO::read_serialized_nodesets(), libMesh::XdrIO::read_serialized_subdomain_names(), libMesh::MeshCommunication::redistribute(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::MeshCommunication::send_coarse_ghosts(), libMesh::Partitioner::set_node_processor_ids(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::LaplaceMeshSmoother::smooth(), libMesh::MeshTools::total_weight(), libMesh::Parallel::Packing< Node * >::unpack(), libMesh::Parallel::Packing< Elem * >::unpack(), libMesh::DistributedMesh::update_parallel_id_counts(), libMesh::MeshTools::weight(), libMesh::NameBasedIO::write(), libMesh::XdrIO::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::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::XdrIO::write_serialized_bc_names(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), libMesh::XdrIO::write_serialized_nodesets(), libMesh::XdrIO::write_serialized_subdomain_names(), libMesh::Nemesis_IO_Helper::write_sidesets(), libMesh::ExodusII_IO_Helper::write_sidesets(), libMesh::ExodusII_IO::write_timestep(), and libMesh::ExodusII_IO_Helper::write_timestep().

100  { return cast_int<processor_id_type>(_communicator.rank()); }
const Parallel::Communicator & _communicator
unsigned int rank() const
Definition: parallel.h:677
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 384 of file sparse_matrix.h.

387  {
388  this->_get_submatrix(submatrix,
389  rows,
390  cols,
391  true); // true means REUSE submatrix
392  }
virtual void _get_submatrix(SparseMatrix< T > &, const std::vector< numeric_index_type > &, const std::vector< numeric_index_type > &, const bool) const
template<typename T >
numeric_index_type libMesh::LaspackMatrix< T >::row_start ( ) const
virtual

return row_start, the index of the first matrix row stored on this processor

Implements libMesh::SparseMatrix< T >.

Definition at line 318 of file laspack_matrix.C.

Referenced by libMesh::LaspackMatrix< T >::close().

319 {
320  return 0;
321 }
template<typename T >
numeric_index_type libMesh::LaspackMatrix< T >::row_stop ( ) const
virtual

return row_stop, the index of the last matrix row (+1) stored on this processor

Implements libMesh::SparseMatrix< T >.

Definition at line 326 of file laspack_matrix.C.

References libMesh::LaspackMatrix< T >::m().

Referenced by libMesh::LaspackMatrix< T >::close().

327 {
328  return this->m();
329 }
virtual numeric_index_type m() const libmesh_override
template<typename T >
void libMesh::LaspackMatrix< T >::set ( const numeric_index_type  i,
const numeric_index_type  j,
const T  value 
)
virtual

Set the element (i,j) to value. Throws an error if the entry does not exist. Still, it is allowed to store zero values in non-existent fields.

Implements libMesh::SparseMatrix< T >.

Definition at line 334 of file laspack_matrix.C.

References libMesh::LaspackMatrix< T >::_QMat, libMesh::LaspackMatrix< T >::_row_start, libMesh::SparseMatrix< T >::initialized(), libMesh::libmesh_assert(), libMesh::LaspackMatrix< T >::m(), libMesh::LaspackMatrix< T >::n(), and libMesh::LaspackMatrix< T >::pos().

337 {
338  libmesh_assert (this->initialized());
339  libmesh_assert_less (i, this->m());
340  libmesh_assert_less (j, this->n());
341 
342  const numeric_index_type position = this->pos(i,j);
343 
344  // Sanity check
345  libmesh_assert_equal_to (*(_row_start[i]+position), j);
346  libmesh_assert_equal_to ((j+1), Q_GetPos (&_QMat, i+1, position));
347 
348  Q_SetEntry (&_QMat, i+1, position, j+1, value);
349 }
virtual numeric_index_type m() const libmesh_override
libmesh_assert(j)
virtual bool initialized() const
dof_id_type numeric_index_type
Definition: id_types.h:92
std::vector< std::vector< numeric_index_type >::const_iterator > _row_start
virtual numeric_index_type n() const libmesh_override
numeric_index_type pos(const numeric_index_type i, const numeric_index_type j) const
template<typename T >
void libMesh::LaspackMatrix< T >::update_sparsity_pattern ( const SparsityPattern::Graph sparsity_pattern)
virtual

Updates the matrix sparsity pattern. This will tell the underlying matrix storage scheme how to map the $ (i,j) $ elements.

Reimplemented from libMesh::SparseMatrix< T >.

Definition at line 39 of file laspack_matrix.C.

References libMesh::libMeshPrivateData::_is_initialized, libMesh::TriangleWrapper::init(), libMesh::LaspackMatrix< T >::init(), libMesh::initialized(), libMesh::libmesh_assert(), and libMesh::libmesh_dbg_var().

Referenced by libMesh::LaspackMatrix< T >::need_full_sparsity_pattern().

40 {
41  // clear data, start over
42  this->clear ();
43 
44  // big trouble if this fails!
45  libmesh_assert(this->_dof_map);
46 
47  const numeric_index_type n_rows = sparsity_pattern.size();
48 
49  // Initialize the _row_start data structure,
50  // allocate storage for the _csr array
51  {
52  std::size_t size = 0;
53 
54  for (numeric_index_type row=0; row<n_rows; row++)
55  size += sparsity_pattern[row].size();
56 
57  _csr.resize (size);
58  _row_start.reserve(n_rows + 1);
59  }
60 
61 
62  // Initize the _csr data structure.
63  {
64  std::vector<numeric_index_type>::iterator position = _csr.begin();
65 
66  _row_start.push_back (position);
67 
68  for (numeric_index_type row=0; row<n_rows; row++)
69  {
70  // insert the row indices
71  for (SparsityPattern::Row::const_iterator col = sparsity_pattern[row].begin();
72  col != sparsity_pattern[row].end(); ++col)
73  {
74  libmesh_assert (position != _csr.end());
75  *position = *col;
76  ++position;
77  }
78 
79  _row_start.push_back (position);
80  }
81  }
82 
83 
84  // Initialize the matrix
85  libmesh_assert (!this->initialized());
86  this->init ();
87  libmesh_assert (this->initialized());
88  //libMesh::out << "n_rows=" << n_rows << std::endl;
89  //libMesh::out << "m()=" << m() << std::endl;
90  libmesh_assert_equal_to (n_rows, this->m());
91 
92  // Tell the matrix about its structure. Initialize it
93  // to zero.
94  for (numeric_index_type i=0; i<n_rows; i++)
95  {
96  const std::vector<numeric_index_type>::const_iterator
97  rs = _row_start[i];
98 
99  const numeric_index_type length = _row_start[i+1] - rs;
100 
101  Q_SetLen (&_QMat, i+1, length);
102 
103  for (numeric_index_type l=0; l<length; l++)
104  {
105  const numeric_index_type j = *(rs+l);
106 
107  // sanity check
108  //libMesh::out << "m()=" << m() << std::endl;
109  //libMesh::out << "(i,j,l) = (" << i
110  // << "," << j
111  // << "," << l
112  // << ")" << std::endl;
113  //libMesh::out << "pos(i,j)=" << pos(i,j)
114  // << std::endl;
115  libmesh_assert_equal_to (this->pos(i,j), l);
116  Q_SetEntry (&_QMat, i+1, l, j+1, 0.);
117  }
118  }
119 
120  // That's it!
121  //libmesh_here();
122 }
std::vector< numeric_index_type > _csr
virtual void init() libmesh_override
virtual numeric_index_type m() const libmesh_override
libmesh_assert(j)
virtual bool initialized() const
dof_id_type numeric_index_type
Definition: id_types.h:92
DofMap const * _dof_map
std::vector< std::vector< numeric_index_type >::const_iterator > _row_start
virtual void clear() libmesh_override
numeric_index_type pos(const numeric_index_type i, const numeric_index_type j) const
template<typename T>
void libMesh::SparseMatrix< T >::vector_mult ( NumericVector< T > &  dest,
const NumericVector< T > &  arg 
) const
inherited

Multiplies the matrix with arg and stores the result in dest.

Definition at line 175 of file sparse_matrix.C.

References libMesh::SparseMatrix< T >::vector_mult_add(), and libMesh::NumericVector< T >::zero().

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

177 {
178  dest.zero();
179  this->vector_mult_add(dest,arg);
180 }
void vector_mult_add(NumericVector< T > &dest, const NumericVector< T > &arg) const
template<typename T>
void libMesh::SparseMatrix< T >::vector_mult_add ( NumericVector< T > &  dest,
const NumericVector< T > &  arg 
) const
inherited

Multiplies the matrix with arg and adds the result to dest.

Definition at line 185 of file sparse_matrix.C.

References libMesh::NumericVector< T >::add_vector().

Referenced by libMesh::SparseMatrix< T >::vector_mult(), and libMesh::ImplicitSystem::weighted_sensitivity_adjoint_solve().

187 {
188  /* This functionality is actually implemented in the \p
189  NumericVector class. */
190  dest.add_vector(arg,*this);
191 }
template<typename T >
void libMesh::LaspackMatrix< T >::zero ( )
virtual

Set all entries to 0.

Implements libMesh::SparseMatrix< T >.

Definition at line 269 of file laspack_matrix.C.

References libMesh::LaspackMatrix< T >::_QMat, libMesh::LaspackMatrix< T >::_row_start, and libMesh::LaspackMatrix< T >::m().

Referenced by libMesh::LaspackMatrix< T >::need_full_sparsity_pattern().

270 {
271  const numeric_index_type n_rows = this->m();
272 
273  for (numeric_index_type row=0; row<n_rows; row++)
274  {
275  const std::vector<numeric_index_type>::const_iterator
276  r_start = _row_start[row];
277 
278  const numeric_index_type len = (_row_start[row+1] - _row_start[row]);
279 
280  // Make sure we agree on the row length
281  libmesh_assert_equal_to (len, Q_GetLen(&_QMat, row+1));
282 
283  for (numeric_index_type l=0; l<len; l++)
284  {
285  const numeric_index_type j = *(r_start + l);
286 
287  // Make sure the data structures are working
288  libmesh_assert_equal_to ((j+1), Q_GetPos (&_QMat, row+1, l));
289 
290  Q_SetEntry (&_QMat, row+1, l, j+1, 0.);
291  }
292  }
293 }
virtual numeric_index_type m() const libmesh_override
dof_id_type numeric_index_type
Definition: id_types.h:92
std::vector< std::vector< numeric_index_type >::const_iterator > _row_start
template<typename T>
void libMesh::SparseMatrix< T >::zero_rows ( std::vector< numeric_index_type > &  rows,
diag_value = 0.0 
)
virtualinherited

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

Reimplemented in libMesh::PetscMatrix< T >.

Definition at line 196 of file sparse_matrix.C.

197 {
198  /* This functionality isn't implemented or stubbed in every subclass yet */
199  libmesh_not_implemented();
200 }

Friends And Related Function Documentation

template<typename T>
friend class LaspackLinearSolver< T >
friend

Definition at line 304 of file laspack_matrix.h.

template<typename T>
friend class LaspackVector< T >
friend

Make other Laspack datatypes friends

Definition at line 303 of file laspack_matrix.h.

Member Data Documentation

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

Flag indicating if the matrix has been closed yet.

Definition at line 298 of file laspack_matrix.h.

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

ReferenceCounter::Counts libMesh::ReferenceCounter::_counts
staticprotectedinherited
template<typename T>
std::vector<numeric_index_type> libMesh::LaspackMatrix< T >::_csr
private

The compressed row indices.

Definition at line 287 of file laspack_matrix.h.

Referenced by libMesh::LaspackMatrix< T >::clear(), and libMesh::LaspackMatrix< T >::pos().

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

The DofMap object associated with this object.

Definition at line 436 of file sparse_matrix.h.

Referenced by libMesh::SparseMatrix< T >::print().

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 134 of file reference_counter.h.

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

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

Flag indicating whether or not the matrix has been initialized.

Definition at line 442 of file sparse_matrix.h.

Referenced by libMesh::PetscMatrix< T >::_get_submatrix(), libMesh::EigenSparseMatrix< T >::clear(), libMesh::LaspackMatrix< T >::clear(), and libMesh::PetscMatrix< T >::get_transpose().

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

Mutual exclusion object to enable thread-safe reference counting.

Definition at line 128 of file reference_counter.h.

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 123 of file reference_counter.h.

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

template<typename T>
std::vector<std::vector<numeric_index_type>::const_iterator> libMesh::LaspackMatrix< T >::_row_start
private

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