MatMPIBAIJSetPreallocation#
Allocates memory for a sparse parallel matrix in MATMPIBAIJ format (block compressed row).
Synopsis#
#include "petscmat.h"  
PetscErrorCode MatMPIBAIJSetPreallocation(Mat B, PetscInt bs, PetscInt d_nz, const PetscInt d_nnz[], PetscInt o_nz, const PetscInt o_nnz[])
Collective
Input Parameters#
- B - the matrix 
- bs - size of block, the blocks are ALWAYS square. One can use - MatSetBlockSizes()to set a different row and column blocksize but the row blocksize always defines the size of the blocks. The column blocksize sets the blocksize of the vectors obtained with- MatCreateVecs()
- d_nz - number of block nonzeros per block row in diagonal portion of local submatrix (same for all local rows) 
- d_nnz - array containing the number of block nonzeros in the various block rows of the in diagonal portion of the local (possibly different for each block row) or - NULL. If you plan to factor the matrix you must leave room for the diagonal entry and set it even if it is zero.
- o_nz - number of block nonzeros per block row in the off-diagonal portion of local submatrix (same for all local rows). 
- o_nnz - array containing the number of nonzeros in the various block rows of the off-diagonal portion of the local submatrix (possibly different for each block row) or - NULL.
If the *_nnz parameter is given then the *_nz parameter is ignored
Options Database Keys#
- -mat_block_size - size of the blocks to use 
- -mat_use_hash_table - - set hash table factor
Notes#
For good matrix assembly performance
the user should preallocate the matrix storage by setting the parameters
d_nz (or d_nnz) and o_nz (or o_nnz).  By setting these parameters accurately,
performance can be increased by more than a factor of 50.
If PETSC_DECIDE or  PETSC_DETERMINE is used for a particular argument on one processor
than it must be used on all processors that share the object for that argument.
Storage Information#
For a square global matrix we define each processor’s diagonal portion to be its local rows and the corresponding columns (a square submatrix); each processor’s off-diagonal portion encompasses the remainder of the local matrix (a rectangular submatrix).
The user can specify preallocated storage for the diagonal part of
the local submatrix with either d_nz or d_nnz (not both).  Set
d_nz = PETSC_DEFAULT and d_nnz = NULL for PETSc to control dynamic
memory allocation.  Likewise, specify preallocated storage for the
off-diagonal part of the local submatrix with o_nz or o_nnz (not both).
Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In the figure below we depict these three local rows and all columns (0-11).
           0 1 2 3 4 5 6 7 8 9 10 11
          --------------------------
   row 3  |o o o d d d o o o o  o  o
   row 4  |o o o d d d o o o o  o  o
   row 5  |o o o d d d o o o o  o  o
          --------------------------
Thus, any entries in the d locations are stored in the d (diagonal)
submatrix, and any entries in the o locations are stored in the
o (off-diagonal) submatrix.  Note that the d and the o submatrices are
stored simply in the MATSEQBAIJ format for compressed row storage.
Now d_nz should indicate the number of block nonzeros per row in the d matrix,
and o_nz should indicate the number of block nonzeros per row in the o matrix.
In general, for PDE problems in which most nonzeros are near the diagonal,
one expects d_nz >> o_nz.
You can call MatGetInfo() to get information on how effective the preallocation was;
for example the fields mallocs,nz_allocated,nz_used,nz_unneeded;
You can also run with the option -info and look for messages with the string
malloc in them to see if additional memory allocation was needed.
See Also#
Mat, MATMPIBAIJ, MatCreate(), MatCreateSeqBAIJ(), MatSetValues(), MatCreateBAIJ(), MatMPIBAIJSetPreallocationCSR(), PetscSplitOwnership()
Level#
intermediate
Location#
Examples#
src/snes/tutorials/ex48.c
src/mat/tutorials/ex17f.F90
src/mat/tutorials/ex17.c
Implementations#
MatMPIBAIJSetPreallocation_MPIBAIJMKL() in src/mat/impls/baij/mpi/baijmkl/mpibaijmkl.c
MatMPIBAIJSetPreallocation_MPIBAIJ() in src/mat/impls/baij/mpi/mpibaij.c
Index of all Mat routines
Table of Contents for all manual pages
Index of all manual pages