:orphan:
# MatSeqBAIJSetPreallocationCSR
Creates a sparse sequential matrix in `MATSEQBAIJ` format using the given nonzero structure and (optional) numerical values 
## Synopsis
```
#include "petscmat.h"  
PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat B, PetscInt bs, const PetscInt i[], const PetscInt j[], const PetscScalar v[])
```
Collective


## Input Parameters

- ***B -*** the matrix
- ***bs -*** the blocksize
- ***i -*** the indices into `j` for the start of each local row (starts with zero)
- ***j -*** the column indices for each local row (starts with zero) these must be sorted for each row
- ***v -*** optional values in the matrix





## Notes
The order of the entries in values is specified by the `MatOption` `MAT_ROW_ORIENTED`.  For example, C programs
may want to use the default `MAT_ROW_ORIENTED` of `PETSC_TRUE` and use an array v[nnz][bs][bs] where the second index is
over rows within a block and the last index is over columns within a block row.  Fortran programs will likely set
`MAT_ROW_ORIENTED` of `PETSC_FALSE` and use a Fortran array v(bs,bs,nnz) in which the first index is over rows within a
block column and the second index is over columns within a block.

Though this routine has Preallocation() in the name it also sets the exact nonzero locations of the matrix entries and usually the numerical values as well


## See Also
 [](ch_matrices), `Mat`, `MatCreate()`, `MatCreateSeqBAIJ()`, `MatSetValues()`, `MatSeqBAIJSetPreallocation()`, `MATSEQBAIJ`

## Level
advanced

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

## Implementations

<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/mat/impls/baij/seq/baij.c.html#MatSeqBAIJSetPreallocationCSR_SeqBAIJ">MatSeqBAIJSetPreallocationCSR_SeqBAIJ in src/mat/impls/baij/seq/baij.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)  
