PETScWrappers::MPI::SparseMatrix Class Reference
[PETScWrappersBasic matrices]

Inheritance diagram for PETScWrappers::MPI::SparseMatrix:

Inheritance graph
[legend]

List of all members.

Classes

class  ExcLocalRowsTooLarge
struct  Traits

Public Member Functions

 SparseMatrix ()
 SparseMatrix (const MPI_Comm &communicator, const unsigned int m, const unsigned int n, const unsigned int local_rows, const unsigned int local_columns, const unsigned int n_nonzero_per_row, const bool is_symmetric=false)
 SparseMatrix (const MPI_Comm &communicator, const unsigned int m, const unsigned int n, const unsigned int local_rows, const unsigned int local_columns, const std::vector< unsigned int > &row_lengths, const bool is_symmetric=false)
template<typename SparsityType >
 SparseMatrix (const MPI_Comm &communicator, const SparsityType &sparsity_pattern, const std::vector< unsigned int > &local_rows_per_process, const std::vector< unsigned int > &local_columns_per_process, const unsigned int this_process, const bool preset_nonzero_locations=true)
SparseMatrixoperator= (const double d)
void reinit (const MPI_Comm &communicator, const unsigned int m, const unsigned int n, const unsigned int local_rows, const unsigned int local_columns, const unsigned int n_nonzero_per_row, const bool is_symmetric=false)
void reinit (const MPI_Comm &communicator, const unsigned int m, const unsigned int n, const unsigned int local_rows, const unsigned int local_columns, const std::vector< unsigned int > &row_lengths, const bool is_symmetric=false)
template<typename SparsityType >
void reinit (const MPI_Comm &communicator, const SparsityType &sparsity_pattern, const std::vector< unsigned int > &local_rows_per_process, const std::vector< unsigned int > &local_columns_per_process, const unsigned int this_process, const bool preset_nonzero_locations=true)
virtual const MPI_Comm & get_mpi_communicator () const

Private Member Functions

void do_reinit (const unsigned int m, const unsigned int n, const unsigned int local_rows, const unsigned int local_columns, const unsigned int n_nonzero_per_row, const bool is_symmetric=false)
void do_reinit (const unsigned int m, const unsigned int n, const unsigned int local_rows, const unsigned int local_columns, const std::vector< unsigned int > &row_lengths, const bool is_symmetric=false)
template<typename SparsityType >
void do_reinit (const SparsityType &sparsity_pattern, const std::vector< unsigned int > &local_rows_per_process, const std::vector< unsigned int > &local_columns_per_process, const unsigned int this_process, const bool preset_nonzero_locations)

Private Attributes

MPI_Comm communicator


Detailed Description

Implementation of a parallel sparse matrix class based on PETSC, with rows of the matrix distributed across an MPI network. All the functionality is actually in the base class, except for the calls to generate a parallel sparse matrix. This is possible since PETSc only works on an abstract matrix type and internally distributes to functions that do the actual work depending on the actual matrix type (much like using virtual functions). Only the functions creating a matrix of specific type differ, and are implemented in this particular class.

There are a number of comments on the communication model as well as access to individual elements in the documentation to the parallel vector class. These comments apply here as well.

Partitioning of matrices

PETSc partitions parallel matrices so that each MPI process "owns" a certain number of rows (i.e. only this process stores the respective entries in these rows). The number of rows each process owns has to be passed to the constructors and reinit() functions via the argument local_rows. The individual values passed as local_rows on all the MPI processes of course have to add up to the global number of rows of the matrix.

In addition to this, PETSc also partitions the rectangular chunk of the matrix it owns (i.e. the local_rows times n() elements in the matrix), so that matrix vector multiplications can be performed efficiently. This column-partitioning therefore has to match the partitioning of the vectors with which the matrix is multiplied, just as the row-partitioning has to match the partitioning of destination vectors. This partitioning is passed to the constructors and reinit() functions through the local_columns variable, which again has to add up to the global number of columns in the matrix. The name local_columns may be named inappropriately since it does not reflect that only these columns are stored locally, but it reflects the fact that these are the columns for which the elements of incoming vectors are stored locally.

To make things even more complicated, PETSc needs a very good estimate of the number of elements to be stored in each row to be efficient. Otherwise it spends most of the time with allocating small chunks of memory, a process that can slow down programs to a crawl if it happens to often. As if a good estimate of the number of entries per row isn't even, it even needs to split this as follows: for each row it owns, it needs an estimate for the number of elements in this row that fall into the columns that are set apart for this process (see above), and the number of elements that are in the rest of the columns.

Since in general this information is not readily available, most of the initializing functions of this class assume that all of the number of elements you give as an argument to n_nonzero_per_row or by row_lengths fall into the columns "owned" by this process, and none into the other ones. This is a fair guess for most of the rows, since in a good domain partitioning, nodes only interact with nodes that are within the same subdomain. It does not hold for nodes on the interfaces of subdomain, however, and for the rows corresponding to these nodes, PETSc will have to allocate additional memory, a costly process.

The only way to avoid this is to tell PETSc where the actual entries of the matrix will be. For this, there are constructors and reinit() functions of this class that take a CompressedSparsityPattern object containing all this information. While in the general case it is sufficient if the constructors and reinit() functions know the number of local rows and columns, the functions getting a sparsity pattern also need to know the number of local rows (local_rows_per_process) and columns (local_columns_per_process) for all other processes, in order to compute which parts of the matrix are which. Thus, it is not sufficient to just count the number of degrees of freedom that belong to a particular process, but you have to have the numbers for all processes available at all processes.

Author:
Wolfgang Bangerth, 2004

Constructor & Destructor Documentation

PETScWrappers::MPI::SparseMatrix::SparseMatrix (  ) 

Default constructor. Create an empty matrix.

PETScWrappers::MPI::SparseMatrix::SparseMatrix ( const MPI_Comm &  communicator,
const unsigned int  m,
const unsigned int  n,
const unsigned int  local_rows,
const unsigned int  local_columns,
const unsigned int  n_nonzero_per_row,
const bool  is_symmetric = false 
)

Create a sparse matrix of dimensions m times n, with an initial guess of n_nonzero_per_row nonzero elements per row. PETSc is able to cope with the situation that more than this number of elements are later allocated for a row, but this involves copying data, and is thus expensive.

For the meaning of the local_row and local_columns parameters, see the class documentation.

The is_symmetric flag determines whether we should tell PETSc that the matrix is going to be symmetric (as indicated by the call MatSetOption(mat, MAT_SYMMETRIC). Note that the PETSc documentation states that one cannot form an ILU decomposition of a matrix for which this flag has been set to true, only an ICC. The default value of this flag is false.

PETScWrappers::MPI::SparseMatrix::SparseMatrix ( const MPI_Comm &  communicator,
const unsigned int  m,
const unsigned int  n,
const unsigned int  local_rows,
const unsigned int  local_columns,
const std::vector< unsigned int > &  row_lengths,
const bool  is_symmetric = false 
)

Initialize a rectangular matrix with m rows and n columns. The maximal number of nonzero entries for each row separately is given by the row_lengths array.

For the meaning of the local_row and local_columns parameters, see the class documentation.

Just as for the other constructors: PETSc is able to cope with the situation that more than this number of elements are later allocated for a row, but this involves copying data, and is thus expensive.

The is_symmetric flag determines whether we should tell PETSc that the matrix is going to be symmetric (as indicated by the call MatSetOption(mat, MAT_SYMMETRIC). Note that the PETSc documentation states that one cannot form an ILU decomposition of a matrix for which this flag has been set to true, only an ICC. The default value of this flag is false.

template<typename SparsityType >
PETScWrappers::MPI::SparseMatrix::SparseMatrix ( const MPI_Comm &  communicator,
const SparsityType &  sparsity_pattern,
const std::vector< unsigned int > &  local_rows_per_process,
const std::vector< unsigned int > &  local_columns_per_process,
const unsigned int  this_process,
const bool  preset_nonzero_locations = true 
) [inline]

Initialize using the given sparsity pattern with communication happening over the provided communicator.

For the meaning of the local_rows_per_process and local_columns_per_process parameters, see the class documentation.

Note that PETSc can be very slow if you do not provide it with a good estimate of the lengths of rows. Using the present function is a very efficient way to do this, as it uses the exact number of nonzero entries for each row of the matrix by using the given sparsity pattern argument. If the preset_nonzero_locations flag is true, this function in addition not only sets the correct row sizes up front, but also pre-allocated the correct nonzero entries in the matrix.

PETsc allows to later add additional nonzero entries to a matrix, by simply writing to these elements. However, this will then lead to additional memory allocations which are very inefficient and will greatly slow down your program. It is therefore significantly more efficient to get memory allocation right from the start.


Member Function Documentation

SparseMatrix& PETScWrappers::MPI::SparseMatrix::operator= ( const double  d  ) 

This operator assigns a scalar to a matrix. Since this does usually not make much sense (should we set all matrix entries to this value? Only the nonzero entries of the sparsity pattern?), this operation is only allowed if the actual value to be assigned is zero. This operator only exists to allow for the obvious notation matrix=0, which sets all elements of the matrix to zero, but keep the sparsity pattern previously used.

Reimplemented from PETScWrappers::MatrixBase.

void PETScWrappers::MPI::SparseMatrix::reinit ( const MPI_Comm &  communicator,
const unsigned int  m,
const unsigned int  n,
const unsigned int  local_rows,
const unsigned int  local_columns,
const unsigned int  n_nonzero_per_row,
const bool  is_symmetric = false 
)

Throw away the present matrix and generate one that has the same properties as if it were created by the constructor of this class with the same argument list as the present function.

void PETScWrappers::MPI::SparseMatrix::reinit ( const MPI_Comm &  communicator,
const unsigned int  m,
const unsigned int  n,
const unsigned int  local_rows,
const unsigned int  local_columns,
const std::vector< unsigned int > &  row_lengths,
const bool  is_symmetric = false 
)

Throw away the present matrix and generate one that has the same properties as if it were created by the constructor of this class with the same argument list as the present function.

template<typename SparsityType >
void PETScWrappers::MPI::SparseMatrix::reinit ( const MPI_Comm &  communicator,
const SparsityType &  sparsity_pattern,
const std::vector< unsigned int > &  local_rows_per_process,
const std::vector< unsigned int > &  local_columns_per_process,
const unsigned int  this_process,
const bool  preset_nonzero_locations = true 
) [inline]

Initialize using the given sparsity pattern with communication happening over the provided communicator.

Note that PETSc can be very slow if you do not provide it with a good estimate of the lengths of rows. Using the present function is a very efficient way to do this, as it uses the exact number of nonzero entries for each row of the matrix by using the given sparsity pattern argument. If the preset_nonzero_locations flag is true, this function in addition not only sets the correct row sizes up front, but also pre-allocated the correct nonzero entries in the matrix.

PETsc allows to later add additional nonzero entries to a matrix, by simply writing to these elements. However, this will then lead to additional memory allocations which are very inefficient and will greatly slow down your program. It is therefore significantly more efficient to get memory allocation right from the start.

const MPI_Comm & SparseMatrix< number >::get_mpi_communicator (  )  const [inline, virtual]

Return a reference to the MPI communicator object in use with this matrix.

Implements PETScWrappers::MatrixBase.

References communicator.

void PETScWrappers::MPI::SparseMatrix::do_reinit ( const unsigned int  m,
const unsigned int  n,
const unsigned int  local_rows,
const unsigned int  local_columns,
const unsigned int  n_nonzero_per_row,
const bool  is_symmetric = false 
) [private]

Do the actual work for the respective reinit() function and the matching constructor, i.e. create a matrix. Getting rid of the previous matrix is left to the caller.

void PETScWrappers::MPI::SparseMatrix::do_reinit ( const unsigned int  m,
const unsigned int  n,
const unsigned int  local_rows,
const unsigned int  local_columns,
const std::vector< unsigned int > &  row_lengths,
const bool  is_symmetric = false 
) [private]

Same as previous function.

template<typename SparsityType >
void PETScWrappers::MPI::SparseMatrix::do_reinit ( const SparsityType &  sparsity_pattern,
const std::vector< unsigned int > &  local_rows_per_process,
const std::vector< unsigned int > &  local_columns_per_process,
const unsigned int  this_process,
const bool  preset_nonzero_locations 
) [inline, private]

Same as previous functions.


Member Data Documentation

Copy of the communicator object to be used for this parallel vector.

Referenced by get_mpi_communicator().


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

deal.II documentation generated on Sat Aug 15 16:52:46 2009 by doxygen 1.5.9