KSPGMRESGetOrthogonalization#

Gets the orthogonalization routine used by KSPGMRES and KSPFGMRES.

Synopsis#

#include "petscksp.h"  
PetscErrorCode KSPGMRESGetOrthogonalization(KSP ksp, PetscErrorCode (**fcn)(KSP ksp, PetscInt it))

Not Collective

Input Parameter#

  • ksp - iterative context obtained from KSPCreate()

Output Parameter#

  • fcn - orthogonalization function

Calling sequence of fcn#

  • ksp - the solver context

  • it - the current iteration

Notes#

Two orthogonalization routines are predefined, 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

Implementations#

KSPGMRESGetOrthogonalization_GMRES() in src/ksp/ksp/impls/gmres/gmres.c


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