petsc-3.14.6 2021-03-30
Report Typos and Errors

KSPGMRESGetOrthogonalization

Gets the orthogonalization routine used by GMRES and FGMRES.

Synopsis

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

Input Parameter

ksp - iterative context obtained from KSPCreate

Output Parameter

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

Notes

Two orthogonalization routines are predefined, including

KSPGMRESModifiedGramSchmidtOrthogonalization()

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

Options Database Keys

-ksp_gmres_classicalgramschmidt - Activates KSPGMRESClassicalGramSchmidtOrthogonalization() (default)
-ksp_gmres_modifiedgramschmidt - Activates KSPGMRESModifiedGramSchmidtOrthogonalization()

See Also

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