libMesh::EpetraMatrix< T > Class Template Reference

#include <trilinos_epetra_matrix.h>

Inheritance diagram for libMesh::EpetraMatrix< T >:

Public Member Functions

 EpetraMatrix (const Parallel::Communicator &comm LIBMESH_CAN_DEFAULT_TO_COMMWORLD)
 
 EpetraMatrix (Epetra_FECrsMatrix *m, const Parallel::Communicator &comm LIBMESH_CAN_DEFAULT_TO_COMMWORLD)
 
virtual ~EpetraMatrix ()
 
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
 
void swap (EpetraMatrix< T > &)
 
Epetra_FECrsMatrix * mat ()
 
const Epetra_FECrsMatrix * mat () const
 
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 Attributes

Epetra_FECrsMatrix * _mat
 
Epetra_Map * _map
 
Epetra_CrsGraph * _graph
 
bool _destroy_mat_on_exit
 
bool _use_transpose
 

Detailed Description

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

Epetra matrix. Provides a nice interface to the Epetra data structures for parallel, sparse matrices.

Author
Benjamin S. Kirk
Date
2008

Definition at line 63 of file trilinos_epetra_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::EpetraMatrix< T >::EpetraMatrix ( 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(...).

Referenced by libMesh::EpetraMatrix< T >::get_transpose().

template<typename T>
libMesh::EpetraMatrix< T >::EpetraMatrix ( Epetra_FECrsMatrix *  m,
const Parallel::Communicator &comm  LIBMESH_CAN_DEFAULT_TO_COMMWORLD 
)

Constructor. Creates a EpetraMatrix assuming you already have a valid Epetra_FECrsMatrix object. In this case, m is NOT destroyed by the EpetraMatrix destructor when this object goes out of scope. This allows ownership of m to remain with the original creator, and to simply provide additional functionality with the EpetraMatrix.

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

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

Definition at line 314 of file trilinos_epetra_matrix.C.

315 {
316  this->clear();
317 }
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::EpetraMatrix< 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 395 of file trilinos_epetra_matrix.C.

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

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

398 {
399  libmesh_assert (this->initialized());
400 
401  int
402  epetra_i = static_cast<int>(i),
403  epetra_j = static_cast<int>(j);
404 
405  T epetra_value = value;
406 
407  _mat->SumIntoGlobalValues (epetra_i, 1, &epetra_value, &epetra_j);
408 }
libmesh_assert(j)
virtual bool initialized() const
Epetra_FECrsMatrix * _mat
template<typename T >
void libMesh::EpetraMatrix< 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 + \texttt{this} $. It is advisable to not only allocate appropriate memory with init() , but also explicitly zero the terms of this whenever you add a non-zero value to X. Note: X will be closed, if not already done, before performing any work.

Implements libMesh::SparseMatrix< T >.

Definition at line 422 of file trilinos_epetra_matrix.C.

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

423 {
424 #ifdef LIBMESH_TRILINOS_HAVE_EPETRAEXT
425  libmesh_assert (this->initialized());
426 
427  // sanity check. but this cannot avoid
428  // crash due to incompatible sparsity structure...
429  libmesh_assert_equal_to (this->m(), X_in.m());
430  libmesh_assert_equal_to (this->n(), X_in.n());
431 
432  EpetraMatrix<T> * X = cast_ptr<EpetraMatrix<T> *> (&X_in);
433 
434  EpetraExt::MatrixMatrix::Add (*X->_mat, false, a_in, *_mat, 1.);
435 #else
436  libmesh_error_msg("ERROR: EpetraExt is required for EpetraMatrix::add()!");
437 #endif
438 }
libmesh_assert(j)
virtual bool initialized() const
virtual numeric_index_type m() const libmesh_override
virtual numeric_index_type n() const libmesh_override
Epetra_FECrsMatrix * _mat
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::EpetraMatrix< 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 Petsc matrix. This is useful for adding an element matrix at assembly time

Implements libMesh::SparseMatrix< T >.

Definition at line 242 of file trilinos_epetra_matrix.C.

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

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

245 {
246  libmesh_assert (this->initialized());
247 
248  const numeric_index_type m = dm.m();
249  const numeric_index_type n = dm.n();
250 
251  libmesh_assert_equal_to (rows.size(), m);
252  libmesh_assert_equal_to (cols.size(), n);
253 
254  _mat->SumIntoGlobalValues(m, (int *)&rows[0], n, (int *)&cols[0], &dm.get_values()[0]);
255 }
libmesh_assert(j)
virtual bool initialized() const
virtual numeric_index_type m() const libmesh_override
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual numeric_index_type n() const libmesh_override
Epetra_FECrsMatrix * _mat
template<typename T >
void libMesh::EpetraMatrix< 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 413 of file trilinos_epetra_matrix.C.

415 {
416  this->add_matrix (dm, dof_indices, dof_indices);
417 }
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::EpetraMatrix< 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 204 of file trilinos_epetra_matrix.C.

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

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

205 {
206  // delete _mat;
207  // delete _map;
208 
209  this->_is_initialized = false;
210 
211  libmesh_assert (!this->initialized());
212 }
libmesh_assert(j)
virtual bool initialized() const
template<typename T >
void libMesh::EpetraMatrix< T >::close ( ) const
virtual

Call the Petsc assemble routines. sends necessary messages to other processors

Implements libMesh::SparseMatrix< T >.

Definition at line 322 of file trilinos_epetra_matrix.C.

References libMesh::libmesh_assert().

Referenced by libMesh::EpetraMatrix< T >::need_full_sparsity_pattern(), and libMesh::AztecLinearSolver< T >::solve().

323 {
325 
326  _mat->GlobalAssemble();
327 }
libmesh_assert(j)
Epetra_FECrsMatrix * _mat
template<typename T >
bool libMesh::EpetraMatrix< T >::closed ( ) const
virtual

see if Petsc matrix has been closed and fully assembled yet

Implements libMesh::SparseMatrix< T >.

Definition at line 479 of file trilinos_epetra_matrix.C.

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

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

480 {
481  libmesh_assert (this->initialized());
482  libmesh_assert(this->_mat);
483 
484  return this->_mat->Filled();
485 }
libmesh_assert(j)
virtual bool initialized() const
Epetra_FECrsMatrix * _mat
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::EpetraMatrix< T >::get_diagonal ( NumericVector< T > &  dest) const
virtual

Copies the diagonal part of the matrix into dest.

Implements libMesh::SparseMatrix< T >.

Definition at line 263 of file trilinos_epetra_matrix.C.

References libMesh::EpetraVector< T >::vec().

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

264 {
265  // Convert vector to EpetraVector.
266  EpetraVector<T> * epetra_dest = cast_ptr<EpetraVector<T> *>(&dest);
267 
268  // Call Epetra function.
269  _mat->ExtractDiagonalCopy(*(epetra_dest->vec()));
270 }
Epetra_FECrsMatrix * _mat
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::EpetraMatrix< 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 275 of file trilinos_epetra_matrix.C.

References libMesh::libMeshPrivateData::_is_initialized, libMesh::EpetraMatrix< T >::_mat, libMesh::EpetraMatrix< T >::_use_transpose, and libMesh::EpetraMatrix< T >::EpetraMatrix().

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

276 {
277  // Make sure the SparseMatrix passed in is really a EpetraMatrix
278  EpetraMatrix<T> & epetra_dest = cast_ref<EpetraMatrix<T> &>(dest);
279 
280  if(&epetra_dest != this)
281  epetra_dest = *this;
282 
283  epetra_dest._use_transpose = !epetra_dest._use_transpose;
284  epetra_dest._mat->SetUseTranspose(epetra_dest._use_transpose);
285 }
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::EpetraMatrix< 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 Petsc 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::EpetraMatrix< T >::init ( )
virtual

Initialize using sparsity structure computed by dof_map.

Implements libMesh::SparseMatrix< T >.

Definition at line 175 of file trilinos_epetra_matrix.C.

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

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

176 {
177  libmesh_assert(this->_dof_map);
178 
179  {
180  // Clear initialized matrices
181  if (this->initialized())
182  this->clear();
183 
184  this->_is_initialized = true;
185  }
186 
187 
188  _mat = new Epetra_FECrsMatrix (Copy, *_graph);
189 }
libmesh_assert(j)
virtual void clear() libmesh_override
virtual bool initialized() const
DofMap const * _dof_map
Epetra_FECrsMatrix * _mat
template<typename T >
Real libMesh::EpetraMatrix< T >::l1_norm ( ) const
virtual

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$. (cf. Haemmerlin-Hoffmann : Numerische Mathematik)

Implements libMesh::SparseMatrix< T >.

Definition at line 217 of file trilinos_epetra_matrix.C.

References libMesh::initialized(), libMesh::libmesh_assert(), and libMesh::Real.

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

218 {
219  libmesh_assert (this->initialized());
220 
222 
223  return static_cast<Real>(_mat->NormOne());
224 }
libmesh_assert(j)
virtual bool initialized() const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Epetra_FECrsMatrix * _mat
template<typename T >
Real libMesh::EpetraMatrix< T >::linfty_norm ( ) const
virtual

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$. (cf. Haemmerlin-Hoffmann : Numerische Mathematik)

Implements libMesh::SparseMatrix< T >.

Definition at line 229 of file trilinos_epetra_matrix.C.

References libMesh::initialized(), libMesh::libmesh_assert(), and libMesh::Real.

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

230 {
231  libmesh_assert (this->initialized());
232 
233 
235 
236  return static_cast<Real>(_mat->NormInf());
237 }
libmesh_assert(j)
virtual bool initialized() const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Epetra_FECrsMatrix * _mat
template<typename T >
numeric_index_type libMesh::EpetraMatrix< T >::m ( ) const
virtual
Returns
m, the row-dimension of the matrix where the marix is $ M \times N $.

Implements libMesh::SparseMatrix< T >.

Definition at line 332 of file trilinos_epetra_matrix.C.

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

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

333 {
334  libmesh_assert (this->initialized());
335 
336  return static_cast<numeric_index_type>(_mat->NumGlobalRows());
337 }
libmesh_assert(j)
virtual bool initialized() const
dof_id_type numeric_index_type
Definition: id_types.h:92
Epetra_FECrsMatrix * _mat
template<typename T>
Epetra_FECrsMatrix* libMesh::EpetraMatrix< T >::mat ( )
inline

Returns the raw PETSc matrix context pointer. Note this is generally not required in user-level code. Just don't do anything crazy like calling LibMeshMatDestroy()!

Definition at line 302 of file trilinos_epetra_matrix.h.

References libMesh::EpetraMatrix< T >::_mat, and libMesh::libmesh_assert().

Referenced by libMesh::TrilinosPreconditioner< T >::init(), libMesh::NoxNonlinearSolver< T >::solve(), and libMesh::AztecLinearSolver< T >::solve().

302 { libmesh_assert(_mat); return _mat; }
libmesh_assert(j)
Epetra_FECrsMatrix * _mat
template<typename T>
const Epetra_FECrsMatrix* libMesh::EpetraMatrix< T >::mat ( ) const
inline

Definition at line 304 of file trilinos_epetra_matrix.h.

References libMesh::EpetraMatrix< T >::_mat, and libMesh::libmesh_assert().

304 { libmesh_assert(_mat); return _mat; }
libmesh_assert(j)
Epetra_FECrsMatrix * _mat
template<typename T >
numeric_index_type libMesh::EpetraMatrix< T >::n ( ) const
virtual
Returns
n, the column-dimension of the matrix where the marix is $ M \times N $.

Implements libMesh::SparseMatrix< T >.

Definition at line 342 of file trilinos_epetra_matrix.C.

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

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

343 {
344  libmesh_assert (this->initialized());
345 
346  return static_cast<numeric_index_type>(_mat->NumGlobalCols());
347 }
libmesh_assert(j)
virtual bool initialized() const
dof_id_type numeric_index_type
Definition: id_types.h:92
Epetra_FECrsMatrix * _mat
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 >
T libMesh::EpetraMatrix< 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 444 of file trilinos_epetra_matrix.C.

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

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

446 {
447  libmesh_assert (this->initialized());
448  libmesh_assert(this->_mat);
449  libmesh_assert (this->_mat->MyGlobalRow(static_cast<int>(i)));
450  libmesh_assert_greater_equal (i, this->row_start());
451  libmesh_assert_less (i, this->row_stop());
452 
453 
454  int row_length;
455  int * row_indices;
456  double * values;
457 
458  _mat->ExtractMyRowView (i-this->row_start(),
459  row_length,
460  values,
461  row_indices);
462 
463  //libMesh::out << "row_length=" << row_length << std::endl;
464 
465  int * index = std::lower_bound (row_indices, row_indices+row_length, j);
466 
467  libmesh_assert_less (*index, row_length);
468  libmesh_assert_equal_to (static_cast<numeric_index_type>(row_indices[*index]), j);
469 
470  //libMesh::out << "val=" << values[*index] << std::endl;
471 
472  return values[*index];
473 }
virtual numeric_index_type row_stop() const libmesh_override
libmesh_assert(j)
virtual bool initialized() const
virtual numeric_index_type row_start() const libmesh_override
Epetra_FECrsMatrix * _mat
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 >
void libMesh::EpetraMatrix< T >::print_personal ( std::ostream &  os = libMesh::out) const
virtual

Print the contents of the matrix, by default to libMesh::out.

Implements libMesh::SparseMatrix< T >.

Definition at line 500 of file trilinos_epetra_matrix.C.

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

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

501 {
502  libmesh_assert (this->initialized());
504 
505  os << *_mat;
506 }
libmesh_assert(j)
virtual bool initialized() const
Epetra_FECrsMatrix * _mat
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::EpetraMatrix< 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 352 of file trilinos_epetra_matrix.C.

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

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

353 {
354  libmesh_assert (this->initialized());
356 
357  return static_cast<numeric_index_type>(_map->MinMyGID());
358 }
libmesh_assert(j)
virtual bool initialized() const
dof_id_type numeric_index_type
Definition: id_types.h:92
template<typename T >
numeric_index_type libMesh::EpetraMatrix< 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 363 of file trilinos_epetra_matrix.C.

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

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

364 {
365  libmesh_assert (this->initialized());
367 
368  return static_cast<numeric_index_type>(_map->MaxMyGID())+1;
369 }
libmesh_assert(j)
virtual bool initialized() const
dof_id_type numeric_index_type
Definition: id_types.h:92
template<typename T >
void libMesh::EpetraMatrix< 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 374 of file trilinos_epetra_matrix.C.

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

377 {
378  libmesh_assert (this->initialized());
379 
380  int
381  epetra_i = static_cast<int>(i),
382  epetra_j = static_cast<int>(j);
383 
384  T epetra_value = value;
385 
386  if (_mat->Filled())
387  _mat->ReplaceGlobalValues (epetra_i, 1, &epetra_value, &epetra_j);
388  else
389  _mat->InsertGlobalValues (epetra_i, 1, &epetra_value, &epetra_j);
390 }
libmesh_assert(j)
virtual bool initialized() const
Epetra_FECrsMatrix * _mat
template<typename T >
void libMesh::EpetraMatrix< T >::swap ( EpetraMatrix< T > &  m)

Swaps the raw PETSc matrix context pointers.

Definition at line 489 of file trilinos_epetra_matrix.C.

References libMesh::EpetraMatrix< T >::_destroy_mat_on_exit, libMesh::EpetraMatrix< T >::_mat, and swap().

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

490 {
491  std::swap(_mat, m._mat);
492  std::swap(_destroy_mat_on_exit, m._destroy_mat_on_exit);
493 }
virtual numeric_index_type m() const libmesh_override
void swap(Iterator &lhs, Iterator &rhs)
Epetra_FECrsMatrix * _mat
template<typename T >
void libMesh::EpetraMatrix< 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 41 of file trilinos_epetra_matrix.C.

References libMesh::libMeshPrivateData::_is_initialized, libMesh::TriangleWrapper::init(), libMesh::EpetraMatrix< T >::init(), libMesh::initialized(), libMesh::libmesh_assert(), libMesh::libmesh_dbg_var(), std::min(), and libMesh::processor_id().

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

42 {
43  // clear data, start over
44  this->clear ();
45 
46  // big trouble if this fails!
47  libmesh_assert(this->_dof_map);
48 
49  const numeric_index_type n_rows = sparsity_pattern.size();
50 
51  const numeric_index_type m = this->_dof_map->n_dofs();
52  const numeric_index_type n = m;
53  const numeric_index_type n_l = this->_dof_map->n_dofs_on_processor(this->processor_id());
54  const numeric_index_type m_l = n_l;
55 
56  // error checking
57 #ifndef NDEBUG
58  {
59  libmesh_assert_equal_to (n, m);
60  libmesh_assert_equal_to (n_l, m_l);
61 
63  summed_m_l = m_l,
64  summed_n_l = n_l;
65 
66  this->comm().sum (summed_m_l);
67  this->comm().sum (summed_n_l);
68 
69  libmesh_assert_equal_to (m, summed_m_l);
70  libmesh_assert_equal_to (n, summed_n_l);
71  }
72 #endif
73 
74  // build a map defining the data distribution
75  _map = new Epetra_Map (static_cast<int>(m),
76  m_l,
77  0,
78  Epetra_MpiComm (this->comm().get()));
79 
80  libmesh_assert_equal_to (static_cast<numeric_index_type>(_map->NumGlobalPoints()), m);
81  libmesh_assert_equal_to (static_cast<numeric_index_type>(_map->MaxAllGID()+1), m);
82 
83  const std::vector<numeric_index_type> & n_nz = this->_dof_map->get_n_nz();
84  const std::vector<numeric_index_type> & n_oz = this->_dof_map->get_n_oz();
85 
86  // Make sure the sparsity pattern isn't empty
87  libmesh_assert_equal_to (n_nz.size(), n_l);
88  libmesh_assert_equal_to (n_oz.size(), n_l);
89 
90  // Epetra wants the total number of nonzeros, both local and remote.
91  std::vector<int> n_nz_tot; n_nz_tot.reserve(n_nz.size());
92 
93  for (numeric_index_type i=0; i<n_nz.size(); i++)
94  n_nz_tot.push_back(std::min(n_nz[i] + n_oz[i], n));
95 
96  if (m==0)
97  return;
98 
99  _graph = new Epetra_CrsGraph(Copy, *_map, &n_nz_tot[0]);
100 
101  // Tell the matrix about its structure. Initialize it
102  // to zero.
103  for (numeric_index_type i=0; i<n_rows; i++)
104  _graph->InsertGlobalIndices(_graph->GRID(i),
105  sparsity_pattern[i].size(),
106  const_cast<int *>((const int *)&sparsity_pattern[i][0]));
107 
108  _graph->FillComplete();
109 
110  //Initialize the matrix
111  libmesh_assert (!this->initialized());
112  this->init ();
113  libmesh_assert (this->initialized());
114 }
dof_id_type n_dofs() const
Definition: dof_map.h:506
dof_id_type n_dofs_on_processor(const processor_id_type proc) const
Definition: dof_map.h:522
libmesh_assert(j)
virtual void clear() libmesh_override
virtual bool initialized() const
virtual void init() libmesh_override
virtual numeric_index_type m() const libmesh_override
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
virtual numeric_index_type n() const libmesh_override
const Parallel::Communicator & comm() const
const std::vector< dof_id_type > & get_n_nz() const
Definition: dof_map.h:400
long double min(long double a, double b)
processor_id_type processor_id() 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::EpetraMatrix< T >::zero ( )
virtual

Set all entries to 0. This method retains sparsity structure.

Implements libMesh::SparseMatrix< T >.

Definition at line 194 of file trilinos_epetra_matrix.C.

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

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

195 {
196  libmesh_assert (this->initialized());
197 
198  _mat->Scale(0.0);
199 }
libmesh_assert(j)
virtual bool initialized() const
Epetra_FECrsMatrix * _mat
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 }

Member Data Documentation

ReferenceCounter::Counts libMesh::ReferenceCounter::_counts
staticprotectedinherited
template<typename T>
bool libMesh::EpetraMatrix< T >::_destroy_mat_on_exit
private

This boolean value should only be set to false for the constructor which takes a PETSc Mat object.

Definition at line 329 of file trilinos_epetra_matrix.h.

Referenced by libMesh::EpetraMatrix< T >::swap().

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>
Epetra_CrsGraph* libMesh::EpetraMatrix< T >::_graph
private

Holds the sparsity pattern

Definition at line 323 of file trilinos_epetra_matrix.h.

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

template<typename T>
Epetra_Map* libMesh::EpetraMatrix< T >::_map
private

Holds the distributed Map

Definition at line 318 of file trilinos_epetra_matrix.h.

template<typename T>
Epetra_FECrsMatrix* libMesh::EpetraMatrix< T >::_mat
private

Actual Epetra datatype to hold matrix entries

Definition at line 313 of file trilinos_epetra_matrix.h.

Referenced by libMesh::EpetraMatrix< T >::get_transpose(), libMesh::EpetraMatrix< T >::mat(), and libMesh::EpetraMatrix< T >::swap().

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>
bool libMesh::EpetraMatrix< T >::_use_transpose
private

Epetra has no GetUseTranspose so we need to keep track of whether we're transposed manually.

Definition at line 335 of file trilinos_epetra_matrix.h.

Referenced by libMesh::EpetraMatrix< T >::get_transpose().


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