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 ()
 
bool need_full_sparsity_pattern () const
 
void update_sparsity_pattern (const SparsityPattern::Graph &)
 
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)
 
void init ()
 
void clear ()
 
void zero ()
 
void close () const
 
numeric_index_type m () const
 
numeric_index_type n () const
 
numeric_index_type row_start () const
 
numeric_index_type row_stop () const
 
void set (const numeric_index_type i, const numeric_index_type j, const T value)
 
void add (const numeric_index_type i, const numeric_index_type j, const T value)
 
void add_matrix (const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols)
 
void add_matrix (const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &dof_indices)
 
void add (const T a, SparseMatrix< T > &X)
 
operator() (const numeric_index_type i, const numeric_index_type j) const
 
Real l1_norm () const
 
Real linfty_norm () const
 
bool closed () const
 
void print_personal (std::ostream &os=libMesh::out) const
 
void print_matlab (const std::string name="NULL") const
 
void get_diagonal (NumericVector< T > &dest) const
 
virtual void get_transpose (SparseMatrix< T > &dest) const
 
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 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 AutoPtr< 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, 2008

Definition at line 57 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 113 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(...).

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 328 of file trilinos_epetra_matrix.C.

329 {
330  this->clear();
331 }

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.

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

429  {
430  libMesh::err << "Error! This function is not yet implemented in the base class!"
431  << std::endl;
432  libmesh_error();
433  }
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 409 of file trilinos_epetra_matrix.C.

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

412 {
413  libmesh_assert (this->initialized());
414 
415  int
416  epetra_i = static_cast<int>(i),
417  epetra_j = static_cast<int>(j);
418 
419  T epetra_value = value;
420 
421  _mat->SumIntoGlobalValues (epetra_i, 1, &epetra_value, &epetra_j);
422 }
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 436 of file trilinos_epetra_matrix.C.

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

437 {
438  libmesh_assert (this->initialized());
439 
440  // sanity check. but this cannot avoid
441  // crash due to incompatible sparsity structure...
442  libmesh_assert_equal_to (this->m(), X_in.m());
443  libmesh_assert_equal_to (this->n(), X_in.n());
444 
445  EpetraMatrix<T>* X = libmesh_cast_ptr<EpetraMatrix<T>*> (&X_in);
446 
447  EpetraExt::MatrixMatrix::Add (*X->_mat, false, a_in, *_mat, 1.);
448 }
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::DenseMatrixBase< T >::m(), and libMesh::DenseMatrixBase< T >::n().

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

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

259 {
260  libmesh_assert (this->initialized());
261 
262  const numeric_index_type m = dm.m();
263  const numeric_index_type n = dm.n();
264 
265  libmesh_assert_equal_to (rows.size(), m);
266  libmesh_assert_equal_to (cols.size(), n);
267 
268  _mat->SumIntoGlobalValues(m, (int *)&rows[0], n, (int *)&cols[0], &dm.get_values()[0]);
269 }
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 427 of file trilinos_epetra_matrix.C.

429 {
430  this->add_matrix (dm, dof_indices, dof_indices);
431 }
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::DofMap::attach_matrix().

113  { _dof_map = &dof_map; }
template<typename T >
AutoPtr< 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 134 of file sparse_matrix.C.

References libMesh::EIGEN_SOLVERS, libMesh::err, libMesh::LASPACK_SOLVERS, libMeshEnums::PETSC_SOLVERS, and libMesh::TRILINOS_SOLVERS.

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

136 {
137  // Build the appropriate vector
138  switch (solver_package)
139  {
140 
141 
142 #ifdef LIBMESH_HAVE_LASPACK
143  case LASPACK_SOLVERS:
144  {
145  AutoPtr<SparseMatrix<T> > ap(new LaspackMatrix<T>(comm));
146  return ap;
147  }
148 #endif
149 
150 
151 #ifdef LIBMESH_HAVE_PETSC
152  case PETSC_SOLVERS:
153  {
154  AutoPtr<SparseMatrix<T> > ap(new PetscMatrix<T>(comm));
155  return ap;
156  }
157 #endif
158 
159 
160 #ifdef LIBMESH_HAVE_TRILINOS
161  case TRILINOS_SOLVERS:
162  {
163  AutoPtr<SparseMatrix<T> > ap(new EpetraMatrix<T>(comm));
164  return ap;
165  }
166 #endif
167 
168 
169 #ifdef LIBMESH_HAVE_EIGEN
170  case EIGEN_SOLVERS:
171  {
172  AutoPtr<SparseMatrix<T> > ap(new EigenSparseMatrix<T>(comm));
173  return ap;
174  }
175 #endif
176 
177 
178  default:
179  libMesh::err << "ERROR: Unrecognized solver package: "
180  << solver_package
181  << std::endl;
182  libmesh_error();
183  }
184 
185  AutoPtr<SparseMatrix<T> > ap(NULL);
186  return ap;
187 }
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().

205 {
206 // delete _mat;
207 // delete _map;
208 
209  this->_is_initialized = false;
210 
211  libmesh_assert (!this->initialized());
212 }
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 336 of file trilinos_epetra_matrix.C.

References libMesh::libmesh_assert().

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

337 {
339 
340  _mat->GlobalAssemble();
341 }
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 488 of file trilinos_epetra_matrix.C.

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

489 {
490  libmesh_assert (this->initialized());
491  libmesh_assert(this->_mat);
492 
493  return this->_mat->Filled();
494 }
const Parallel::Communicator& libMesh::ParallelObject::comm ( ) const
inlineinherited
Returns
a reference to the Parallel::Communicator object used by this mesh.

Definition at line 86 of file parallel_object.h.

References libMesh::ParallelObject::_communicator.

Referenced by libMesh::__libmesh_petsc_diff_solver_jacobian(), libMesh::__libmesh_petsc_diff_solver_monitor(), libMesh::__libmesh_petsc_diff_solver_residual(), libMesh::__libmesh_petsc_snes_jacobian(), libMesh::__libmesh_petsc_snes_residual(), libMesh::MeshRefinement::_coarsen_elements(), libMesh::ExactSolution::_compute_error(), libMesh::MetisPartitioner::_do_partition(), libMesh::ParmetisPartitioner::_do_repartition(), libMesh::UniformRefinementEstimator::_estimate_error(), libMesh::SlepcEigenSolver< T >::_petsc_shell_matrix_get_diagonal(), libMesh::PetscLinearSolver< T >::_petsc_shell_matrix_get_diagonal(), libMesh::SlepcEigenSolver< T >::_petsc_shell_matrix_mult(), libMesh::PetscLinearSolver< T >::_petsc_shell_matrix_mult(), libMesh::PetscLinearSolver< T >::_petsc_shell_matrix_mult_add(), libMesh::EquationSystems::_read_impl(), libMesh::MeshRefinement::_refine_elements(), libMesh::ParallelMesh::add_elem(), libMesh::ImplicitSystem::add_matrix(), libMesh::ParallelMesh::add_node(), libMesh::System::add_vector(), libMesh::UnstructuredMesh::all_second_order(), libMesh::LaplaceMeshSmoother::allgather_graph(), libMesh::FEMSystem::assemble_qoi(), libMesh::MeshCommunication::assign_global_indices(), libMesh::ParmetisPartitioner::assign_partitioning(), libMesh::DofMap::attach_matrix(), libMesh::MeshTools::bounding_box(), libMesh::System::calculate_norm(), libMesh::MeshRefinement::coarsen_elements(), 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::MeshCommunication::delete_remote_elements(), libMesh::DofMap::distribute_dofs(), DMlibMeshFunction(), DMlibMeshJacobian(), DMLibMeshSetSystem(), DMVariableBounds_libMesh(), libMesh::MeshRefinement::eliminate_unrefined_patches(), libMesh::WeightedPatchRecoveryErrorEstimator::estimate_error(), libMesh::PatchRecoveryErrorEstimator::estimate_error(), libMesh::JumpErrorEstimator::estimate_error(), libMesh::AdjointRefinementEstimator::estimate_error(), libMesh::MeshRefinement::flag_elements_by_elem_fraction(), libMesh::MeshRefinement::flag_elements_by_error_fraction(), libMesh::MeshRefinement::flag_elements_by_nelem_target(), libMesh::for(), libMesh::CondensedEigenSystem::get_eigenpair(), libMesh::ImplicitSystem::get_linear_solver(), libMesh::LocationMap< T >::init(), libMesh::PetscDiffSolver::init(), libMesh::TimeSolver::init(), libMesh::SystemSubsetBySubdomain::init(), libMesh::EigenSystem::init_data(), libMesh::EigenSystem::init_matrices(), libMesh::ParmetisPartitioner::initialize(), libMesh::MeshTools::libmesh_assert_valid_dof_ids(), libMesh::ParallelMesh::libmesh_assert_valid_parallel_flags(), libMesh::MeshTools::libmesh_assert_valid_procids< Elem >(), libMesh::MeshTools::libmesh_assert_valid_procids< Node >(), libMesh::MeshTools::libmesh_assert_valid_refinement_flags(), libMesh::MeshRefinement::limit_level_mismatch_at_edge(), libMesh::MeshRefinement::limit_level_mismatch_at_node(), libMesh::MeshRefinement::make_coarsening_compatible(), libMesh::MeshCommunication::make_elems_parallel_consistent(), libMesh::MeshRefinement::make_flags_parallel_consistent(), libMesh::MeshCommunication::make_node_ids_parallel_consistent(), libMesh::MeshCommunication::make_node_proc_ids_parallel_consistent(), libMesh::MeshCommunication::make_nodes_parallel_consistent(), libMesh::MeshRefinement::make_refinement_compatible(), libMesh::FEMSystem::mesh_position_set(), libMesh::MeshSerializer::MeshSerializer(), libMesh::ParallelMesh::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::ParallelMesh::parallel_max_elem_id(), libMesh::ParallelMesh::parallel_max_node_id(), libMesh::ParallelMesh::parallel_n_elem(), libMesh::ParallelMesh::parallel_n_nodes(), libMesh::Partitioner::partition(), libMesh::Partitioner::partition_unpartitioned_elements(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::MeshBase::prepare_for_use(), libMesh::System::project_vector(), libMesh::Nemesis_IO::read(), libMesh::XdrIO::read(), libMesh::System::read_header(), libMesh::System::read_legacy_data(), libMesh::System::read_SCALAR_dofs(), libMesh::XdrIO::read_serialized_bc_names(), libMesh::XdrIO::read_serialized_bcs(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::XdrIO::read_serialized_nodes(), libMesh::XdrIO::read_serialized_nodesets(), libMesh::XdrIO::read_serialized_subdomain_names(), libMesh::System::read_serialized_vector(), libMesh::MeshBase::recalculate_n_partitions(), libMesh::MeshRefinement::refine_and_coarsen_elements(), libMesh::MeshRefinement::refine_elements(), libMesh::Partitioner::set_node_processor_ids(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::LaplaceMeshSmoother::smooth(), libMesh::MeshBase::subdomain_ids(), libMesh::BoundaryInfo::sync(), libMesh::Parallel::sync_element_data_by_parent_id(), libMesh::MeshRefinement::test_level_one(), libMesh::MeshRefinement::test_unflagged(), libMesh::MeshTools::total_weight(), libMesh::CheckpointIO::write(), libMesh::XdrIO::write(), libMesh::UnstructuredMesh::write(), libMesh::LegacyXdrIO::write_mesh(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bcs(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), libMesh::XdrIO::write_serialized_nodesets(), and libMesh::DivaIO::write_stream().

87  { return _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.

Referenced by libMesh::CondensedEigenSystem::solve().

371  {
372  this->_get_submatrix(submatrix,
373  rows,
374  cols,
375  false); // false means DO NOT REUSE submatrix
376  }
void libMesh::ReferenceCounter::disable_print_counter_info ( )
staticinherited

Definition at line 106 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

107 {
108  _enable_print_counter = false;
109  return;
110 }
void libMesh::ReferenceCounter::enable_print_counter_info ( )
staticinherited

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

Definition at line 100 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

101 {
102  _enable_print_counter = true;
103  return;
104 }
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 277 of file trilinos_epetra_matrix.C.

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

278 {
279  // Convert vector to EpetraVector.
280  EpetraVector<T>* epetra_dest = libmesh_cast_ptr<EpetraVector<T>*>(&dest);
281 
282  // Call Epetra function.
283  _mat->ExtractDiagonalCopy(*(epetra_dest->vec()));
284 }
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 }
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 289 of file trilinos_epetra_matrix.C.

References libMesh::EpetraMatrix< T >::_mat, and libMesh::EpetraMatrix< T >::_use_transpose.

290 {
291  // Make sure the SparseMatrix passed in is really a EpetraMatrix
292  EpetraMatrix<T>& epetra_dest = libmesh_cast_ref<EpetraMatrix<T>&>(dest);
293 
294  if(&epetra_dest != this)
295  epetra_dest = *this;
296 
297  epetra_dest._use_transpose = !epetra_dest._use_transpose;
298  epetra_dest._mat->SetUseTranspose(epetra_dest._use_transpose);
299 }
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 163 of file reference_counter.h.

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

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

164 {
165  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
166  std::pair<unsigned int, unsigned int>& p = _counts[name];
167 
168  p.first++;
169 }
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 176 of file reference_counter.h.

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

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

177 {
178  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
179  std::pair<unsigned int, unsigned int>& p = _counts[name];
180 
181  p.second++;
182 }
template<typename T >
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 >.

Definition at line 119 of file trilinos_epetra_matrix.C.

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

126 {
127  if ((m==0) || (n==0))
128  return;
129 
130  {
131  // Clear initialized matrices
132  if (this->initialized())
133  this->clear();
134 
135  libmesh_assert (!this->_mat);
136  libmesh_assert (!this->_map);
137 
138  this->_is_initialized = true;
139  }
140 
141  // error checking
142 #ifndef NDEBUG
143  {
144  libmesh_assert_equal_to (n, m);
145  libmesh_assert_equal_to (n_l, m_l);
146 
148  summed_m_l = m_l,
149  summed_n_l = n_l;
150 
151  this->comm().sum (summed_m_l);
152  this->comm().sum (summed_n_l);
153 
154  libmesh_assert_equal_to (m, summed_m_l);
155  libmesh_assert_equal_to (n, summed_n_l);
156  }
157 #endif
158 
159  // build a map defining the data distribution
160  _map = new Epetra_Map (static_cast<int>(m),
161  m_l,
162  0,
163  Epetra_MpiComm (this->comm().get()));
164 
165  libmesh_assert_equal_to (static_cast<numeric_index_type>(_map->NumGlobalPoints()), m);
166  libmesh_assert_equal_to (static_cast<numeric_index_type>(_map->MaxAllGID()+1), m);
167 
168  _mat = new Epetra_FECrsMatrix (Copy, *_map, nnz + noz);
169 }
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().

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 }
template<typename T>
virtual bool libMesh::SparseMatrix< T >::initialized ( ) const
inlinevirtualinherited
Returns
true if the matrix has been initialized, false otherwise.

Definition at line 107 of file sparse_matrix.h.

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

107 { return _is_initialized; }
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.

218 {
219  libmesh_assert (this->initialized());
220 
222 
223  return static_cast<Real>(_mat->NormOne());
224 }
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.

230 {
231  libmesh_assert (this->initialized());
232 
233 
235 
236  return static_cast<Real>(_mat->NormInf());
237 }
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 346 of file trilinos_epetra_matrix.C.

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

347 {
348  libmesh_assert (this->initialized());
349 
350  return static_cast<numeric_index_type>(_mat->NumGlobalRows());
351 }
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 305 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().

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

Definition at line 307 of file trilinos_epetra_matrix.h.

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

307 { libmesh_assert(_mat); return _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 356 of file trilinos_epetra_matrix.C.

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

357 {
358  libmesh_assert (this->initialized());
359 
360  return static_cast<numeric_index_type>(_mat->NumGlobalCols());
361 }
static unsigned int libMesh::ReferenceCounter::n_objects ( )
inlinestaticinherited

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

Definition at line 79 of file reference_counter.h.

References libMesh::ReferenceCounter::_n_objects.

80  { return _n_objects; }
processor_id_type libMesh::ParallelObject::n_processors ( ) const
inlineinherited
Returns
the number of processors in the group.

Definition at line 92 of file parallel_object.h.

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

Referenced by libMesh::ParmetisPartitioner::_do_repartition(), libMesh::ParallelMesh::add_elem(), libMesh::ParallelMesh::add_node(), libMesh::LaplaceMeshSmoother::allgather_graph(), libMesh::ParmetisPartitioner::assign_partitioning(), libMesh::ParallelMesh::assign_unique_ids(), libMesh::AztecLinearSolver< T >::AztecLinearSolver(), libMesh::ParallelMesh::clear(), libMesh::Nemesis_IO_Helper::compute_border_node_ids(), libMesh::Nemesis_IO_Helper::construct_nemesis_filename(), libMesh::UnstructuredMesh::create_pid_mesh(), libMesh::DofMap::distribute_dofs(), libMesh::DofMap::distribute_local_dofs_node_major(), libMesh::DofMap::distribute_local_dofs_var_major(), libMesh::EnsightIO::EnsightIO(), libMesh::MeshBase::get_info(), libMesh::EquationSystems::init(), libMesh::SystemSubsetBySubdomain::init(), libMesh::ParmetisPartitioner::initialize(), libMesh::Nemesis_IO_Helper::initialize(), libMesh::MeshTools::libmesh_assert_valid_dof_ids(), libMesh::MeshTools::libmesh_assert_valid_procids< Elem >(), libMesh::MeshTools::libmesh_assert_valid_procids< Node >(), libMesh::MeshTools::libmesh_assert_valid_refinement_flags(), libMesh::MeshBase::n_active_elem_on_proc(), libMesh::MeshBase::n_elem_on_proc(), libMesh::MeshBase::n_nodes_on_proc(), 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::MeshTools::processor_bounding_box(), libMesh::System::project_vector(), libMesh::Nemesis_IO::read(), libMesh::CheckpointIO::read(), libMesh::UnstructuredMesh::read(), libMesh::System::read_parallel_data(), libMesh::System::read_SCALAR_dofs(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::System::read_serialized_vector(), libMesh::Partitioner::repartition(), libMesh::Partitioner::set_node_processor_ids(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::BoundaryInfo::sync(), libMesh::ParallelMesh::update_parallel_id_counts(), libMesh::CheckpointIO::write(), libMesh::GMVIO::write_binary(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::System::write_parallel_data(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bcs(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), and libMesh::XdrIO::write_serialized_nodesets().

93  { return libmesh_cast_int<processor_id_type>(_communicator.size()); }
template<typename T>
bool libMesh::EpetraMatrix< T >::need_full_sparsity_pattern ( ) const
inlinevirtual

The EpetraMatrix needs the full sparsity pattern.

Reimplemented from libMesh::SparseMatrix< T >.

Definition at line 99 of file trilinos_epetra_matrix.h.

100  { return true; }
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 454 of file trilinos_epetra_matrix.C.

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

456 {
457  libmesh_assert (this->initialized());
458  libmesh_assert(this->_mat);
459  libmesh_assert (this->_mat->MyGlobalRow(static_cast<int>(i)));
460  libmesh_assert_greater_equal (i, this->row_start());
461  libmesh_assert_less (i, this->row_stop());
462 
463 
464  int row_length, *row_indices;
465  double *values;
466 
467  _mat->ExtractMyRowView (i-this->row_start(),
468  row_length,
469  values,
470  row_indices);
471 
472  //libMesh::out << "row_length=" << row_length << std::endl;
473 
474  int *index = std::lower_bound (row_indices, row_indices+row_length, j);
475 
476  libmesh_assert_less (*index, row_length);
477  libmesh_assert_equal_to (static_cast<numeric_index_type>(row_indices[*index]), j);
478 
479  //libMesh::out << "val=" << values[*index] << std::endl;
480 
481  return values[*index];
482 }
template<>
void libMesh::SparseMatrix< Complex >::print ( std::ostream &  os,
const bool  sparse 
) const
inherited

Definition at line 100 of file sparse_matrix.C.

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

References libMesh::comm, libMesh::initialized(), libMesh::libmesh_assert(), libMesh::n_processors(), and libMesh::processor_id().

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

222 {
223  parallel_object_only();
224 
225  libmesh_assert (this->initialized());
226 
227  if(!this->_dof_map)
228  {
229  os << std::endl << "Error! Trying to print a matrix with no dof_map set!" << std::endl << std::endl;
230  libmesh_error();
231  }
232 
233  // We'll print the matrix from processor 0 to make sure
234  // it's serialized properly
235  if (this->processor_id() == 0)
236  {
237  libmesh_assert_equal_to (this->_dof_map->first_dof(), 0);
238  for (numeric_index_type i=this->_dof_map->first_dof();
239  i!=this->_dof_map->end_dof(); ++i)
240  {
241  if(sparse)
242  {
243  for (numeric_index_type j=0; j<this->n(); j++)
244  {
245  T c = (*this)(i,j);
246  if (c != static_cast<T>(0.0))
247  {
248  os << i << " " << j << " " << c << std::endl;
249  }
250  }
251  }
252  else
253  {
254  for (numeric_index_type j=0; j<this->n(); j++)
255  os << (*this)(i,j) << " ";
256  os << std::endl;
257  }
258  }
259 
260  std::vector<numeric_index_type> ibuf, jbuf;
261  std::vector<T> cbuf;
262  numeric_index_type currenti = this->_dof_map->end_dof();
263  for (processor_id_type p=1; p < this->n_processors(); ++p)
264  {
265  this->comm().receive(p, ibuf);
266  this->comm().receive(p, jbuf);
267  this->comm().receive(p, cbuf);
268  libmesh_assert_equal_to (ibuf.size(), jbuf.size());
269  libmesh_assert_equal_to (ibuf.size(), cbuf.size());
270 
271  if (ibuf.empty())
272  continue;
273  libmesh_assert_greater_equal (ibuf.front(), currenti);
274  libmesh_assert_greater_equal (ibuf.back(), ibuf.front());
275 
276  std::size_t currentb = 0;
277  for (;currenti <= ibuf.back(); ++currenti)
278  {
279  if(sparse)
280  {
281  for (numeric_index_type j=0; j<this->n(); j++)
282  {
283  if (currentb < ibuf.size() &&
284  ibuf[currentb] == currenti &&
285  jbuf[currentb] == j)
286  {
287  os << currenti << " " << j << " " << cbuf[currentb] << std::endl;
288  currentb++;
289  }
290  }
291  }
292  else
293  {
294  for (numeric_index_type j=0; j<this->n(); j++)
295  {
296  if (currentb < ibuf.size() &&
297  ibuf[currentb] == currenti &&
298  jbuf[currentb] == j)
299  {
300  os << cbuf[currentb] << " ";
301  currentb++;
302  }
303  else
304  os << static_cast<T>(0.0) << " ";
305  }
306  os << std::endl;
307  }
308  }
309  }
310  if(!sparse)
311  {
312  for (; currenti != this->m(); ++currenti)
313  {
314  for (numeric_index_type j=0; j<this->n(); j++)
315  os << static_cast<T>(0.0) << " ";
316  os << std::endl;
317  }
318  }
319  }
320  else
321  {
322  std::vector<numeric_index_type> ibuf, jbuf;
323  std::vector<T> cbuf;
324 
325  // We'll assume each processor has access to entire
326  // matrix rows, so (*this)(i,j) is valid if i is a local index.
327  for (numeric_index_type i=this->_dof_map->first_dof();
328  i!=this->_dof_map->end_dof(); ++i)
329  {
330  for (numeric_index_type j=0; j<this->n(); j++)
331  {
332  T c = (*this)(i,j);
333  if (c != static_cast<T>(0.0))
334  {
335  ibuf.push_back(i);
336  jbuf.push_back(j);
337  cbuf.push_back(c);
338  }
339  }
340  }
341  this->comm().send(0,ibuf);
342  this->comm().send(0,jbuf);
343  this->comm().send(0,cbuf);
344  }
345 }
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().

89 {
91 }
template<typename T >
void libMesh::EpetraMatrix< T >::print_matlab ( const std::string  name = "NULL") const
virtual

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 from libMesh::SparseMatrix< T >.

Definition at line 242 of file trilinos_epetra_matrix.C.

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

243 {
244  libmesh_assert (this->initialized());
245 
246  // libmesh_assert (this->closed());
247  this->close();
248 
249  libmesh_not_implemented();
250 }
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 509 of file trilinos_epetra_matrix.C.

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

510 {
511  libmesh_assert (this->initialized());
513 
514  os << *_mat;
515 }
processor_id_type libMesh::ParallelObject::processor_id ( ) const
inlineinherited
Returns
the rank of this processor in the group.

Definition at line 98 of file parallel_object.h.

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

Referenced by libMesh::MetisPartitioner::_do_partition(), libMesh::EquationSystems::_read_impl(), libMesh::SerialMesh::active_local_elements_begin(), libMesh::ParallelMesh::active_local_elements_begin(), libMesh::SerialMesh::active_local_elements_end(), libMesh::ParallelMesh::active_local_elements_end(), libMesh::SerialMesh::active_local_subdomain_elements_begin(), libMesh::ParallelMesh::active_local_subdomain_elements_begin(), libMesh::SerialMesh::active_local_subdomain_elements_end(), libMesh::ParallelMesh::active_local_subdomain_elements_end(), libMesh::SerialMesh::active_not_local_elements_begin(), libMesh::ParallelMesh::active_not_local_elements_begin(), libMesh::SerialMesh::active_not_local_elements_end(), libMesh::ParallelMesh::active_not_local_elements_end(), libMesh::ParallelMesh::add_elem(), libMesh::DofMap::add_neighbors_to_send_list(), libMesh::ParallelMesh::add_node(), libMesh::UnstructuredMesh::all_second_order(), libMesh::ParmetisPartitioner::assign_partitioning(), 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::DofMap::build_sparsity(), libMesh::ParallelMesh::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::ExodusII_IO_Helper::create(), libMesh::DofMap::distribute_dofs(), libMesh::DofMap::distribute_local_dofs_node_major(), libMesh::DofMap::distribute_local_dofs_var_major(), libMesh::DofMap::end_dof(), libMesh::DofMap::end_old_dof(), libMesh::EnsightIO::EnsightIO(), libMesh::UnstructuredMesh::find_neighbors(), libMesh::DofMap::first_dof(), libMesh::DofMap::first_old_dof(), libMesh::Nemesis_IO_Helper::get_cmap_params(), libMesh::Nemesis_IO_Helper::get_eb_info_global(), libMesh::Nemesis_IO_Helper::get_elem_cmap(), libMesh::Nemesis_IO_Helper::get_elem_map(), libMesh::MeshBase::get_info(), libMesh::Nemesis_IO_Helper::get_init_global(), libMesh::Nemesis_IO_Helper::get_init_info(), libMesh::Nemesis_IO_Helper::get_loadbal_param(), libMesh::Nemesis_IO_Helper::get_node_cmap(), libMesh::Nemesis_IO_Helper::get_node_map(), libMesh::Nemesis_IO_Helper::get_ns_param_global(), libMesh::Nemesis_IO_Helper::get_ss_param_global(), libMesh::MeshFunction::gradient(), libMesh::MeshFunction::hessian(), libMesh::SystemSubsetBySubdomain::init(), libMesh::ParmetisPartitioner::initialize(), libMesh::ExodusII_IO_Helper::initialize(), libMesh::ExodusII_IO_Helper::initialize_discontinuous(), libMesh::ExodusII_IO_Helper::initialize_element_variables(), libMesh::ExodusII_IO_Helper::initialize_global_variables(), libMesh::ExodusII_IO_Helper::initialize_nodal_variables(), libMesh::SparsityPattern::Build::join(), libMesh::DofMap::last_dof(), libMesh::MeshTools::libmesh_assert_valid_procids< Elem >(), libMesh::MeshTools::libmesh_assert_valid_procids< Node >(), libMesh::SerialMesh::local_elements_begin(), libMesh::ParallelMesh::local_elements_begin(), libMesh::SerialMesh::local_elements_end(), libMesh::ParallelMesh::local_elements_end(), libMesh::SerialMesh::local_level_elements_begin(), libMesh::ParallelMesh::local_level_elements_begin(), libMesh::SerialMesh::local_level_elements_end(), libMesh::ParallelMesh::local_level_elements_end(), libMesh::SerialMesh::local_nodes_begin(), libMesh::ParallelMesh::local_nodes_begin(), libMesh::SerialMesh::local_nodes_end(), libMesh::ParallelMesh::local_nodes_end(), libMesh::SerialMesh::local_not_level_elements_begin(), libMesh::ParallelMesh::local_not_level_elements_begin(), libMesh::SerialMesh::local_not_level_elements_end(), libMesh::ParallelMesh::local_not_level_elements_end(), libMesh::MeshRefinement::make_coarsening_compatible(), libMesh::MeshBase::n_active_local_elem(), libMesh::BoundaryInfo::n_boundary_conds(), libMesh::BoundaryInfo::n_edge_conds(), libMesh::DofMap::n_local_dofs(), libMesh::System::n_local_dofs(), libMesh::MeshBase::n_local_elem(), libMesh::MeshBase::n_local_nodes(), libMesh::BoundaryInfo::n_nodeset_conds(), libMesh::SerialMesh::not_local_elements_begin(), libMesh::ParallelMesh::not_local_elements_begin(), libMesh::SerialMesh::not_local_elements_end(), libMesh::ParallelMesh::not_local_elements_end(), libMesh::WeightedPatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::SparsityPattern::Build::operator()(), libMesh::PatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::MeshFunction::operator()(), libMesh::ParallelMesh::ParallelMesh(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::System::project_vector(), 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::Nemesis_IO::read(), libMesh::CheckpointIO::read(), libMesh::XdrIO::read(), libMesh::UnstructuredMesh::read(), libMesh::CheckpointIO::read_connectivity(), libMesh::ExodusII_IO_Helper::read_elem_num_map(), libMesh::System::read_header(), libMesh::System::read_legacy_data(), libMesh::ExodusII_IO_Helper::read_node_num_map(), libMesh::System::read_parallel_data(), libMesh::System::read_SCALAR_dofs(), libMesh::XdrIO::read_serialized_bc_names(), libMesh::XdrIO::read_serialized_bcs(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::System::read_serialized_data(), libMesh::XdrIO::read_serialized_nodes(), libMesh::XdrIO::read_serialized_nodesets(), libMesh::XdrIO::read_serialized_subdomain_names(), libMesh::System::read_serialized_vector(), libMesh::System::read_serialized_vectors(), libMesh::MeshData::read_xdr(), libMesh::Partitioner::set_node_processor_ids(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::LaplaceMeshSmoother::smooth(), libMesh::BoundaryInfo::sync(), libMesh::MeshTools::total_weight(), libMesh::ParallelMesh::update_parallel_id_counts(), libMesh::MeshTools::weight(), libMesh::ExodusII_IO::write(), libMesh::CheckpointIO::write(), libMesh::XdrIO::write(), libMesh::UnstructuredMesh::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_Helper::write_elements_discontinuous(), libMesh::ExodusII_IO::write_global_data(), libMesh::ExodusII_IO_Helper::write_global_values(), libMesh::System::write_header(), libMesh::ExodusII_IO::write_information_records(), libMesh::ExodusII_IO_Helper::write_information_records(), libMesh::ExodusII_IO_Helper::write_nodal_coordinates(), libMesh::ExodusII_IO_Helper::write_nodal_coordinates_discontinuous(), 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::ExodusII_IO_Helper::write_nodesets(), libMesh::Nemesis_IO_Helper::write_nodesets(), libMesh::System::write_parallel_data(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bc_names(), libMesh::XdrIO::write_serialized_bcs(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::System::write_serialized_data(), libMesh::XdrIO::write_serialized_nodes(), libMesh::XdrIO::write_serialized_nodesets(), libMesh::XdrIO::write_serialized_subdomain_names(), libMesh::System::write_serialized_vector(), libMesh::System::write_serialized_vectors(), libMesh::ExodusII_IO_Helper::write_sidesets(), libMesh::Nemesis_IO_Helper::write_sidesets(), libMesh::ExodusII_IO::write_timestep(), and libMesh::ExodusII_IO_Helper::write_timestep().

99  { return libmesh_cast_int<processor_id_type>(_communicator.rank()); }
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  }
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 366 of file trilinos_epetra_matrix.C.

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

367 {
368  libmesh_assert (this->initialized());
370 
371  return static_cast<numeric_index_type>(_map->MinMyGID());
372 }
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 377 of file trilinos_epetra_matrix.C.

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

378 {
379  libmesh_assert (this->initialized());
381 
382  return static_cast<numeric_index_type>(_map->MaxMyGID())+1;
383 }
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 388 of file trilinos_epetra_matrix.C.

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

391 {
392  libmesh_assert (this->initialized());
393 
394  int
395  epetra_i = static_cast<int>(i),
396  epetra_j = static_cast<int>(j);
397 
398  T epetra_value = value;
399 
400  if (_mat->Filled())
401  _mat->ReplaceGlobalValues (epetra_i, 1, &epetra_value, &epetra_j);
402  else
403  _mat->InsertGlobalValues (epetra_i, 1, &epetra_value, &epetra_j);
404 }
template<typename T >
void libMesh::EpetraMatrix< T >::swap ( EpetraMatrix< T > &  m)

Swaps the raw PETSc matrix context pointers.

Definition at line 498 of file trilinos_epetra_matrix.C.

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

499 {
500  std::swap(_mat, m._mat);
501  std::swap(_destroy_mat_on_exit, m._destroy_mat_on_exit);
502 }
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::comm, libMesh::TriangleWrapper::init(), libMesh::initialized(), libMesh::libmesh_assert(), std::min(), and libMesh::processor_id().

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 }
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 191 of file sparse_matrix.C.

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

193 {
194  dest.zero();
195  this->vector_mult_add(dest,arg);
196 }
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 201 of file sparse_matrix.C.

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

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

203 {
204  /* This functionality is actually implemented in the \p
205  NumericVector class. */
206  dest.add_vector(arg,*this);
207 }
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().

195 {
196  libmesh_assert (this->initialized());
197 
198  _mat->Scale(0.0);
199 }
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 212 of file sparse_matrix.C.

213 {
214  /* This functionality isn't implemented or stubbed in every subclass yet */
215  libmesh_not_implemented();
216 }

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 332 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 438 of file sparse_matrix.h.

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

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 137 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 326 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 444 of file sparse_matrix.h.

Referenced by libMesh::PetscMatrix< T >::_get_submatrix(), libMesh::PetscMatrix< T >::get_transpose(), and libMesh::SparseMatrix< Number >::initialized().

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

Holds the distributed Map

Definition at line 321 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 316 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 131 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 126 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 338 of file trilinos_epetra_matrix.h.

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


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

Site Created By: libMesh Developers
Last modified: February 07 2014 16:58:02 UTC

Hosted By:
SourceForge.net Logo