KSPGMRESSetOrthogonalization#

Sets the orthogonalization routine used by KSPGMRES and KSPFGMRES.

Synopsis#

Logically Collective

Input Parameters#

  • ksp - iterative context obtained from KSPCreate()

  • fcn - orthogonalization function

Calling Sequence of function#

errorcode = PetscErrorCode fcn(KSP ksp,PetscInt it);
it is one minus the number of GMRES iterations since last restart;
i.e. the size of Krylov space minus one

Options Database Keys#

  • -ksp_gmres_classicalgramschmidt - Activates KSPGMRESClassicalGramSchmidtOrthogonalization() (default)

  • -ksp_gmres_modifiedgramschmidt - Activates KSPGMRESModifiedGramSchmidtOrthogonalization()

Notes#

Two orthogonalization routines are predefined, including KSPGMRESModifiedGramSchmidtOrthogonalization() and the default KSPGMRESClassicalGramSchmidtOrthogonalization().

Use KSPGMRESSetCGSRefinementType() to determine if iterative refinement is used to increase stability.

See Also#

KSP: Linear System Solvers, KSPGMRESSetRestart(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetCGSRefinementType(), KSPGMRESSetOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(), KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESGetCGSRefinementType()

Level#

intermediate

Location#

src/ksp/ksp/impls/gmres/gmres2.c

Examples#

src/ksp/ksp/tutorials/ex5f.F90.html

Implementations#

KSPGMRESSetOrthogonalization_GMRES in src/ksp/ksp/impls/gmres/gmres.c


Edit on GitLab

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