:orphan:
# MatCreateSeqBAIJ
Creates a sparse matrix in `MATSEQAIJ` (block compressed row) 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 MatCreateSeqBAIJ(MPI_Comm comm, PetscInt bs, PetscInt m, PetscInt n, PetscInt nz, const PetscInt nnz[], Mat *A)
```
Collective


## Input Parameters

- ***comm -*** MPI communicator, set to `PETSC_COMM_SELF`
- ***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()`
- ***m -*** number of rows
- ***n -*** number of columns
- ***nz -*** number of nonzero blocks  per block row (same for all rows)
- ***nnz -*** array containing the number of nonzero blocks in the various block rows
(possibly different for each block row) or `NULL`



## Output Parameter

- ***A -*** the matrix



## 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
It is recommended that one use the `MatCreate()`, `MatSetType()` and/or `MatSetFromOptions()`,
MatXXXXSetPreallocation() paradigm instead of this routine directly.
[MatXXXXSetPreallocation() is, for example, `MatSeqAIJSetPreallocation()`]

The number of rows and columns must be divisible by blocksize.

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

A nonzero block is any block that as 1 or more nonzeros in it

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


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

## Level
intermediate

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

## Examples
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/ksp/tutorials/ex74.c.html">src/ksp/ksp/tutorials/ex74.c</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/tao/unconstrained/tutorials/rosenbrock1.c.html">src/tao/unconstrained/tutorials/rosenbrock1.c</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/tao/unconstrained/tutorials/rosenbrock1f.F90.html">src/tao/unconstrained/tutorials/rosenbrock1f.F90</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/tao/unconstrained/tutorials/rosenbrock2.c.html">src/tao/unconstrained/tutorials/rosenbrock2.c</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/tao/unconstrained/tutorials/rosenbrock3.c.html">src/tao/unconstrained/tutorials/rosenbrock3.c</A><BR>


---
[Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/mat/impls/baij/seq/baij.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)  
