KSPPIPEGCR#

Implements a Pipelined Generalized Conjugate Residual method [SSM16]. Flexible Krylov Methods. Pipelined Krylov Methods

Options Database Keys#

  • -ksp_pipegcr_mmax - 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 - the number of vectors to preallocated for storing Krylov directions. Once exhausted new directions are allocated blockwise (default: 5)

  • -ksp_pipegcr_truncation_type <standard,notay> - which previous search directions to orthogonalize against

Notes#

Compare to KSPGCR

The KSPPIPEGCR Krylov 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 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 KSPPIPEGCR algorithm. The method implemented requires at most the storage of 4 x mmax + 5 vectors, roughly twice as much as KSPGCR.

Only supports left preconditioning.

The natural “norm” for this method is \((u,Au)\), where \(u\) is the preconditioned residual. This norm is available at no additional computational cost, as with standard KSPCG. Choosing preconditioned or unpreconditioned norm types involves a blocking reduction which prevents any benefit from pipelining.

MPI configuration may be necessary for reductions to make asynchronous progress, which is important for performance of pipelined methods. See What steps are necessary to make the pipelined solvers execute efficiently?

Contributed by#

Sascha M. Schnepp and Patrick Sanan

References#

[SSM16]

P. Sanan, S. M. Schnepp, and D. A. May. Pipelined, flexible Krylov subspace methods. SIAM Journal on Scientific Computing, 38(5):C441–C470, 2016. doi:10.1137/15M1049130.

See Also#

KSP: Linear System Solvers, Flexible Krylov Methods, Pipelined Krylov Methods, What steps are necessary to make the pipelined solvers execute efficiently?, KSPCreate(), KSPSetType(), KSPType, KSP, KSPPIPEFGMRES, KSPPIPECG, KSPPIPECR, KSPPIPEFCG, KSPPIPEGCRSetTruncationType(), KSPPIPEGCRSetNprealloc(), KSPPIPEGCRSetUnrollW(), KSPPIPEGCRSetMmax(), KSPPIPEGCRGetTruncationType(), KSPPIPEGCRGetNprealloc(), KSPPIPEGCRGetMmax(), KSPGCR, `KSPPIPEGCRGetUnrollW()

Level#

intermediate

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