KSPFGMRESSetModifyPC#
Sets the routine used by KSPFGMRES to modify the preconditioner. Flexible Krylov Methods
Synopsis#
#include "petscksp.h" 
PetscErrorCode KSPFGMRESSetModifyPC(KSP ksp, PetscErrorCode (*fcn)(KSP ksp, PetscInt total_its, PetscInt local_its, PetscReal res_norm, void *ctx), void *ctx, PetscErrorCode (*destroy)(void *ctx))
Logically Collective
Input Parameters#
- ksp - iterative context obtained from - KSPCreate()
- fcn - modifypc function 
- ctx - optional context 
- destroy - optional context destroy routine 
Calling sequence of fcn#
- ksp - the ksp context being used. 
- total_its - the total number of FGMRES iterations that have occurred. 
- local_its - the number of FGMRES iterations since last restart. 
- res_norm - the current residual norm. 
- ctx - optional context variable 
Calling sequence of destroy#
- ctx - the function context 
Options Database Keys#
- -ksp_fgmres_modifypcnochange - do not change the - PC
- -ksp_fgmres_modifypcksp - changes the inner KSP solver tolerances 
Note#
Several fcn routines are predefined, including KSPFGMRESModifyPCNoChange(), and  KSPFGMRESModifyPCKSP()
See Also#
KSP: Linear System Solvers, Flexible Krylov Methods, KSPFGMRES, KSPFGMRESModifyPCNoChange(), KSPFGMRESModifyPCKSP()
Level#
intermediate
Location#
Implementations#
KSPFGMRESSetModifyPC_FGMRES() in src/ksp/ksp/impls/gmres/fgmres/fgmres.c
Index of all KSP routines
Table of Contents for all manual pages
Index of all manual pages