$extrastylesheet
libMesh::EpetraMatrix< T > Class Template Reference

#include <trilinos_epetra_matrix.h>

Inheritance diagram for libMesh::EpetraMatrix< T >:

List of all members.

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

Friends

std::ostream & operator<< (std::ostream &os, const SparseMatrix< T > &m)

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 [protected, inherited]

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

{
  this->clear();
}

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 [inline, protected, virtual, inherited]

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

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

  {
    libmesh_not_implemented();
  }
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().

{
  libmesh_assert (this->initialized());

  int
    epetra_i = static_cast<int>(i),
    epetra_j = static_cast<int>(j);

  T epetra_value = value;

  _mat->SumIntoGlobalValues (epetra_i, 1, &epetra_value, &epetra_j);
}
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().

{
  libmesh_assert (this->initialized());

  // sanity check. but this cannot avoid
  // crash due to incompatible sparsity structure...
  libmesh_assert_equal_to (this->m(), X_in.m());
  libmesh_assert_equal_to (this->n(), X_in.n());

  EpetraMatrix<T>* X = cast_ptr<EpetraMatrix<T>*> (&X_in);

  EpetraExt::MatrixMatrix::Add (*X->_mat, false, a_in, *_mat, 1.);
}
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 
) [virtual, inherited]

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

{
  libmesh_assert_equal_to (dm.m() / brows.size(), dm.n() / bcols.size());

  const numeric_index_type blocksize = cast_int<numeric_index_type>
    (dm.m() / brows.size());

  libmesh_assert_equal_to (dm.m()%blocksize, 0);
  libmesh_assert_equal_to (dm.n()%blocksize, 0);

  std::vector<numeric_index_type> rows, cols;

  rows.reserve(blocksize*brows.size());
  cols.reserve(blocksize*bcols.size());

  for (unsigned int ib=0; ib<brows.size(); ib++)
    {
      numeric_index_type i=brows[ib]*blocksize;

      for (unsigned int v=0; v<blocksize; v++)
        rows.push_back(i++);
    }

  for (unsigned int jb=0; jb<bcols.size(); jb++)
    {
      numeric_index_type j=bcols[jb]*blocksize;

      for (unsigned int v=0; v<blocksize; v++)
        cols.push_back(j++);
    }

  this->add_matrix (dm, rows, cols);
}
template<typename T>
virtual void libMesh::SparseMatrix< T >::add_block_matrix ( const DenseMatrix< T > &  dm,
const std::vector< numeric_index_type > &  dof_indices 
) [inline, virtual, inherited]

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.

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

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

{
  libmesh_assert (this->initialized());

  const numeric_index_type m = dm.m();
  const numeric_index_type n = dm.n();

  libmesh_assert_equal_to (rows.size(), m);
  libmesh_assert_equal_to (cols.size(), n);

  _mat->SumIntoGlobalValues(m, (int *)&rows[0], n, (int *)&cols[0], &dm.get_values()[0]);
}
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.

{
  this->add_matrix (dm, dof_indices, dof_indices);
}
template<typename T>
void libMesh::SparseMatrix< T >::attach_dof_map ( const DofMap dof_map) [inline, inherited]

Get a pointer to the DofMap to use.

Definition at line 112 of file sparse_matrix.h.

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

  { _dof_map = &dof_map; }
template<typename T >
UniquePtr< SparseMatrix< T > > libMesh::SparseMatrix< T >::build ( const Parallel::Communicator comm,
const SolverPackage  solver_package = libMesh::default_solver_package() 
) [static, inherited]

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

Definition at line 135 of file sparse_matrix.C.

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

{
  // Build the appropriate vector
  switch (solver_package)
    {

#ifdef LIBMESH_HAVE_LASPACK
    case LASPACK_SOLVERS:
      return UniquePtr<SparseMatrix<T> >(new LaspackMatrix<T>(comm));
#endif


#ifdef LIBMESH_HAVE_PETSC
    case PETSC_SOLVERS:
      return UniquePtr<SparseMatrix<T> >(new PetscMatrix<T>(comm));
#endif


#ifdef LIBMESH_HAVE_TRILINOS
    case TRILINOS_SOLVERS:
      return UniquePtr<SparseMatrix<T> >(new EpetraMatrix<T>(comm));
#endif


#ifdef LIBMESH_HAVE_EIGEN
    case EIGEN_SOLVERS:
      return UniquePtr<SparseMatrix<T> >(new EigenSparseMatrix<T>(comm));
#endif

    default:
      libmesh_error_msg("ERROR:  Unrecognized solver package: " << solver_package);
    }

  libmesh_error_msg("We'll never get here!");
  return UniquePtr<SparseMatrix<T> >();
}
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().

{
  //  delete _mat;
  //  delete _map;

  this->_is_initialized = false;

  libmesh_assert (!this->initialized());
}
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::AztecLinearSolver< T >::solve().

{
  libmesh_assert(_mat);

  _mat->GlobalAssemble();
}
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 474 of file trilinos_epetra_matrix.C.

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

{
  libmesh_assert (this->initialized());
  libmesh_assert(this->_mat);

  return this->_mat->Filled();
}
const Parallel::Communicator& libMesh::ParallelObject::comm ( ) const [inline, inherited]
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_monitor(), libMesh::__libmesh_petsc_diff_solver_residual(), 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::ImplicitSystem::add_matrix(), 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::MeshBase::cache_elem_dims(), 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(), DMlibMeshSetSystem_libMesh(), libMesh::MeshRefinement::eliminate_unrefined_patches(), libMesh::WeightedPatchRecoveryErrorEstimator::estimate_error(), libMesh::PatchRecoveryErrorEstimator::estimate_error(), libMesh::JumpErrorEstimator::estimate_error(), libMesh::AdjointRefinementEstimator::estimate_error(), libMesh::ExactErrorEstimator::estimate_error(), libMesh::MeshRefinement::flag_elements_by_elem_fraction(), libMesh::MeshRefinement::flag_elements_by_error_fraction(), libMesh::MeshRefinement::flag_elements_by_nelem_target(), libMesh::CondensedEigenSystem::get_eigenpair(), libMesh::ImplicitSystem::get_linear_solver(), libMesh::LocationMap< T >::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::petsc_auto_fieldsplit(), 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::Parallel::sync_node_data_by_element_id(), libMesh::MeshRefinement::test_level_one(), libMesh::MeshRefinement::test_unflagged(), libMesh::MeshTools::total_weight(), libMesh::NameBasedIO::write(), libMesh::CheckpointIO::write(), libMesh::XdrIO::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().

  { 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 [inline, virtual, inherited]

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

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

  {
    this->_get_submatrix(submatrix,
                         rows,
                         cols,
                         false); // false means DO NOT REUSE submatrix
  }

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.

{
  _enable_print_counter = true;
  return;
}
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().

{
  // Convert vector to EpetraVector.
  EpetraVector<T>* epetra_dest = cast_ptr<EpetraVector<T>*>(&dest);

  // Call Epetra function.
  _mat->ExtractDiagonalCopy(*(epetra_dest->vec()));
}
std::string libMesh::ReferenceCounter::get_info ( ) [static, inherited]

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

{
#if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)

  std::ostringstream oss;

  oss << '\n'
      << " ---------------------------------------------------------------------------- \n"
      << "| Reference count information                                                |\n"
      << " ---------------------------------------------------------------------------- \n";

  for (Counts::iterator it = _counts.begin();
       it != _counts.end(); ++it)
    {
      const std::string name(it->first);
      const unsigned int creations    = it->second.first;
      const unsigned int destructions = it->second.second;

      oss << "| " << name << " reference count information:\n"
          << "|  Creations:    " << creations    << '\n'
          << "|  Destructions: " << destructions << '\n';
    }

  oss << " ---------------------------------------------------------------------------- \n";

  return oss.str();

#else

  return "";

#endif
}
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::EpetraMatrix< T >::_mat, and libMesh::EpetraMatrix< T >::_use_transpose.

{
  // Make sure the SparseMatrix passed in is really a EpetraMatrix
  EpetraMatrix<T>& epetra_dest = cast_ref<EpetraMatrix<T>&>(dest);

  if(&epetra_dest != this)
    epetra_dest = *this;

  epetra_dest._use_transpose = !epetra_dest._use_transpose;
  epetra_dest._mat->SetUseTranspose(epetra_dest._use_transpose);
}
void libMesh::ReferenceCounter::increment_constructor_count ( const std::string &  name) [inline, protected, inherited]

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

{
  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
  std::pair<unsigned int, unsigned int>& p = _counts[name];

  p.first++;
}
void libMesh::ReferenceCounter::increment_destructor_count ( const std::string &  name) [inline, protected, inherited]

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

{
  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
  std::pair<unsigned int, unsigned int>& p = _counts[name];

  p.second++;
}
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::initialized(), and libMesh::libmesh_assert().

{
  if ((m==0) || (n==0))
    return;

  {
    // Clear initialized matrices
    if (this->initialized())
      this->clear();

    libmesh_assert (!this->_mat);
    libmesh_assert (!this->_map);

    this->_is_initialized = true;
  }

  // error checking
#ifndef NDEBUG
  {
    libmesh_assert_equal_to (n, m);
    libmesh_assert_equal_to (n_l, m_l);

    numeric_index_type
      summed_m_l = m_l,
      summed_n_l = n_l;

    this->comm().sum (summed_m_l);
    this->comm().sum (summed_n_l);

    libmesh_assert_equal_to (m, summed_m_l);
    libmesh_assert_equal_to (n, summed_n_l);
  }
#endif

  // build a map defining the data distribution
  _map = new Epetra_Map (static_cast<int>(m),
                         m_l,
                         0,
                         Epetra_MpiComm (this->comm().get()));

  libmesh_assert_equal_to (static_cast<numeric_index_type>(_map->NumGlobalPoints()), m);
  libmesh_assert_equal_to (static_cast<numeric_index_type>(_map->MaxAllGID()+1), m);

  _mat = new Epetra_FECrsMatrix (Copy, *_map, nnz + noz);
}
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().

{
  libmesh_assert(this->_dof_map);

  {
    // Clear initialized matrices
    if (this->initialized())
      this->clear();

    this->_is_initialized = true;
  }


  _mat = new Epetra_FECrsMatrix (Copy, *_graph);
}
template<typename T>
virtual bool libMesh::SparseMatrix< T >::initialized ( ) const [inline, virtual, inherited]
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().

{ 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.

{
  libmesh_assert (this->initialized());

  libmesh_assert(_mat);

  return static_cast<Real>(_mat->NormOne());
}
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.

{
  libmesh_assert (this->initialized());


  libmesh_assert(_mat);

  return static_cast<Real>(_mat->NormInf());
}
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().

{
  libmesh_assert (this->initialized());

  return static_cast<numeric_index_type>(_mat->NumGlobalRows());
}
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 297 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().

{ libmesh_assert(_mat); return _mat; }
template<typename T>
const Epetra_FECrsMatrix* libMesh::EpetraMatrix< T >::mat ( ) const [inline]
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().

{
  libmesh_assert (this->initialized());

  return static_cast<numeric_index_type>(_mat->NumGlobalCols());
}
static unsigned int libMesh::ReferenceCounter::n_objects ( ) [inline, static, inherited]

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.

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

  { return _n_objects; }
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::FEMSystem::assembly(), 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::DofMap::local_variable_indices(), 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::NameBasedIO::read(), libMesh::Nemesis_IO::read(), libMesh::CheckpointIO::read(), libMesh::XdrIO::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().

  { return cast_int<processor_id_type>(_communicator.size()); }
template<typename T>
bool libMesh::EpetraMatrix< T >::need_full_sparsity_pattern ( ) const [inline, virtual]

The EpetraMatrix needs the full sparsity pattern.

Reimplemented from libMesh::SparseMatrix< T >.

Definition at line 99 of file trilinos_epetra_matrix.h.

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

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

{
  libmesh_assert (this->initialized());
  libmesh_assert(this->_mat);
  libmesh_assert (this->_mat->MyGlobalRow(static_cast<int>(i)));
  libmesh_assert_greater_equal (i, this->row_start());
  libmesh_assert_less (i, this->row_stop());


  int row_length, *row_indices;
  double *values;

  _mat->ExtractMyRowView (i-this->row_start(),
                          row_length,
                          values,
                          row_indices);

  //libMesh::out << "row_length=" << row_length << std::endl;

  int *index = std::lower_bound (row_indices, row_indices+row_length, j);

  libmesh_assert_less (*index, row_length);
  libmesh_assert_equal_to (static_cast<numeric_index_type>(row_indices[*index]), j);

  //libMesh::out << "val=" << values[*index] << std::endl;

  return values[*index];
}
template<>
void libMesh::SparseMatrix< Complex >::print ( std::ostream &  os,
const bool  sparse 
) const [inherited]

Definition at line 101 of file sparse_matrix.C.

{
  // std::complex<>::operator<<() is defined, but use this form

  if(sparse)
    {
      libmesh_not_implemented();
    }

  os << "Real part:" << std::endl;
  for (numeric_index_type i=0; i<this->m(); i++)
    {
      for (numeric_index_type j=0; j<this->n(); j++)
        os << std::setw(8) << (*this)(i,j).real() << " ";
      os << std::endl;
    }

  os << std::endl << "Imaginary part:" << std::endl;
  for (numeric_index_type i=0; i<this->m(); i++)
    {
      for (numeric_index_type j=0; j<this->n(); j++)
        os << std::setw(8) << (*this)(i,j).imag() << " ";
      os << std::endl;
    }
}
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::initialized(), libMesh::libmesh_assert(), libMesh::n_processors(), and libMesh::processor_id().

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

{
  parallel_object_only();

  libmesh_assert (this->initialized());

  if(!this->_dof_map)
    libmesh_error_msg("Error!  Trying to print a matrix with no dof_map set!");

  // We'll print the matrix from processor 0 to make sure
  // it's serialized properly
  if (this->processor_id() == 0)
    {
      libmesh_assert_equal_to (this->_dof_map->first_dof(), 0);
      for (numeric_index_type i=this->_dof_map->first_dof();
           i!=this->_dof_map->end_dof(); ++i)
        {
          if(sparse)
            {
              for (numeric_index_type j=0; j<this->n(); j++)
                {
                  T c = (*this)(i,j);
                  if (c != static_cast<T>(0.0))
                    {
                      os << i << " " << j << " " << c << std::endl;
                    }
                }
            }
          else
            {
              for (numeric_index_type j=0; j<this->n(); j++)
                os << (*this)(i,j) << " ";
              os << std::endl;
            }
        }

      std::vector<numeric_index_type> ibuf, jbuf;
      std::vector<T> cbuf;
      numeric_index_type currenti = this->_dof_map->end_dof();
      for (processor_id_type p=1; p < this->n_processors(); ++p)
        {
          this->comm().receive(p, ibuf);
          this->comm().receive(p, jbuf);
          this->comm().receive(p, cbuf);
          libmesh_assert_equal_to (ibuf.size(), jbuf.size());
          libmesh_assert_equal_to (ibuf.size(), cbuf.size());

          if (ibuf.empty())
            continue;
          libmesh_assert_greater_equal (ibuf.front(), currenti);
          libmesh_assert_greater_equal (ibuf.back(), ibuf.front());

          std::size_t currentb = 0;
          for (;currenti <= ibuf.back(); ++currenti)
            {
              if(sparse)
                {
                  for (numeric_index_type j=0; j<this->n(); j++)
                    {
                      if (currentb < ibuf.size() &&
                          ibuf[currentb] == currenti &&
                          jbuf[currentb] == j)
                        {
                          os << currenti << " " << j << " " << cbuf[currentb] << std::endl;
                          currentb++;
                        }
                    }
                }
              else
                {
                  for (numeric_index_type j=0; j<this->n(); j++)
                    {
                      if (currentb < ibuf.size() &&
                          ibuf[currentb] == currenti &&
                          jbuf[currentb] == j)
                        {
                          os << cbuf[currentb] << " ";
                          currentb++;
                        }
                      else
                        os << static_cast<T>(0.0) << " ";
                    }
                  os << std::endl;
                }
            }
        }
      if(!sparse)
        {
          for (; currenti != this->m(); ++currenti)
            {
              for (numeric_index_type j=0; j<this->n(); j++)
                os << static_cast<T>(0.0) << " ";
              os << std::endl;
            }
        }
    }
  else
    {
      std::vector<numeric_index_type> ibuf, jbuf;
      std::vector<T> cbuf;

      // We'll assume each processor has access to entire
      // matrix rows, so (*this)(i,j) is valid if i is a local index.
      for (numeric_index_type i=this->_dof_map->first_dof();
           i!=this->_dof_map->end_dof(); ++i)
        {
          for (numeric_index_type j=0; j<this->n(); j++)
            {
              T c = (*this)(i,j);
              if (c != static_cast<T>(0.0))
                {
                  ibuf.push_back(i);
                  jbuf.push_back(j);
                  cbuf.push_back(c);
                }
            }
        }
      this->comm().send(0,ibuf);
      this->comm().send(0,jbuf);
      this->comm().send(0,cbuf);
    }
}
void libMesh::ReferenceCounter::print_info ( std::ostream &  out = libMesh::out) [static, inherited]

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

template<typename T>
virtual void libMesh::SparseMatrix< T >::print_matlab ( const std::string &  = "") const [inline, virtual, inherited]

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

  {
    libmesh_not_implemented();
  }
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 495 of file trilinos_epetra_matrix.C.

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

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::FEMSystem::assembly(), libMesh::ParmetisPartitioner::assign_partitioning(), libMesh::ParallelMesh::assign_unique_ids(), 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::ParallelMesh::delete_elem(), libMesh::ParallelMesh::delete_node(), 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::SerialMesh::facelocal_elements_begin(), libMesh::ParallelMesh::facelocal_elements_begin(), libMesh::SerialMesh::facelocal_elements_end(), libMesh::ParallelMesh::facelocal_elements_end(), libMesh::MeshFunction::find_element(), 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::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::ParallelMesh::insert_elem(), 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::DofMap::local_variable_indices(), libMesh::MeshRefinement::make_coarsening_compatible(), libMesh::MeshBase::n_active_local_elem(), libMesh::BoundaryInfo::n_boundary_conds(), libMesh::BoundaryInfo::n_edge_conds(), libMesh::DofMap::n_local_dofs(), libMesh::System::n_local_dofs(), libMesh::MeshBase::n_local_elem(), libMesh::MeshBase::n_local_nodes(), libMesh::BoundaryInfo::n_nodeset_conds(), libMesh::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::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::NameBasedIO::read(), libMesh::Nemesis_IO::read(), libMesh::CheckpointIO::read(), libMesh::XdrIO::read(), 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::SerialMesh::semilocal_elements_begin(), libMesh::ParallelMesh::semilocal_elements_begin(), libMesh::SerialMesh::semilocal_elements_end(), libMesh::ParallelMesh::semilocal_elements_end(), 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::NameBasedIO::write(), libMesh::ExodusII_IO::write(), libMesh::CheckpointIO::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::System::write_header(), libMesh::ExodusII_IO::write_information_records(), libMesh::ExodusII_IO_Helper::write_information_records(), libMesh::ExodusII_IO_Helper::write_nodal_coordinates(), libMesh::UCDIO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data_discontinuous(), libMesh::ExodusII_IO_Helper::write_nodal_values(), libMesh::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().

  { return 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 [inline, virtual, inherited]

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

  {
    this->_get_submatrix(submatrix,
                         rows,
                         cols,
                         true); // true means REUSE submatrix
  }
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().

{
  libmesh_assert (this->initialized());
  libmesh_assert(_map);

  return static_cast<numeric_index_type>(_map->MinMyGID());
}
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().

{
  libmesh_assert (this->initialized());
  libmesh_assert(_map);

  return static_cast<numeric_index_type>(_map->MaxMyGID())+1;
}
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().

{
  libmesh_assert (this->initialized());

  int
    epetra_i = static_cast<int>(i),
    epetra_j = static_cast<int>(j);

  T epetra_value = value;

  if (_mat->Filled())
    _mat->ReplaceGlobalValues (epetra_i, 1, &epetra_value, &epetra_j);
  else
    _mat->InsertGlobalValues (epetra_i, 1, &epetra_value, &epetra_j);
}
template<typename T >
void libMesh::EpetraMatrix< T >::swap ( EpetraMatrix< T > &  m)

Swaps the raw PETSc matrix context pointers.

Definition at line 484 of file trilinos_epetra_matrix.C.

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

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

{
  // clear data, start over
  this->clear ();

  // big trouble if this fails!
  libmesh_assert(this->_dof_map);

  const numeric_index_type n_rows = sparsity_pattern.size();

  const numeric_index_type m   = this->_dof_map->n_dofs();
  const numeric_index_type n   = m;
  const numeric_index_type n_l = this->_dof_map->n_dofs_on_processor(this->processor_id());
  const numeric_index_type m_l = n_l;

  // error checking
#ifndef NDEBUG
  {
    libmesh_assert_equal_to (n, m);
    libmesh_assert_equal_to (n_l, m_l);

    numeric_index_type
      summed_m_l = m_l,
      summed_n_l = n_l;

    this->comm().sum (summed_m_l);
    this->comm().sum (summed_n_l);

    libmesh_assert_equal_to (m, summed_m_l);
    libmesh_assert_equal_to (n, summed_n_l);
  }
#endif

  // build a map defining the data distribution
  _map = new Epetra_Map (static_cast<int>(m),
                         m_l,
                         0,
                         Epetra_MpiComm (this->comm().get()));

  libmesh_assert_equal_to (static_cast<numeric_index_type>(_map->NumGlobalPoints()), m);
  libmesh_assert_equal_to (static_cast<numeric_index_type>(_map->MaxAllGID()+1), m);

  const std::vector<numeric_index_type>& n_nz = this->_dof_map->get_n_nz();
  const std::vector<numeric_index_type>& n_oz = this->_dof_map->get_n_oz();

  // Make sure the sparsity pattern isn't empty
  libmesh_assert_equal_to (n_nz.size(), n_l);
  libmesh_assert_equal_to (n_oz.size(), n_l);

  // Epetra wants the total number of nonzeros, both local and remote.
  std::vector<int> n_nz_tot;  n_nz_tot.reserve(n_nz.size());

  for (numeric_index_type i=0; i<n_nz.size(); i++)
    n_nz_tot.push_back(std::min(n_nz[i] + n_oz[i], n));

  if (m==0)
    return;

  _graph = new Epetra_CrsGraph(Copy, *_map, &n_nz_tot[0]);

  // Tell the matrix about its structure.  Initialize it
  // to zero.
  for (numeric_index_type i=0; i<n_rows; i++)
    _graph->InsertGlobalIndices(_graph->GRID(i),
                                sparsity_pattern[i].size(),
                                const_cast<int *>((const int *)&sparsity_pattern[i][0]));

  _graph->FillComplete();

  //Initialize the matrix
  libmesh_assert (!this->initialized());
  this->init ();
  libmesh_assert (this->initialized());
}
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::NumericVector< T >::zero().

{
  dest.zero();
  this->vector_mult_add(dest,arg);
}
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::ImplicitSystem::weighted_sensitivity_adjoint_solve().

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

{
  libmesh_assert (this->initialized());

  _mat->Scale(0.0);
}
template<typename T>
void libMesh::SparseMatrix< T >::zero_rows ( std::vector< numeric_index_type > &  rows,
diag_value = 0.0 
) [virtual, inherited]

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.

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

Friends And Related Function Documentation

template<typename T>
std::ostream& operator<< ( std::ostream &  os,
const SparseMatrix< T > &  m 
) [friend, inherited]

Same as the print method above, but allows you to print to a stream in the standard syntax.

template <typename u>=""> friend std::ostream& operator << (std::ostream& os, const SparseMatrix<U>& m);

Obscure C++ note 1: the above syntax, which does not require any prior declaration of operator<<, declares *any* instantiation of SparseMatrix<X> is friend to *any* instantiation of operator<<(ostream&, SparseMatrix<Y>&). It would not happen in practice, but in principle it means that SparseMatrix<Complex> would be friend to operator<<(ostream&, SparseMatrix<Real>).

Obscure C++ note 2: The form below, which requires a previous declaration of the operator<<(stream&, SparseMatrix<T>&) function (see top of this file), means that any instantiation of SparseMatrix<T> is friend to the specialization operator<<(ostream&, SparseMatrix<T>&), but e.g. SparseMatrix<U> is *not* friend to the same function. So this is slightly different to the form above...

This method seems to be the "preferred" technique, see http://www.parashift.com/c++-faq-lite/template-friends.html


Member Data Documentation

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 324 of file trilinos_epetra_matrix.h.

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

template<typename T>
DofMap const* libMesh::SparseMatrix< T >::_dof_map [protected, inherited]

The DofMap object associated with this object.

Definition at line 434 of file sparse_matrix.h.

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

bool libMesh::ReferenceCounter::_enable_print_counter = true [static, protected, inherited]

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 318 of file trilinos_epetra_matrix.h.

template<typename T>
bool libMesh::SparseMatrix< T >::_is_initialized [protected, inherited]

Flag indicating whether or not the matrix has been initialized.

Definition at line 440 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 313 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 308 of file trilinos_epetra_matrix.h.

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

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 [static, protected, inherited]

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 330 of file trilinos_epetra_matrix.h.

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


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