KSPLGMRES#

Augments the standard GMRES approximation space with approximations to the error from previous restart cycles as in [BJM05].

Options Database Keys#

  • -ksp_gmres_restart - total approximation space size (Krylov directions + error approximations)

  • -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_lgmres_augment - number of error approximations to augment the Krylov space with

  • -ksp_lgmres_constant - use a constant approximate space size (only affects restart cycles < num. error approx.(k), i.e. the first k restarts)

Notes#

Supports both left and right preconditioning, but not symmetric.

Augmenting with 1,2, or 3 approximations is generally optimal.

This method is an accelerator for KSPGMRES - it is not useful for problems that stall. When gmres(m) stalls then lgmres with a size m approximation space will also generally stall.

If gmres(m) converges in a small number of restart cycles, then lgmres also tends not to be very helpful.

Developer Notes#

To run LGMRES(m, k) as described in [BJM05], use:

   -ksp_gmres_restart <m+k>
  -ksp_lgmres_augment <k>

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

References#

BJM05(1,2)

AH Baker, ER Jessup, and T Manteuffel. A technique for accelerating the convergence of restarted gmres. SIAM JOURNAL ON MATRIX ANALYSIS AND APPLICATIONS, 2005.

See Also#

KSP: Linear System Solvers, KSPCreate(), KSPSetType(), KSPType, KSP, KSPFGMRES, KSPGMRES, KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(), KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESMonitorKrylov(), KSPLGMRESSetAugDim(), KSPGMRESSetConstant(), KSPLGMRESSetConstant(), KSPLGMRESSetAugDim()

Level#

beginner

Location#

src/ksp/ksp/impls/gmres/lgmres/lgmres.c


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