# PCBJacobiGetSubKSP
Gets the local `KSP` contexts for all blocks on this processor. 
## Synopsis
```
#include "petscpc.h" 
PetscErrorCode PCBJacobiGetSubKSP(PC pc, PetscInt *n_local, PetscInt *first_local, KSP *ksp[])
```
Not Collective


## Input Parameter

- ***pc -*** the preconditioner context



## Output Parameters

- ***n_local -*** the number of blocks on this processor, or NULL
- ***first_local -*** the global number of the first block on this processor, or NULL
- ***ksp -*** the array of KSP contexts



## Notes
After `PCBJacobiGetSubKSP()` the array of `KSP` contexts is not to be freed.

Currently for some matrix implementations only 1 block per processor
is supported.

You must call `KSPSetUp()` or `PCSetUp()` before calling `PCBJacobiGetSubKSP()`.


## Fortran Usage
You must pass in a `KSP` array that is large enough to contain all the local `KSP`s.

You can call `PCBJacobiGetSubKSP`(pc,nlocal,firstlocal,`PETSC_NULL_KSP`,ierr) to determine how large the
`KSP` array must be.




## See Also
 `PCBJACOBI`, `PCASM`, `PCASMGetSubKSP()`

## Level
advanced

## Location
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/pc/impls/bjacobi/bjacobi.c.html#PCBJacobiGetSubKSP">src/ksp/pc/impls/bjacobi/bjacobi.c</A>

## Examples
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/ksp/tutorials/ex7.c.html">src/ksp/ksp/tutorials/ex7.c.html</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/ksp/tutorials/ex7f.F90.html">src/ksp/ksp/tutorials/ex7f.F90.html</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/snes/tutorials/ex3.c.html">src/snes/tutorials/ex3.c.html</A><BR>

## Implementations

<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/pc/impls/bjacobi/bjacobi.c.html#PCBJacobiGetSubKSP_BJacobi">PCBJacobiGetSubKSP_BJacobi in src/ksp/pc/impls/bjacobi/bjacobi.c</A><BR>


---
[Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/ksp/pc/impls/bjacobi/bjacobi.c)


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