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, 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
#
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
#
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#
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