KSPFCG#
Implements the Flexible Conjugate Gradient method (FCG) [Not00], [AV91]. Unlike most KSP
methods this allows the preconditioner to be nonlinear. Flexible Krylov Methods
Options Database Keys#
-ksp_fcg_mmax
- maximum number of search directions-ksp_fcg_nprealloc
- number of directions to preallocate-ksp_fcg_truncation_type <standard,notay> - truncation approach for directions
Notes#
Compare to KSPFCG
Supports left preconditioning only.
Contributed by#
Patrick Sanan
References#
- AV91
O. Axelsson and P. S. Vassilevski. A black box generalized Conjugate Gradient solver with inner iterations and variable-step preconditioning. SIAM Journal on Matrix Analysis and Applications, 12(4):625–644, 1991.
- Not00
Yvan Notay. Flexible conjugate gradients. SIAM Journal on Scientific Computing, 22(4):1444–1460, 2000.
See Also#
KSP: Linear System Solvers, Flexible Krylov Methods, KSPGCR
, KSPFGMRES
, KSPCG
, KSPFCGSetMmax()
, KSPFCGGetMmax()
, KSPFCGSetNprealloc()
, KSPFCGGetNprealloc()
, KSPFCGSetTruncationType()
, KSPFCGGetTruncationType()
,
KSPFCGGetTruncationType
Level#
beginner
Location#
Index of all KSP routines
Table of Contents for all manual pages
Index of all manual pages