petsc-3.6.1 2015-08-06
Report Typos and Errors

KSPGMRESSetOrthogonalization

Sets the orthogonalization routine used by GMRES and FGMRES.

Synopsis

#include "petscksp.h"  
PetscErrorCode  KSPGMRESSetOrthogonalization(KSP ksp,PetscErrorCode (*fcn)(KSP,PetscInt))
Logically Collective on KSP

Input Parameters

ksp - iterative context obtained from KSPCreate
fcn - orthogonalization function

Calling Sequence of function

  errorcode = int fcn(KSP ksp,int 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()

Keywords

KSP, GMRES, set, orthogonalization, Gram-Schmidt, iterative refinement

See Also

KSPGMRESSetRestart(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetCGSRefinementType(), KSPGMRESSetOrthogonalization(),
KSPGMRESModifiedGramSchmidtOrthogonalization(), KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESGetCGSRefinementType()

Level:intermediate
Location:
src/ksp/ksp/impls/gmres/gmres2.c
Index of all KSP routines
Table of Contents for all manual pages
Index of all manual pages