KSPGMRESSetOrthogonalization#
Sets the orthogonalization routine used by KSPGMRES
and KSPFGMRES
.
Synopsis#
#include "petscksp.h"
PetscErrorCode KSPGMRESSetOrthogonalization(KSP ksp, PetscErrorCode (*fcn)(KSP ksp, PetscInt it))
Logically Collective
Input Parameters#
ksp - iterative context obtained from
KSPCreate()
fcn - orthogonalization function
Calling sequence of fcn
#
ksp - the solver context
it - the current iteration
Options Database Keys#
-ksp_gmres_classicalgramschmidt - Activates KSPGMRESClassicalGramSchmidtOrthogonalization() (default)
-ksp_gmres_modifiedgramschmidt - Activates KSPGMRESModifiedGramSchmidtOrthogonalization()
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()
, KSPGMRESModifiedGramSchmidtOrthogonalization()
,
KSPGMRESClassicalGramSchmidtOrthogonalization()
, KSPGMRESGetCGSRefinementType()
Level#
intermediate
Location#
Implementations#
KSPGMRESSetOrthogonalization_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