KSPGCR#

Implements the preconditioned flexible Generalized Conjugate Residual (GCR) method [EES83]. Flexible Krylov Methods

Options Database Key#

  • -ksp_gcr_restart - the number of stored vectors to orthogonalize against

Notes#

This method supports non-symmetric matrices and permits the use of a preconditioner which may vary from one iteration to the next.

Users can define a method to vary the preconditioner between iterates via KSPFlexibleSetModifyPC().

When a restart occurs, the initial starting solution is given by the current estimate for x.

Unlike KSPGMRES and KSPFGMRES, when using GCR, the solution and residual vector can be directly accessed at any iterate, with zero computational cost, via a call to KSPBuildSolution() and KSPBuildResidual() respectively.

The stopping condition test is only applied after the iteration count exceeds the value set by KSPSetCheckNormIteration() (also available as -ksp_check_norm_iteration). The residual norm reported by the monitor and stored in the residual history will be listed as 0.0 before that iteration; the norm is not actually zero, it is simply not computed until then.

The method implemented requires the storage of 2 x restart + 1 vectors, twice as much as KSPGMRES.

Supports only right preconditioning.

Contributed by#

Dave May

References#

[EES83]

S.C. Eisenstat, H.C. Elman, and M.H. Schultz. Variational iterative methods for nonsymmetric systems of linear equations. SIAM Journal on Numerical Analysis, 20(2):345–357, 1983.

See Also#

KSP: Linear System Solvers, Flexible Krylov Methods, KSPFCG, KSPPIPEGCR, KSPPIPEFCG, KSPFGMRES, KSPCG, KSPCreate(), KSPSetType(), KSPType, KSP, KSPGCRSetRestart(), KSPGCRGetRestart(), KSPFlexibleSetModifyPC(), KSPGMRES

Level#

beginner

Location#

src/ksp/ksp/impls/gcr/gcr.c

Examples#

src/dm/impls/stag/tutorials/ex4.c


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