KSPHPDDM#
Interface with the HPDDM library. This KSP may be used to further select methods that are currently not implemented natively in PETSc, e.g., GCRODR [PdSM+06], a recycled Krylov method which is similar to KSPLGMRES, see [JT16] for a comparison. ex75.c shows how to reproduce the results from the aforementioned paper [PdSM+06]. A chronological bibliography of relevant publications linked with KSP available in HPDDM through KSPHPDDM, and not available directly in PETSc, may be found below. The interface is explained in details in [JRZ21]. See also [OLeary80], [JL17] and [CGL+13]
Options Database Keys#
-ksp_gmres_restart restart - see
KSPGMRES, default is 30-ksp_hpddm_type type - see
KSPHPDDMType, default isgmres-ksp_hpddm_precision (__fp16|single|double|__float128) - see
PetscPrecision, default is the same asPetscScalar-ksp_hpddm_deflation_tol eps - tolerance when deflating right-hand sides inside block methods (only relevant with block methods), default is -1.0
-ksp_hpddm_enlarge_krylov_subspace p - split the initial right-hand side into multiple vectors (only relevant with nonblock methods), default is 1
-ksp_hpddm_orthogonalization (cgs|mgs) - see
KSPGMRES, default iscgs-ksp_hpddm_qr (cholqr|cgs|mgs) - distributed QR factorizations, only relevant with block methods, default is
cholqr-ksp_hpddm_variant (left|right|flexible) - this option is superseded by
KSPSetPCSide(), default isleft-ksp_hpddm_recycle n - number of harmonic Ritz vectors to compute (only relevant with GCRODR or BGCRODR), default is 0
-ksp_hpddm_recycle_target (SM|LM|SR|LR|SI|LI) - criterion to select harmonic Ritz vectors, default is
SM(only relevant with GCRODR or BGCRODR). For BGCRODR, if PETSc is compiled with SLEPc, this option is not relevant, since SLEPc is used instead. Options are set with the prefix-ksp_hpddm_recycle_eps_-ksp_hpddm_recycle_strategy (A|B) - generalized eigenvalue problem to solve for recycling (only relevant with flexible GCRODR or BGCRODR), default is
A-ksp_hpddm_recycle_symmetric (true|false) - symmetric generalized eigenproblems in BGCRODR, useful to switch to distributed solvers like
EPSELEMENTALorEPSSCALAPACK(only relevant when PETSc is compiled with SLEPc), default isfalse
References#
Henri Calandra, Serge Gratton, Rafael Lago, Xavier Vasseur, and Luiz Mariano Carvalho. A modified block flexible GMRES method with deflation at each iteration for the solution of non-Hermitian linear systems with multiple right-hand sides. SIAM Journal on Scientific Computing, 35(5):S345–S367, 2013.
Hao Ji and Yaohang Li. A breakdown-free block conjugate gradient method. BIT Numerical Mathematics, 57:379–403, 2017.
Pierre Jolivet, Jose Roman, and Stefano Zampini. KSPHPDDM and PCHPDDM: extending PETSc with advanced Krylov methods and robust multilevel overlapping Schwarz preconditioners. Computers and Mathematics with Applications, 84:277–295, 2021.
Pierre Jolivet and Pierre-Henri Tournier. Block iterative methods and recycling for improved scalability of linear solvers. In SC'16: Proceedings of the International Conference for High Performance Computing, Networking, Storage and Analysis, 190–203. IEEE, 2016.
Dianne P O'Leary. The block conjugate gradient algorithm and related methods. Linear algebra and its applications, 29:293–322, 1980.
See Also#
KSP: Linear System Solvers, Flexible Krylov Methods, KSPCreate(), KSPSetType(), KSPType, KSP, KSPGMRES, KSPCG, KSPLGMRES, KSPDGMRES
Level#
intermediate
Location#
Examples#
src/ksp/ksp/tutorials/ex75.c
src/ksp/ksp/tutorials/ex75f.F90
src/ksp/ksp/tutorials/ex79.c
src/ksp/ksp/tutorials/ex76.c
src/ksp/ksp/tutorials/ex76f.F90
src/ksp/ksp/tutorials/ex78.c
Index of all KSP routines
Table of Contents for all manual pages
Index of all manual pages