:orphan: # KSPPIPEGCRSetModifyPC Sets the routine used by `KSPPIPEGCR` to modify the preconditioner at each iteration ## Synopsis ``` #include "petscksp.h" PetscErrorCode KSPPIPEGCRSetModifyPC(KSP ksp, PetscErrorCode (*function)(KSP, PetscInt, PetscReal, void *), void *data, PetscErrorCode (*destroy)(void *)) ``` Logically Collective ## Input Parameters + ksp - iterative context obtained from KSPCreate() . function - user defined function to modify the preconditioner . ctx - user provided context for the modify preconditioner function - destroy - the function to use to destroy the user provided application context. ## Calling Sequence of `function` ```none PetscErrorCode function(KSP ksp, PetscInt n, PetscReal rnorm, void *ctx) ``` - ***ksp -*** iterative context - ***n -*** the total number of PIPEGCR iterations that have occurred - ***rnorm -*** 2-norm residual value - ***ctx -*** the user provided application context ## Calling Sequence of `destroy` ```none PetscErrorCode destroy(void *ctx) ``` ## Notes The default modifypc routine is `KSPPIPEGCRModifyPCNoChange()` ## See Also [](ch_ksp), `KSPPIPEGCR`, `KSPPIPEGCRModifyPCNoChange()` ## Level intermediate ## Location src/ksp/ksp/impls/gcr/pipegcr/pipegcr.c ## Implementations KSPPIPEGCRSetModifyPC_PIPEGCR in src/ksp/ksp/impls/gcr/pipegcr/pipegcr.c
--- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/ksp/ksp/impls/gcr/pipegcr/pipegcr.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)