KSPGMRESClassicalGramSchmidtOrthogonalization#

This is the basic orthogonalization routine using classical Gram-Schmidt with possible iterative refinement to improve the stability

Synopsis#

Collective, No Fortran Support

Input Parameters#

  • ksp - KSP object, must be associated with KSPGMRES, KSPFGMRES, or KSPLGMRES Krylov method

  • it - one less than the current GMRES restart iteration, i.e. the size of the Krylov space

Options Database Keys#

  • -ksp_gmres_classicalgramschmidt - Activates KSPGMRESClassicalGramSchmidtOrthogonalization()

  • -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.

Note#

Use KSPGMRESSetCGSRefinementType() to determine if iterative refinement is to be used. This is much faster than KSPGMRESModifiedGramSchmidtOrthogonalization() but has the small possibility of stability issues that can usually be handled by using a single step of iterative refinement with KSPGMRESSetCGSRefinementType()

See Also#

KSP: Linear System Solvers, KSPGMRESCGSRefinementType, KSPGMRESSetOrthogonalization(), KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESGetOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization()

Level#

intermediate

Location#

src/ksp/ksp/impls/gmres/borthog2.c


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