:orphan:
# MatSeqSBAIJSetPreallocation
Creates a sparse symmetric matrix in block AIJ (block compressed row) `MATSEQSBAIJ` format.  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 MatSeqSBAIJSetPreallocation(Mat B, PetscInt bs, PetscInt nz, const PetscInt nnz[])
```
Collective


## Input Parameters

- ***B -*** the symmetric 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 upper triangular plus
diagonal portion of each block (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 (only works if a negative bs is passed in





## Notes
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](sec_matsparse) for details.

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.

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


## See Also
 [](ch_matrices), `Mat`, [Sparse Matrices](sec_matsparse), `MATSEQSBAIJ`, `MatCreate()`, `MatCreateSeqAIJ()`, `MatSetValues()`, `MatCreateSBAIJ()`

## Level
intermediate

## Location
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/mat/impls/sbaij/seq/sbaij.c.html#MatSeqSBAIJSetPreallocation">src/mat/impls/sbaij/seq/sbaij.c</A>

## Examples
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/ksp/tutorials/ex2.c.html">src/ksp/ksp/tutorials/ex2.c</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/ksp/tutorials/ex59.c.html">src/ksp/ksp/tutorials/ex59.c</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/mat/tutorials/ex17.c.html">src/mat/tutorials/ex17.c</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/mat/tutorials/ex17f.F90.html">src/mat/tutorials/ex17f.F90</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/snes/tutorials/ex48.c.html">src/snes/tutorials/ex48.c</A><BR>

## Implementations

<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/mat/impls/sbaij/seq/sbaij.c.html#MatSeqSBAIJSetPreallocation_SeqSBAIJ">MatSeqSBAIJSetPreallocation_SeqSBAIJ in src/mat/impls/sbaij/seq/sbaij.c</A><BR>


---
[Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/mat/impls/sbaij/seq/sbaij.c)


[Index of all Mat routines](index.md)  
[Table of Contents for all manual pages](/manualpages/index.md)  
[Index of all manual pages](/manualpages/singleindex.md)  
