libMesh::PetscMatrix< T > Class Template Reference

#include <petsc_linear_solver.h>

Inheritance diagram for libMesh::PetscMatrix< T >:

Public Member Functions

 PetscMatrix (const Parallel::Communicator &comm LIBMESH_CAN_DEFAULT_TO_COMMWORLD)
 
 PetscMatrix (Mat m, const Parallel::Communicator &comm LIBMESH_CAN_DEFAULT_TO_COMMWORLD)
 
 ~PetscMatrix ()
 
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 (const numeric_index_type m, const numeric_index_type n, const numeric_index_type m_l, const numeric_index_type n_l, const std::vector< numeric_index_type > &n_nz, const std::vector< numeric_index_type > &n_oz, const numeric_index_type blocksize=1)
 
void init ()
 
void clear ()
 
void zero ()
 
void zero_rows (std::vector< numeric_index_type > &rows, T diag_value=0.0)
 
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)
 
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 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
 
virtual void get_diagonal (NumericVector< T > &dest) const
 
virtual void get_transpose (SparseMatrix< T > &dest) const
 
void swap (PetscMatrix< T > &)
 
Mat mat ()
 
virtual bool initialized () const
 
void attach_dof_map (const DofMap &dof_map)
 
virtual bool need_full_sparsity_pattern () const
 
virtual void update_sparsity_pattern (const SparsityPattern::Graph &)
 
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 > &submatrix, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols, const bool reuse_submatrix) 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

Mat _mat
 
bool _destroy_mat_on_exit
 

Detailed Description

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

Petsc matrix. Provides a nice interface to the Petsc C-based data structures for parallel, sparse matrices.

Author
Benjamin S. Kirk, 2002

Definition at line 87 of file petsc_linear_solver.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::PetscMatrix< T >::PetscMatrix ( const Parallel::Communicator &comm  LIBMESH_CAN_DEFAULT_TO_COMMWORLD)
explicit

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::PetscMatrix< T >::PetscMatrix ( Mat  m,
const Parallel::Communicator &comm  LIBMESH_CAN_DEFAULT_TO_COMMWORLD 
)
explicit

Constructor. Creates a PetscMatrix assuming you already have a valid Mat object. In this case, m is NOT destroyed by the PetscMatrix 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 PetscMatrix.

template<typename T >
libMesh::PetscMatrix< T >::~PetscMatrix ( )

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

Definition at line 115 of file petsc_matrix.C.

116 {
117  this->clear();
118 }

Member Function Documentation

template<typename T >
void libMesh::PetscMatrix< T >::_get_submatrix ( SparseMatrix< T > &  submatrix,
const std::vector< numeric_index_type > &  rows,
const std::vector< numeric_index_type > &  cols,
const bool  reuse_submatrix 
) const
protectedvirtual

This function either creates or re-initializes a matrix called "submatrix" which is defined by the row and column indices given in the "rows" and "cols" entries. This function is implemented in terms of the MatGetSubMatrix() routine of PETSc. The boolean reuse_submatrix parameter determines whether or not PETSc will treat "submatrix" as one which has already been used (had memory allocated) or as a new matrix.

Reimplemented from libMesh::SparseMatrix< T >.

Definition at line 732 of file petsc_matrix.C.

References libMesh::SparseMatrix< T >::_is_initialized, libMesh::PetscMatrix< T >::_mat, libMesh::SparseMatrix< T >::clear(), libMesh::PetscMatrix< T >::close(), libMesh::comm, ierr, libMesh::SparseMatrix< T >::initialized(), and PETSC_USE_POINTER.

736 {
737  // Can only extract submatrices from closed matrices
738  this->close();
739 
740  // Make sure the SparseMatrix passed in is really a PetscMatrix
741  PetscMatrix<T>* petsc_submatrix = libmesh_cast_ptr<PetscMatrix<T>*>(&submatrix);
742 
743  // If we're not reusing submatrix and submatrix is already initialized
744  // then we need to clear it, otherwise we get a memory leak.
745  if( !reuse_submatrix && submatrix.initialized() )
746  submatrix.clear();
747 
748  // Construct row and column index sets.
750  IS isrow, iscol;
751 
752  ierr = ISCreateLibMesh(this->comm().get(),
753  rows.size(),
754  (PetscInt*) &rows[0],
756  &isrow); LIBMESH_CHKERRABORT(ierr);
757 
758  ierr = ISCreateLibMesh(this->comm().get(),
759  cols.size(),
760  (PetscInt*) &cols[0],
762  &iscol); LIBMESH_CHKERRABORT(ierr);
763 
764  // Extract submatrix
765  ierr = MatGetSubMatrix(_mat,
766  isrow,
767  iscol,
768 #if PETSC_RELEASE_LESS_THAN(3,0,1)
769  PETSC_DECIDE,
770 #endif
771  (reuse_submatrix ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX),
772  &(petsc_submatrix->_mat)); LIBMESH_CHKERRABORT(ierr);
773 
774  // Specify that the new submatrix is initialized and close it.
775  petsc_submatrix->_is_initialized = true;
776  petsc_submatrix->close();
777 
778  // Clean up PETSc data structures
779  ierr = LibMeshISDestroy(&isrow); LIBMESH_CHKERRABORT(ierr);
780  ierr = LibMeshISDestroy(&iscol); LIBMESH_CHKERRABORT(ierr);
781 }
template<typename T >
void libMesh::PetscMatrix< 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 954 of file petsc_matrix.C.

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

957 {
958  libmesh_assert (this->initialized());
959 
961  PetscInt i_val=i, j_val=j;
962 
963  PetscScalar petsc_value = static_cast<PetscScalar>(value);
964  ierr = MatSetValues(_mat, 1, &i_val, 1, &j_val,
965  &petsc_value, ADD_VALUES);
966  LIBMESH_CHKERRABORT(ierr);
967 }
template<typename T >
void libMesh::PetscMatrix< 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} $. Use this with caution, the sparse matrices need to have the same nonzero pattern, otherwise PETSc will crash! 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 985 of file petsc_matrix.C.

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

986 {
987  libmesh_assert (this->initialized());
988 
989  // sanity check. but this cannot avoid
990  // crash due to incompatible sparsity structure...
991  libmesh_assert_equal_to (this->m(), X_in.m());
992  libmesh_assert_equal_to (this->n(), X_in.n());
993 
994  PetscScalar a = static_cast<PetscScalar> (a_in);
995  PetscMatrix<T>* X = libmesh_cast_ptr<PetscMatrix<T>*> (&X_in);
996 
997  libmesh_assert (X);
998 
1000 
1001  // the matrix from which we copy the values has to be assembled/closed
1002  // X->close ();
1003  libmesh_assert(X->closed());
1004 
1005  semiparallel_only();
1006 
1007 // 2.2.x & earlier style
1008 #if PETSC_VERSION_LESS_THAN(2,3,0)
1009 
1010  ierr = MatAXPY(&a, X->_mat, _mat, SAME_NONZERO_PATTERN);
1011  LIBMESH_CHKERRABORT(ierr);
1012 
1013 // 2.3.x & newer
1014 #else
1015 
1016  ierr = MatAXPY(_mat, a, X->_mat, DIFFERENT_NONZERO_PATTERN);
1017  LIBMESH_CHKERRABORT(ierr);
1018 
1019 #endif
1020 }
template<typename T >
void libMesh::PetscMatrix< T >::add_block_matrix ( const DenseMatrix< T > &  dm,
const std::vector< numeric_index_type > &  brows,
const std::vector< numeric_index_type > &  bcols 
)
virtual

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

Definition at line 693 of file petsc_matrix.C.

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

Referenced by libMesh::PetscMatrix< T >::add_block_matrix().

696 {
697  libmesh_assert (this->initialized());
698 
699  const numeric_index_type n_rows = dm.m();
700  const numeric_index_type n_cols = dm.n();
701  const numeric_index_type n_brows = brows.size();
702  const numeric_index_type n_bcols = bcols.size();
703  const numeric_index_type blocksize = n_rows / n_brows;
704 
705  libmesh_assert_equal_to (n_cols / n_bcols, blocksize);
706  libmesh_assert_equal_to (blocksize*n_brows, n_rows);
707  libmesh_assert_equal_to (blocksize*n_bcols, n_cols);
708 
710 
711 #ifndef NDEBUG
712  PetscInt petsc_blocksize;
713  ierr = MatGetBlockSize(_mat, &petsc_blocksize);
714  LIBMESH_CHKERRABORT(ierr);
715  libmesh_assert_equal_to (blocksize, static_cast<numeric_index_type>(petsc_blocksize));
716 #endif
717 
718  // These casts are required for PETSc <= 2.1.5
719  ierr = MatSetValuesBlocked(_mat,
720  n_brows, (PetscInt*) &brows[0],
721  n_bcols, (PetscInt*) &bcols[0],
722  (PetscScalar*) &dm.get_values()[0],
723  ADD_VALUES);
724  LIBMESH_CHKERRABORT(ierr);
725 }
template<typename T>
virtual void libMesh::PetscMatrix< T >::add_block_matrix ( const DenseMatrix< T > &  dm,
const std::vector< numeric_index_type > &  dof_indices 
)
inlinevirtual

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

Reimplemented from libMesh::SparseMatrix< T >.

Definition at line 264 of file petsc_matrix.h.

References libMesh::PetscMatrix< T >::add_block_matrix().

266  { this->add_block_matrix (dm, dof_indices, dof_indices); }
template<typename T >
void libMesh::PetscMatrix< 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 664 of file petsc_matrix.C.

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

667 {
668  libmesh_assert (this->initialized());
669 
670  const numeric_index_type n_rows = dm.m();
671  const numeric_index_type n_cols = dm.n();
672 
673  libmesh_assert_equal_to (rows.size(), n_rows);
674  libmesh_assert_equal_to (cols.size(), n_cols);
675 
677 
678  // These casts are required for PETSc <= 2.1.5
679  ierr = MatSetValues(_mat,
680  n_rows, (PetscInt*) &rows[0],
681  n_cols, (PetscInt*) &cols[0],
682  (PetscScalar*) &dm.get_values()[0],
683  ADD_VALUES);
684  LIBMESH_CHKERRABORT(ierr);
685 }
template<typename T >
void libMesh::PetscMatrix< T >::add_matrix ( const DenseMatrix< T > &  dm,
const std::vector< numeric_index_type > &  dof_indices 
)
virtual

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

Implements libMesh::SparseMatrix< T >.

Definition at line 972 of file petsc_matrix.C.

974 {
975  this->add_matrix (dm, dof_indices, dof_indices);
976 }
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::PetscMatrix< 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 437 of file petsc_matrix.C.

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

438 {
440 
441  if ((this->initialized()) && (this->_destroy_mat_on_exit))
442  {
443  semiparallel_only();
444 
445  ierr = LibMeshMatDestroy (&_mat);
446  LIBMESH_CHKERRABORT(ierr);
447 
448  this->_is_initialized = false;
449  }
450 }
template<typename T >
void libMesh::PetscMatrix< T >::close ( ) const
virtual

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

Implements libMesh::SparseMatrix< T >.

Definition at line 850 of file petsc_matrix.C.

References ierr.

Referenced by libMesh::PetscMatrix< T >::_get_submatrix(), libMesh::PetscLinearSolver< T >::adjoint_solve(), libMesh::PetscMatrix< T >::get_transpose(), libMesh::PetscLinearSolver< T >::solve(), libMesh::SlepcEigenSolver< T >::solve_generalized(), and libMesh::SlepcEigenSolver< T >::solve_standard().

851 {
852  semiparallel_only();
853 
854  // BSK - 1/19/2004
855  // strictly this check should be OK, but it seems to
856  // fail on matrix-free matrices. Do they falsely
857  // state they are assembled? Check with the developers...
858 // if (this->closed())
859 // return;
860 
862 
863  ierr = MatAssemblyBegin (_mat, MAT_FINAL_ASSEMBLY);
864  LIBMESH_CHKERRABORT(ierr);
865  ierr = MatAssemblyEnd (_mat, MAT_FINAL_ASSEMBLY);
866  LIBMESH_CHKERRABORT(ierr);
867 }
template<typename T >
bool libMesh::PetscMatrix< T >::closed ( ) const
virtual

see if Petsc matrix has been closed and fully assembled yet

Implements libMesh::SparseMatrix< T >.

Definition at line 1097 of file petsc_matrix.C.

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

Referenced by libMesh::PetscMatrix< T >::add().

1098 {
1099  libmesh_assert (this->initialized());
1100 
1101  PetscErrorCode ierr=0;
1102  PetscBool assembled;
1103 
1104  ierr = MatAssembled(_mat, &assembled);
1105  LIBMESH_CHKERRABORT(ierr);
1106 
1107  return (assembled == PETSC_TRUE);
1108 }
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::PetscMatrix< T >::get_diagonal ( NumericVector< T > &  dest) const
virtual

Copies the diagonal part of the matrix into dest.

Implements libMesh::SparseMatrix< T >.

Definition at line 786 of file petsc_matrix.C.

References ierr, libMesh::out, and libMesh::PetscVector< T >::vec().

787 {
788  // Make sure the NumericVector passed in is really a PetscVector
789  PetscVector<T>& petsc_dest = libmesh_cast_ref<PetscVector<T>&>(dest);
790 
791  // Call PETSc function.
792 
793 #if PETSC_VERSION_LESS_THAN(2,3,1)
794 
795  libMesh::out << "This method has been developed with PETSc 2.3.1. "
796  << "No one has made it backwards compatible with older "
797  << "versions of PETSc so far; however, it might work "
798  << "without any change with some older version." << std::endl;
799  libmesh_error();
800 
801 #else
802 
803  // Needs a const_cast since PETSc does not work with const.
805  MatGetDiagonal(const_cast<PetscMatrix<T>*>(this)->mat(),petsc_dest.vec()); LIBMESH_CHKERRABORT(ierr);
806 
807 #endif
808 
809 }
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::PetscMatrix< 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 814 of file petsc_matrix.C.

References libMesh::SparseMatrix< T >::_is_initialized, libMesh::PetscMatrix< T >::_mat, libMesh::SparseMatrix< T >::clear(), libMesh::PetscMatrix< T >::close(), and ierr.

815 {
816  // Make sure the SparseMatrix passed in is really a PetscMatrix
817  PetscMatrix<T>& petsc_dest = libmesh_cast_ref<PetscMatrix<T>&>(dest);
818 
819  // If we aren't reusing the matrix then need to clear dest,
820  // otherwise we get a memory leak
821  if(&petsc_dest != this)
822  dest.clear();
823 
825 #if PETSC_VERSION_LESS_THAN(3,0,0)
826  if (&petsc_dest == this)
827  ierr = MatTranspose(_mat,PETSC_NULL);
828  else
829  ierr = MatTranspose(_mat,&petsc_dest._mat);
830  LIBMESH_CHKERRABORT(ierr);
831 #else
832  // FIXME - we can probably use MAT_REUSE_MATRIX in more situations
833  if (&petsc_dest == this)
834  ierr = MatTranspose(_mat,MAT_REUSE_MATRIX,&petsc_dest._mat);
835  else
836  ierr = MatTranspose(_mat,MAT_INITIAL_MATRIX,&petsc_dest._mat);
837  LIBMESH_CHKERRABORT(ierr);
838 #endif
839 
840  // Specify that the transposed matrix is initialized and close it.
841  petsc_dest._is_initialized = true;
842  petsc_dest.close();
843 }
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::PetscMatrix< 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 122 of file petsc_matrix.C.

References libMesh::libMeshPrivateData::_is_initialized, libMesh::comm, ierr, libMesh::initialized(), libMesh::libmesh_ignore(), libMesh::n_local, and libMesh::zero.

Referenced by libMesh::PetscLinearSolver< T >::adjoint_solve(), and libMesh::PetscLinearSolver< T >::solve().

129 {
130  // So compilers don't warn when !LIBMESH_ENABLE_BLOCKED_STORAGE
131  libmesh_ignore(blocksize_in);
132 
133  // Clear initialized matrices
134  if (this->initialized())
135  this->clear();
136 
137  this->_is_initialized = true;
138 
139 
140  PetscErrorCode ierr = 0;
141  PetscInt m_global = static_cast<PetscInt>(m_in);
142  PetscInt n_global = static_cast<PetscInt>(n_in);
143  PetscInt m_local = static_cast<PetscInt>(m_l);
144  PetscInt n_local = static_cast<PetscInt>(n_l);
145  PetscInt n_nz = static_cast<PetscInt>(nnz);
146  PetscInt n_oz = static_cast<PetscInt>(noz);
147 
148  ierr = MatCreate(this->comm().get(), &_mat);
149  LIBMESH_CHKERRABORT(ierr);
150  ierr = MatSetSizes(_mat, m_local, n_local, m_global, n_global);
151  LIBMESH_CHKERRABORT(ierr);
152 
153 #ifdef LIBMESH_ENABLE_BLOCKED_STORAGE
154  PetscInt blocksize = static_cast<PetscInt>(blocksize_in);
155  if (blocksize > 1)
156  {
157  // specified blocksize, bs>1.
158  // double check sizes.
159  libmesh_assert_equal_to (m_local % blocksize, 0);
160  libmesh_assert_equal_to (n_local % blocksize, 0);
161  libmesh_assert_equal_to (m_global % blocksize, 0);
162  libmesh_assert_equal_to (n_global % blocksize, 0);
163  libmesh_assert_equal_to (n_nz % blocksize, 0);
164  libmesh_assert_equal_to (n_oz % blocksize, 0);
165 
166  ierr = MatSetType(_mat, MATBAIJ); // Automatically chooses seqbaij or mpibaij
167  LIBMESH_CHKERRABORT(ierr);
168  ierr = MatSetBlockSize(_mat, blocksize);
169  LIBMESH_CHKERRABORT(ierr);
170  ierr = MatSeqBAIJSetPreallocation(_mat, blocksize, n_nz/blocksize, PETSC_NULL);
171  LIBMESH_CHKERRABORT(ierr);
172  ierr = MatMPIBAIJSetPreallocation(_mat, blocksize,
173  n_nz/blocksize, PETSC_NULL,
174  n_oz/blocksize, PETSC_NULL);
175  LIBMESH_CHKERRABORT(ierr);
176  }
177  else
178 #endif
179  {
180  ierr = MatSetType(_mat, MATAIJ); // Automatically chooses seqaij or mpiaij
181  LIBMESH_CHKERRABORT(ierr);
182  ierr = MatSeqAIJSetPreallocation(_mat, n_nz, PETSC_NULL);
183  LIBMESH_CHKERRABORT(ierr);
184  ierr = MatMPIAIJSetPreallocation(_mat, n_nz, PETSC_NULL, n_oz, PETSC_NULL);
185  LIBMESH_CHKERRABORT(ierr);
186  }
187 
188  // Make it an error for PETSc to allocate new nonzero entries during assembly
189 #if PETSC_VERSION_LESS_THAN(3,0,0)
190  ierr = MatSetOption(_mat, MAT_NEW_NONZERO_ALLOCATION_ERR);
191 #else
192  ierr = MatSetOption(_mat, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE);
193 #endif
194  LIBMESH_CHKERRABORT(ierr);
195 
196  // Is prefix information available somewhere? Perhaps pass in the system name?
197  ierr = MatSetOptionsPrefix(_mat, "");
198  LIBMESH_CHKERRABORT(ierr);
199  ierr = MatSetFromOptions(_mat);
200  LIBMESH_CHKERRABORT(ierr);
201 
202  this->zero ();
203 }
template<typename T >
void libMesh::PetscMatrix< 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 std::vector< numeric_index_type > &  n_nz,
const std::vector< numeric_index_type > &  n_oz,
const numeric_index_type  blocksize = 1 
)

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. noz is the number of off-processor nonzeros per row. Optionally supports a block size, which indicates dense coupled blocks for systems with multiple variables all of the same type.

Definition at line 207 of file petsc_matrix.C.

References libMesh::libMeshPrivateData::_is_initialized, libMesh::comm, ierr, libMesh::initialized(), libMesh::libmesh_ignore(), libMesh::n_local, and libMesh::zero.

214 {
215  // So compilers don't warn when !LIBMESH_ENABLE_BLOCKED_STORAGE
216  libmesh_ignore(blocksize_in);
217 
218  // Clear initialized matrices
219  if (this->initialized())
220  this->clear();
221 
222  this->_is_initialized = true;
223 
224  // Make sure the sparsity pattern isn't empty unless the matrix is 0x0
225  libmesh_assert_equal_to (n_nz.size(), m_l);
226  libmesh_assert_equal_to (n_oz.size(), m_l);
227 
228  PetscErrorCode ierr = 0;
229  PetscInt m_global = static_cast<PetscInt>(m_in);
230  PetscInt n_global = static_cast<PetscInt>(n_in);
231  PetscInt m_local = static_cast<PetscInt>(m_l);
232  PetscInt n_local = static_cast<PetscInt>(n_l);
233 
234  ierr = MatCreate(this->comm().get(), &_mat);
235  LIBMESH_CHKERRABORT(ierr);
236  ierr = MatSetSizes(_mat, m_local, n_local, m_global, n_global);
237  LIBMESH_CHKERRABORT(ierr);
238 
239 #ifdef LIBMESH_ENABLE_BLOCKED_STORAGE
240  PetscInt blocksize = static_cast<PetscInt>(blocksize_in);
241  if (blocksize > 1)
242  {
243  // specified blocksize, bs>1.
244  // double check sizes.
245  libmesh_assert_equal_to (m_local % blocksize, 0);
246  libmesh_assert_equal_to (n_local % blocksize, 0);
247  libmesh_assert_equal_to (m_global % blocksize, 0);
248  libmesh_assert_equal_to (n_global % blocksize, 0);
249 
250  ierr = MatSetType(_mat, MATBAIJ); // Automatically chooses seqbaij or mpibaij
251  LIBMESH_CHKERRABORT(ierr);
252  ierr = MatSetBlockSize(_mat, blocksize);
253  LIBMESH_CHKERRABORT(ierr);
254 
255  // transform the per-entry n_nz and n_oz arrays into their block counterparts.
256  std::vector<numeric_index_type> b_n_nz, b_n_oz;
257 
258  transform_preallocation_arrays (blocksize,
259  n_nz, n_oz,
260  b_n_nz, b_n_oz);
261 
262  ierr = MatSeqBAIJSetPreallocation(_mat, blocksize, 0, (PetscInt*)(b_n_nz.empty()?NULL:&b_n_nz[0]));
263  LIBMESH_CHKERRABORT(ierr);
264 
265  ierr = MatMPIBAIJSetPreallocation(_mat, blocksize,
266  0, (PetscInt*)(b_n_nz.empty()?NULL:&b_n_nz[0]),
267  0, (PetscInt*)(b_n_oz.empty()?NULL:&b_n_oz[0]));
268  LIBMESH_CHKERRABORT(ierr);
269  }
270  else
271 #endif
272  {
273  ierr = MatSetType(_mat, MATAIJ); // Automatically chooses seqaij or mpiaij
274  LIBMESH_CHKERRABORT(ierr);
275  ierr = MatSeqAIJSetPreallocation(_mat, 0, (PetscInt*)(n_nz.empty()?NULL:&n_nz[0]));
276  LIBMESH_CHKERRABORT(ierr);
277  ierr = MatMPIAIJSetPreallocation(_mat,
278  0, (PetscInt*)(n_nz.empty()?NULL:&n_nz[0]),
279  0, (PetscInt*)(n_oz.empty()?NULL:&n_oz[0]));
280  LIBMESH_CHKERRABORT(ierr);
281  }
282 
283  // Is prefix information available somewhere? Perhaps pass in the system name?
284  ierr = MatSetOptionsPrefix(_mat, "");
285  LIBMESH_CHKERRABORT(ierr);
286  ierr = MatSetFromOptions(_mat);
287  LIBMESH_CHKERRABORT(ierr);
288 
289 
290  this->zero();
291 }
template<typename T >
void libMesh::PetscMatrix< T >::init ( )
virtual

Initialize using sparsity structure computed by dof_map.

Implements libMesh::SparseMatrix< T >.

Definition at line 295 of file petsc_matrix.C.

References libMesh::libMeshPrivateData::_is_initialized, libMesh::comm, ierr, libMesh::initialized(), libMesh::libmesh_assert(), libMesh::n_local, libMesh::processor_id(), and libMesh::zero.

296 {
297  libmesh_assert(this->_dof_map);
298 
299  // Clear initialized matrices
300  if (this->initialized())
301  this->clear();
302 
303  this->_is_initialized = true;
304 
305 
306  const numeric_index_type my_m = this->_dof_map->n_dofs();
307  const numeric_index_type my_n = my_m;
308  const numeric_index_type n_l = this->_dof_map->n_dofs_on_processor(this->processor_id());
309  const numeric_index_type m_l = n_l;
310 
311 
312  const std::vector<numeric_index_type>& n_nz = this->_dof_map->get_n_nz();
313  const std::vector<numeric_index_type>& n_oz = this->_dof_map->get_n_oz();
314 
315  // Make sure the sparsity pattern isn't empty unless the matrix is 0x0
316  libmesh_assert_equal_to (n_nz.size(), m_l);
317  libmesh_assert_equal_to (n_oz.size(), m_l);
318 
319  PetscErrorCode ierr = 0;
320  PetscInt m_global = static_cast<PetscInt>(my_m);
321  PetscInt n_global = static_cast<PetscInt>(my_n);
322  PetscInt m_local = static_cast<PetscInt>(m_l);
323  PetscInt n_local = static_cast<PetscInt>(n_l);
324 
325  ierr = MatCreate(this->comm().get(), &_mat);
326  LIBMESH_CHKERRABORT(ierr);
327  ierr = MatSetSizes(_mat, m_local, n_local, m_global, n_global);
328  LIBMESH_CHKERRABORT(ierr);
329 
330 #ifdef LIBMESH_ENABLE_BLOCKED_STORAGE
331  PetscInt blocksize = static_cast<PetscInt>(this->_dof_map->block_size());
332  if (blocksize > 1)
333  {
334  // specified blocksize, bs>1.
335  // double check sizes.
336  libmesh_assert_equal_to (m_local % blocksize, 0);
337  libmesh_assert_equal_to (n_local % blocksize, 0);
338  libmesh_assert_equal_to (m_global % blocksize, 0);
339  libmesh_assert_equal_to (n_global % blocksize, 0);
340 
341  ierr = MatSetType(_mat, MATBAIJ); // Automatically chooses seqbaij or mpibaij
342  LIBMESH_CHKERRABORT(ierr);
343  ierr = MatSetBlockSize(_mat, blocksize);
344  LIBMESH_CHKERRABORT(ierr);
345 
346  // transform the per-entry n_nz and n_oz arrays into their block counterparts.
347  std::vector<numeric_index_type> b_n_nz, b_n_oz;
348 
349  transform_preallocation_arrays (blocksize,
350  n_nz, n_oz,
351  b_n_nz, b_n_oz);
352 
353  ierr = MatSeqBAIJSetPreallocation(_mat, blocksize, 0, (PetscInt*)(b_n_nz.empty()?NULL:&b_n_nz[0]));
354  LIBMESH_CHKERRABORT(ierr);
355 
356  ierr = MatMPIBAIJSetPreallocation(_mat, blocksize,
357  0, (PetscInt*)(b_n_nz.empty()?NULL:&b_n_nz[0]),
358  0, (PetscInt*)(b_n_oz.empty()?NULL:&b_n_oz[0]));
359  LIBMESH_CHKERRABORT(ierr);
360  }
361  else
362 #endif
363  {
364  // no block storage case
365  ierr = MatSetType(_mat, MATAIJ); // Automatically chooses seqaij or mpiaij
366  LIBMESH_CHKERRABORT(ierr);
367 
368  ierr = MatSeqAIJSetPreallocation(_mat, 0, (PetscInt*)(n_nz.empty()?NULL:&n_nz[0]));
369  LIBMESH_CHKERRABORT(ierr);
370  ierr = MatMPIAIJSetPreallocation(_mat,
371  0, (PetscInt*)(n_nz.empty()?NULL:&n_nz[0]),
372  0, (PetscInt*)(n_oz.empty()?NULL:&n_oz[0]));
373  LIBMESH_CHKERRABORT(ierr);
374  }
375 
376  // Is prefix information available somewhere? Perhaps pass in the system name?
377  ierr = MatSetOptionsPrefix(_mat, "");
378  LIBMESH_CHKERRABORT(ierr);
379  ierr = MatSetFromOptions(_mat);
380  LIBMESH_CHKERRABORT(ierr);
381 
382  this->zero();
383 }
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::PetscMatrix< 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 455 of file petsc_matrix.C.

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

456 {
457  libmesh_assert (this->initialized());
458 
459  semiparallel_only();
460 
462  PetscReal petsc_value;
463  Real value;
464 
465  libmesh_assert (this->closed());
466 
467  ierr = MatNorm(_mat, NORM_1, &petsc_value);
468  LIBMESH_CHKERRABORT(ierr);
469 
470  value = static_cast<Real>(petsc_value);
471 
472  return value;
473 }
template<typename T >
Real libMesh::PetscMatrix< 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 478 of file petsc_matrix.C.

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

479 {
480  libmesh_assert (this->initialized());
481 
482  semiparallel_only();
483 
485  PetscReal petsc_value;
486  Real value;
487 
488  libmesh_assert (this->closed());
489 
490  ierr = MatNorm(_mat, NORM_INFINITY, &petsc_value);
491  LIBMESH_CHKERRABORT(ierr);
492 
493  value = static_cast<Real>(petsc_value);
494 
495  return value;
496 }
template<typename T >
numeric_index_type libMesh::PetscMatrix< 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 872 of file petsc_matrix.C.

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

873 {
874  libmesh_assert (this->initialized());
875 
876  PetscInt petsc_m=0, petsc_n=0;
878 
879  ierr = MatGetSize (_mat, &petsc_m, &petsc_n);
880  LIBMESH_CHKERRABORT(ierr);
881 
882  return static_cast<numeric_index_type>(petsc_m);
883 }
template<typename T>
Mat libMesh::PetscMatrix< T >::mat ( )
inline
template<typename T >
numeric_index_type libMesh::PetscMatrix< 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 888 of file petsc_matrix.C.

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

889 {
890  libmesh_assert (this->initialized());
891 
892  PetscInt petsc_m=0, petsc_n=0;
894 
895  ierr = MatGetSize (_mat, &petsc_m, &petsc_n);
896  LIBMESH_CHKERRABORT(ierr);
897 
898  return static_cast<numeric_index_type>(petsc_n);
899 }
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>
virtual bool libMesh::SparseMatrix< T >::need_full_sparsity_pattern ( ) const
inlinevirtualinherited

returns true if this sparse matrix format needs to be fed the graph of the sparse matrix. This is true in the case of the LaspackMatrix, but not for the PetscMatrix. In the case where the full graph is not required we can efficiently approximate it to provide a good estimate of the required size of the sparse matrix.

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

Definition at line 122 of file sparse_matrix.h.

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

123  { return false; }
template<typename T >
T libMesh::PetscMatrix< 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 1026 of file petsc_matrix.C.

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

1028 {
1029  libmesh_assert (this->initialized());
1030 
1031 #if PETSC_VERSION_LESS_THAN(2,2,1)
1032 
1033  // PETSc 2.2.0 & older
1034  PetscScalar *petsc_row;
1035  int* petsc_cols;
1036 
1037 #else
1038 
1039  // PETSc 2.2.1 & newer
1040  const PetscScalar *petsc_row;
1041  const PetscInt *petsc_cols;
1042 
1043 #endif
1044 
1045 
1046  // If the entry is not in the sparse matrix, it is 0.
1047  T value=0.;
1048 
1050  ierr=0;
1051  PetscInt
1052  ncols=0,
1053  i_val=static_cast<PetscInt>(i_in),
1054  j_val=static_cast<PetscInt>(j_in);
1055 
1056 
1057  // the matrix needs to be closed for this to work
1058  // this->close();
1059  // but closing it is a semiparallel operation; we want operator()
1060  // to run on one processor.
1061  libmesh_assert(this->closed());
1062 
1063  ierr = MatGetRow(_mat, i_val, &ncols, &petsc_cols, &petsc_row);
1064  LIBMESH_CHKERRABORT(ierr);
1065 
1066  // Perform a binary search to find the contiguous index in
1067  // petsc_cols (resp. petsc_row) corresponding to global index j_val
1068  std::pair<const PetscInt*, const PetscInt*> p =
1069  std::equal_range (&petsc_cols[0], &petsc_cols[0] + ncols, j_val);
1070 
1071  // Found an entry for j_val
1072  if (p.first != p.second)
1073  {
1074  // The entry in the contiguous row corresponding
1075  // to the j_val column of interest
1076  const std::size_t j =
1077  std::distance (const_cast<PetscInt*>(&petsc_cols[0]),
1078  const_cast<PetscInt*>(p.first));
1079 
1080  libmesh_assert_less (static_cast<PetscInt>(j), ncols);
1081  libmesh_assert_equal_to (petsc_cols[j], j_val);
1082 
1083  value = static_cast<T> (petsc_row[j]);
1084  }
1085 
1086  ierr = MatRestoreRow(_mat, i_val,
1087  &ncols, &petsc_cols, &petsc_row);
1088  LIBMESH_CHKERRABORT(ierr);
1089 
1090  return value;
1091 }
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::PetscMatrix< 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.

Create an ASCII file containing the matrix if a filename was provided.

Otherwise the matrix will be dumped to the screen.

Destroy the viewer.

Reimplemented from libMesh::SparseMatrix< T >.

Definition at line 501 of file petsc_matrix.C.

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

502 {
503  libmesh_assert (this->initialized());
504 
505  semiparallel_only();
506 
507  // libmesh_assert (this->closed());
508  this->close();
509 
511  PetscViewer petsc_viewer;
512 
513 
514  ierr = PetscViewerCreate (this->comm().get(),
515  &petsc_viewer);
516  LIBMESH_CHKERRABORT(ierr);
517 
522  if (name != "NULL")
523  {
524  ierr = PetscViewerASCIIOpen( this->comm().get(),
525  name.c_str(),
526  &petsc_viewer);
527  LIBMESH_CHKERRABORT(ierr);
528 
529  ierr = PetscViewerSetFormat (petsc_viewer,
530  PETSC_VIEWER_ASCII_MATLAB);
531  LIBMESH_CHKERRABORT(ierr);
532 
533  ierr = MatView (_mat, petsc_viewer);
534  LIBMESH_CHKERRABORT(ierr);
535  }
536 
540  else
541  {
542  ierr = PetscViewerSetFormat (PETSC_VIEWER_STDOUT_WORLD,
543  PETSC_VIEWER_ASCII_MATLAB);
544  LIBMESH_CHKERRABORT(ierr);
545 
546  ierr = MatView (_mat, PETSC_VIEWER_STDOUT_WORLD);
547  LIBMESH_CHKERRABORT(ierr);
548  }
549 
550 
554  ierr = LibMeshPetscViewerDestroy (&petsc_viewer);
555  LIBMESH_CHKERRABORT(ierr);
556 }
template<typename T >
void libMesh::PetscMatrix< T >::print_personal ( std::ostream &  os = libMesh::out) const
virtual

Print the contents of the matrix to the screen with the PETSc viewer. This function only allows printing to standard out, this is because we have limited ourselves to one PETSc implementation for writing.

Implements libMesh::SparseMatrix< T >.

Definition at line 563 of file petsc_matrix.C.

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

564 {
565  libmesh_assert (this->initialized());
566 
567  // Routine must be called in parallel on parallel matrices
568  // and serial on serial matrices.
569  semiparallel_only();
570 
571 // #ifndef NDEBUG
572 // if (os != std::cout)
573 // libMesh::err << "Warning! PETSc can only print to std::cout!" << std::endl;
574 // #endif
575 
576  // Matrix must be in an assembled state to be printed
577  this->close();
578 
580 
581  // Print to screen if ostream is stdout
582  if (os == std::cout)
583  {
584  ierr = MatView(_mat, PETSC_VIEWER_STDOUT_SELF);
585  LIBMESH_CHKERRABORT(ierr);
586  }
587 
588  // Otherwise, print to the requested file, in a roundabout way...
589  else
590  {
591  // We will create a temporary filename, and file, for PETSc to
592  // write to.
593  std::string temp_filename;
594 
595  {
596  // Template for temporary filename
597  char c[] = "temp_petsc_matrix.XXXXXX";
598 
599  // Generate temporary, unique filename only on processor 0. We will
600  // use this filename for PetscViewerASCIIOpen, before copying it into
601  // the user's stream
602  if (this->processor_id() == 0)
603  {
604  int fd = mkstemp(c);
605 
606  // Check to see that mkstemp did not fail.
607  if (fd == -1)
608  libmesh_error();
609 
610  // mkstemp returns a file descriptor for an open file,
611  // so let's close it before we hand it to PETSc!
612  ::close (fd);
613  }
614 
615  // Store temporary filename as string, makes it easier to broadcast
616  temp_filename = c;
617  }
618 
619  // Now broadcast the filename from processor 0 to all processors.
620  this->comm().broadcast(temp_filename);
621 
622  // PetscViewer object for passing to MatView
623  PetscViewer petsc_viewer;
624 
625  // This PETSc function only takes a string and handles the opening/closing
626  // of the file internally. Since print_personal() takes a reference to
627  // an ostream, we have to do an extra step... print_personal() should probably
628  // have a version that takes a string to get rid of this problem.
629  ierr = PetscViewerASCIIOpen( this->comm().get(),
630  temp_filename.c_str(),
631  &petsc_viewer);
632  LIBMESH_CHKERRABORT(ierr);
633 
634  // Probably don't need to set the format if it's default...
635  // ierr = PetscViewerSetFormat (petsc_viewer,
636  // PETSC_VIEWER_DEFAULT);
637  // LIBMESH_CHKERRABORT(ierr);
638 
639  // Finally print the matrix using the viewer
640  ierr = MatView (_mat, petsc_viewer);
641  LIBMESH_CHKERRABORT(ierr);
642 
643  if (this->processor_id() == 0)
644  {
645  // Now the inefficient bit: open temp_filename as an ostream and copy the contents
646  // into the user's desired ostream. We can't just do a direct file copy, we don't even have the filename!
647  std::ifstream input_stream(temp_filename.c_str());
648  os << input_stream.rdbuf(); // The "most elegant" way to copy one stream into another.
649  // os.close(); // close not defined in ostream
650 
651  // Now remove the temporary file
652  input_stream.close();
653  std::remove(temp_filename.c_str());
654  }
655  }
656 }
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::PetscMatrix< 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 904 of file petsc_matrix.C.

References ierr, libMesh::initialized(), libMesh::libmesh_assert(), and libMesh::MacroFunctions::stop().

905 {
906  libmesh_assert (this->initialized());
907 
908  PetscInt start=0, stop=0;
910 
911  ierr = MatGetOwnershipRange(_mat, &start, &stop);
912  LIBMESH_CHKERRABORT(ierr);
913 
914  return static_cast<numeric_index_type>(start);
915 }
template<typename T >
numeric_index_type libMesh::PetscMatrix< 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 920 of file petsc_matrix.C.

References ierr, libMesh::initialized(), libMesh::libmesh_assert(), and libMesh::MacroFunctions::stop().

921 {
922  libmesh_assert (this->initialized());
923 
924  PetscInt start=0, stop=0;
926 
927  ierr = MatGetOwnershipRange(_mat, &start, &stop);
928  LIBMESH_CHKERRABORT(ierr);
929 
930  return static_cast<numeric_index_type>(stop);
931 }
template<typename T >
void libMesh::PetscMatrix< 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 936 of file petsc_matrix.C.

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

939 {
940  libmesh_assert (this->initialized());
941 
943  PetscInt i_val=i, j_val=j;
944 
945  PetscScalar petsc_value = static_cast<PetscScalar>(value);
946  ierr = MatSetValues(_mat, 1, &i_val, 1, &j_val,
947  &petsc_value, INSERT_VALUES);
948  LIBMESH_CHKERRABORT(ierr);
949 }
template<typename T >
void libMesh::PetscMatrix< T >::swap ( PetscMatrix< T > &  m_in)

Swaps the raw PETSc matrix context pointers.

Definition at line 1113 of file petsc_matrix.C.

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

Referenced by libMesh::__libmesh_petsc_diff_solver_jacobian().

1114 {
1115  std::swap(_mat, m_in._mat);
1117 }
template<typename T>
virtual void libMesh::SparseMatrix< T >::update_sparsity_pattern ( const SparsityPattern::Graph )
inlinevirtualinherited

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

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

Definition at line 130 of file sparse_matrix.h.

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

130 {}
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::PetscMatrix< T >::zero ( )
virtual

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

Implements libMesh::SparseMatrix< T >.

Definition at line 388 of file petsc_matrix.C.

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

389 {
390  libmesh_assert (this->initialized());
391 
392  semiparallel_only();
393 
395 
396  PetscInt m_l, n_l;
397 
398  ierr = MatGetLocalSize(_mat,&m_l,&n_l);
399  LIBMESH_CHKERRABORT(ierr);
400 
401  if (n_l)
402  {
403  ierr = MatZeroEntries(_mat);
404  LIBMESH_CHKERRABORT(ierr);
405  }
406 }
template<typename T >
void libMesh::PetscMatrix< T >::zero_rows ( std::vector< numeric_index_type > &  rows,
diag_value = 0.0 
)
virtual

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

Reimplemented from libMesh::SparseMatrix< T >.

Definition at line 409 of file petsc_matrix.C.

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

410 {
411  libmesh_assert (this->initialized());
412 
413  semiparallel_only();
414 
416 
417 #if PETSC_RELEASE_LESS_THAN(3,1,1)
418  if(!rows.empty())
419  ierr = MatZeroRows(_mat, rows.size(), (PetscInt*)&rows[0], diag_value);
420  else
421  ierr = MatZeroRows(_mat, 0, PETSC_NULL, diag_value);
422 #else
423  // As of petsc-dev at the time of 3.1.0, MatZeroRows now takes two additional
424  // optional arguments. The optional arguments (x,b) can be used to specify the
425  // solutions for the zeroed rows (x) and right hand side (b) to update.
426  // Could be useful for setting boundary conditions...
427  if(!rows.empty())
428  ierr = MatZeroRows(_mat, rows.size(), (PetscInt*)&rows[0], diag_value, PETSC_NULL, PETSC_NULL);
429  else
430  ierr = MatZeroRows(_mat, 0, PETSC_NULL, diag_value, PETSC_NULL, PETSC_NULL);
431 #endif
432 
433  LIBMESH_CHKERRABORT(ierr);
434 }

Member Data Documentation

ReferenceCounter::Counts libMesh::ReferenceCounter::_counts
staticprotectedinherited
template<typename T>
bool libMesh::PetscMatrix< 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 391 of file petsc_matrix.h.

Referenced by libMesh::PetscMatrix< 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>
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().

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


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

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

Hosted By:
SourceForge.net Logo