# MatSeqAIJSetPreallocation
For good matrix assembly performance the user should preallocate the matrix storage by setting the parameter nz (or the array nnz).  By setting these parameters accurately, performance during matrix assembly can be increased by more than a factor of 50. 
## Synopsis
```
#include "petscmat.h" 
PetscErrorCode MatSeqAIJSetPreallocation(Mat B, PetscInt nz, const PetscInt nnz[])
```
Collective


## Input Parameters

- ***B -*** The matrix
- ***nz -*** number of nonzeros per row (same for all rows)
- ***nnz -*** array containing the number of nonzeros in the various rows
(possibly different for each row) or NULL



## Notes
If nnz is given then nz is ignored

The `MATSEQAIJ` format also called
compressed row storage, is fully compatible with standard Fortran 77
storage.  That is, the stored row and column indices can begin at
either one (as in Fortran) or zero.  See the users' manual for details.

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.  For large problems you MUST preallocate memory or you
will get TERRIBLE performance, see the users' manual chapter on matrices.

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.


## Developer Notes
Use nz of `MAT_SKIP_ALLOCATION` to not allocate any space for the matrix
entries or columns indices

By default, this format uses inodes (identical nodes) when possible, to
improve numerical efficiency of matrix-vector products and solves. We
search for consecutive rows with the same nonzero structure, thereby
reusing matrix information to achieve increased efficiency.


## Options Database Keys

- ***-mat_no_inode  -*** Do not use inodes
- ***-mat_inode_limit <limit> -*** Sets inode limit (max limit=5)





## See Also
 `MatCreate()`, `MatCreateAIJ()`, `MatSetValues()`, `MatSeqAIJSetColumnIndices()`, `MatCreateSeqAIJWithArrays()`, `MatGetInfo()`,
`MatSeqAIJSetTotalPreallocation()`

## Level
intermediate

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

## Examples
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/mat/tutorials/ex15.c.html">src/mat/tutorials/ex15.c.html</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/mat/tutorials/ex15f.F90.html">src/mat/tutorials/ex15f.F90.html</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/mat/tutorials/ex16.c.html">src/mat/tutorials/ex16.c.html</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/mat/tutorials/ex17.c.html">src/mat/tutorials/ex17.c.html</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/mat/tutorials/ex17f.F90.html">src/mat/tutorials/ex17f.F90.html</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/mat/tutorials/ex5k.kokkos.cxx.html">src/mat/tutorials/ex5k.kokkos.cxx.html</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/pc/tutorials/ex3.c.html">src/ksp/pc/tutorials/ex3.c.html</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/ksp/tutorials/ex15f.F90.html">src/ksp/ksp/tutorials/ex15f.F90.html</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/ksp/tutorials/ex18.c.html">src/ksp/ksp/tutorials/ex18.c.html</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/ksp/tutorials/ex2.c.html">src/ksp/ksp/tutorials/ex2.c.html</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/ksp/tutorials/ex3.c.html">src/ksp/ksp/tutorials/ex3.c.html</A><BR>

## Implementations

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


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


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