petsc-3.7.7 2017-09-25
Report Typos and Errors

KSPPIPEGCR

Implements the preconditioned Generalized Conjugate Residual method with pipelining. The PIPEGCR Krylov method supports non-symmetric matrices and permits the use of a preconditioner which may vary from one iteration to the next. Users can can define a method to vary the preconditioner between iterates via KSPPIPEGCRSetModifyPC(). Restarts are solves with x0 not equal to zero. When a restart occurs, the initial starting solution is given by the current estimate for x which was obtained by the last restart iterations of the PIPEGCR algorithm. The method implemented requires at most the storage of 4 x mmax + 5 vectors, roughly twice as much as GCR.

Only supports left preconditioning.

The natural norm for this method is (u,Au). This norm is available at no computational costs. Choosing norm types preconditioned or unpreconditioned involves a blocking reduction which prevents any benefit from pipelining.

Options Database Keys

-ksp_pipegcr_mmax <N> - the max number of Krylov directions to orthogonalize against
-ksp_pipegcr_unroll_w - unroll w at the storage cost of a maximum of (mmax+1) extra vectors with the benefit of better pipelining (default: PETSC_TRUE)
-ksp_pipegcr_nprealloc <N> - the number of vectors to preallocated for storing Krylov directions. Once exhausted new directions are allocated blockwise (default: 5)
-ksp_pipegcr_truncation - Truncate number of previous Krylov directions
-ksp_pipegcr_trancation_restart - Truncation-restart strategy: Keep at most mmax Krylov directions then restart (the default)

Reference

Pipelined, Flexible Krylov Subspace Methods Patrick Sanan, Sascha M. Schnepp, and Dave A. May

See Also

KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP,
KSPPIPEFGMRES, KSPPIPECG, KSPPIPECR, KSPPIPEFCG

Level:beginner
Location:
src/ksp/ksp/impls/gcr/pipegcr/pipegcr.c
Index of all KSP routines
Table of Contents for all manual pages
Index of all manual pages