petsc-3.7.1 2016-05-15
Report Typos and Errors

KSPFGMRES

Implements the Flexible Generalized Minimal Residual method. developed by Saad with restart

Options Database Keys

-ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against Many br
-ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence) Many br
-ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of Many brvectors are allocated as needed) Many br
-ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default) Many br
-ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower) Many br
-ksp_gmres_cgs_refinement_type <never,ifneeded,always> - determine if iterative refinement is used to increase the Many brstability of the classical Gram-Schmidt orthogonalization. Many br
-ksp_gmres_krylov_monitor - plot the Krylov space generated Many br
-ksp_fgmres_modifypcnochange - do not change the preconditioner between iterations Many br
-ksp_fgmres_modifypcksp - modify the preconditioner using KSPFGMRESModifyPCKSP() Many br

Many br

Notes: See KSPFGMRESSetModifyPC() for how to vary the preconditioner between iterations Many brOnly right preconditioning is supported. Many br

Notes: 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) Many brbe bi-CG-stab with a preconditioner of Jacobi. Many br

Developer Notes: This object is subclassed off of KSPGMRES Many br

See Also

KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPGMRES, KSPLGMRES,
KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), Many brKSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(), Many brKSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESMonitorKrylov(), KSPFGMRESSetModifyPC(), Many brKSPFGMRESModifyPCKSP() Many br

Level:beginner
Location:
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