KSPFGMRES#

Implements the Flexible Generalized Minimal Residual method. Flexible Krylov Methods

Options Database Keys#

  • -ksp_gmres_restart - the number of Krylov directions to orthogonalize against

  • -ksp_gmres_haptol - sets the tolerance for “happy ending” (exact convergence)

  • -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of vectors are allocated as needed)

  • -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)

  • -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)

  • -ksp_gmres_cgs_refinement_type <refine_never,refine_ifneeded,refine_always> - determine if iterative refinement is used to increase the stability of the classical Gram-Schmidt orthogonalization.

  • -ksp_gmres_krylov_monitor - plot the Krylov space generated

  • -ksp_fgmres_modifypcnochange - do not change the preconditioner between iterations

  • -ksp_fgmres_modifypcksp - modify the preconditioner using KSPFGMRESModifyPCKSP()

Notes#

See KSPFGMRESSetModifyPC() for how to vary the preconditioner between iterations

Only right preconditioning is supported.

The following options -ksp_type fgmres -pc_type ksp -ksp_ksp_type bcgs -ksp_view -ksp_pc_type jacobi make the preconditioner (or inner solver) be bi-CG-stab with a preconditioner of PCJACOBI

Developer Note#

This object is subclassed off of KSPGMRES, see the source code in src/ksp/ksp/impls/gmres for comments on the structure of the code

Contributed by#

Allison Baker

See Also#

KSP: Linear System Solvers, Flexible Krylov Methods, KSPCreate(), KSPSetType(), KSPType, KSP, KSPGMRES, KSPLGMRES, KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(), KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESMonitorKrylov(), KSPFGMRESSetModifyPC(), KSPFGMRESModifyPCKSP()

Level#

beginner

Location#

src/ksp/ksp/impls/gmres/fgmres/fgmres.c

Examples#

src/ts/tutorials/ex30.c
src/dm/impls/stag/tutorials/ex8.c
src/dm/impls/stag/tutorials/ex3.c
src/dm/impls/stag/tutorials/ex2.c
src/ksp/pc/tutorials/ex4.c
src/dm/impls/stag/tutorials/ex4.c
src/ksp/ksp/tutorials/ex7.c


Index of all KSP routines
Table of Contents for all manual pages
Index of all manual pages