:orphan:
# KSPFGMRESSetModifyPC
Sets the routine used by `KSPFGMRES` to modify the preconditioner. [](sec_flexibleksp) 
## Synopsis
```
#include "petscksp.h" 
PetscErrorCode KSPFGMRESSetModifyPC(KSP ksp, PetscErrorCode (*fcn)(KSP, PetscInt, PetscInt, PetscReal, void *), void *ctx, PetscErrorCode (*d)(void *))
```
Logically Collective


## Input Parameters

- ***ksp -*** iterative context obtained from `KSPCreate()`
- ***fcn -*** modifypc function
- ***ctx -*** optional context
- ***d -*** optional context destroy routine



## Calling Sequence of `function`
```none
PetscErrorCode fcn(KSP ksp, PetscInt total_its, PetscInt loc_its, PetscReal res_norm, void *ctx);
```

- ***ksp -*** the ksp context being used.
- ***total_its     -*** the total number of FGMRES iterations that have occurred.
- ***loc_its       -*** the number of FGMRES iterations since last restart.
- ***res_norm      -*** the current residual norm.
- ***ctx           -*** optional context variable



## Calling Sequence of `d`
```none
PetscErrorCode d(void *ctx)
```


## Options Database Keys

- ***-ksp_fgmres_modifypcnochange -*** do not change the `PC`
- ***-ksp_fgmres_modifypcksp -*** changes the inner KSP solver tolerances





## Note
Several modifypc routines are predefined, including  `KSPFGMRESModifyPCNoChange()`, and  `KSPFGMRESModifyPCKSP()`


## See Also
 [](ch_ksp), [](sec_flexibleksp), `KSPFGMRES`, `KSPFGMRESModifyPCNoChange()`, `KSPFGMRESModifyPCKSP()`

## Level
intermediate

## Location
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/ksp/impls/gmres/fgmres/modpcf.c.html#KSPFGMRESSetModifyPC">src/ksp/ksp/impls/gmres/fgmres/modpcf.c</A>

## Implementations

<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/ksp/impls/gmres/fgmres/fgmres.c.html#KSPFGMRESSetModifyPC_FGMRES">KSPFGMRESSetModifyPC_FGMRES in src/ksp/ksp/impls/gmres/fgmres/fgmres.c</A><BR>


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


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