:orphan:
# PCSetModifySubMatrices
Sets a user-defined routine for modifying the submatrices that arise within certain subdomain-based preconditioners. The basic submatrices are extracted from the preconditioner matrix as usual; the user can then alter these (for example, to set different boundary conditions for each submatrix) before they are used for the local solves. 
## Synopsis
```
#include "petscksp.h" 
PetscErrorCode PCSetModifySubMatrices(PC pc, PetscErrorCode (*func)(PC, PetscInt, const IS[], const IS[], Mat[], void *), void *ctx)
```
Logically Collective


## Input Parameters

- ***pc -*** the preconditioner context
- ***func -*** routine for modifying the submatrices
- ***ctx -*** optional user-defined context (may be null)



## Calling sequence of `func`
```none
PetscErrorCode func(PC pc, PetscInt nsub, IS *row, IS *col, Mat *submat, void *ctx);
```

- ***pc -*** the preconditioner context
- ***row -*** an array of index sets that contain the global row numbers
that comprise each local submatrix
- ***col -*** an array of index sets that contain the global column numbers
that comprise each local submatrix
- ***submat -*** array of local submatrices
- ***ctx -*** optional user-defined context for private data for the
user-defined func routine (may be null)





## Notes
`PCSetModifySubMatrices()` MUST be called before `KSPSetUp()` and
`KSPSolve()`.

A routine set by `PCSetModifySubMatrices()` is currently called within
the block Jacobi (`PCBJACOBI`) and additive Schwarz (`PCASM`)
preconditioners.  All other preconditioners ignore this routine.


## See Also
 `PC`, `PCBJACOBI`, `PCASM`, `PCModifySubMatrices()`

## Level
advanced

## Location
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/pc/interface/precon.c.html#PCSetModifySubMatrices">src/ksp/pc/interface/precon.c</A>


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


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