MatSeqBAIJSetPreallocation#

Sets the block size and expected nonzeros per row in the matrix. For good matrix assembly performance the user should preallocate the matrix storage by setting the parameter nz (or the array nnz).

Synopsis#

#include "petscmat.h"  
PetscErrorCode MatSeqBAIJSetPreallocation(Mat B, PetscInt bs, PetscInt nz, const PetscInt 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()

  • nz - number of block nonzeros per block row (same for all rows)

  • nnz - array containing the number of block nonzeros in the various block rows (possibly different for each block row) or NULL

Options Database Keys#

  • -mat_no_unroll - uses code that does not unroll the loops in the block calculations (much slower)

  • -mat_block_size - size of the blocks to use

Notes#

If the nnz parameter is given then the nz parameter is ignored

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.

The MATSEQBAIJ format is fully compatible with standard Fortran storage. That is, the stored row and column indices can begin at either one (as in Fortran) or zero.

Specify the preallocated storage with either nz or nnz (not both). Set nz = PETSC_DEFAULT and nnz = NULL for PETSc to control dynamic memory allocation. See Sparse Matrices for details.

See Also#

Matrices, Mat, Sparse Matrices, MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateBAIJ(), MatGetInfo()

Level#

intermediate

Location#

src/mat/impls/baij/seq/baij.c

Examples#

src/snes/tutorials/ex48.c
src/mat/tutorials/ex17.c
src/mat/tutorials/ex17f.F90
src/ksp/ksp/tutorials/ex59.c

Implementations#

MatSeqBAIJSetPreallocation_SeqBAIJ() in src/mat/impls/baij/seq/baij.c


Index of all Mat routines
Table of Contents for all manual pages
Index of all manual pages