PCHPDDM#
Interface with the HPDDM library. This PC
may be used to build multilevel spectral domain decomposition methods based on the GenEO framework [SDH+11] [ADGJT21].
It may be viewed as an alternative to spectral
AMGe or PCBDDC
with adaptive selection of constraints. The interface is explained in details in [JRZ21]
The matrix used for building the preconditioner (Pmat) may be unassembled (MATIS
), assembled (MATAIJ
, MATBAIJ
, or MATSBAIJ
), hierarchical (MATHTOOL
), MATNORMAL
, MATNORMALHERMITIAN
, or MATSCHURCOMPLEMENT
(when PCHPDDM
is used as part of an outer PCFIELDSPLIT
).
For multilevel preconditioning, when using an assembled or hierarchical Pmat, one must provide an auxiliary local Mat
(unassembled local operator for GenEO) using
PCHPDDMSetAuxiliaryMat()
. Calling this routine is not needed when using a MATIS
Pmat, assembly is done internally using MatConvert()
.
Options Database Keys#
-pc_hpddm_define_subdomains <true, default=false> - on the finest level, calls
PCASMSetLocalSubdomains()
with theIS
supplied inPCHPDDMSetAuxiliaryMat()
(not relevant with an unassembled Pmat)-pc_hpddm_has_neumann <true, default=false> - on the finest level, informs the
PC
that the local Neumann matrix is supplied inPCHPDDMSetAuxiliaryMat()
-pc_hpddm_coarse_correction <type, default=deflated> - determines the
PCHPDDMCoarseCorrectionType
when callingPCApply()
Options for subdomain solvers, subdomain eigensolvers (for computing deflation vectors), and the coarse solver can be set using the following options database prefixes.
-pc_hpddm_levels_%d_pc_
-pc_hpddm_levels_%d_ksp_
-pc_hpddm_levels_%d_eps_
-pc_hpddm_levels_%d_p
-pc_hpddm_levels_%d_mat_type
-pc_hpddm_coarse_
-pc_hpddm_coarse_p
-pc_hpddm_coarse_mat_type
-pc_hpddm_coarse_mat_filter
E.g., -pc_hpddm_levels_1_sub_pc_type lu -pc_hpddm_levels_1_eps_nev 10 -pc_hpddm_levels_2_p 4 -pc_hpddm_levels_2_sub_pc_type lu -pc_hpddm_levels_2_eps_nev 10
-pc_hpddm_coarse_p 2 -pc_hpddm_coarse_mat_type baij will use 10 deflation vectors per subdomain on the fine “level 1”,
aggregate the fine subdomains into 4 “level 2” subdomains, then use 10 deflation vectors per subdomain on “level 2”,
and assemble the coarse matrix (of dimension 4 x 10 = 40) on two processes as a MATBAIJ
(default is MATSBAIJ
).
In order to activate a “level N+1” coarse correction, it is mandatory to call -pc_hpddm_levels_N_eps_nev
Notes#
This preconditioner requires that PETSc is built with SLEPc (--download-slepc
).
By default, the underlying concurrent eigenproblems
are solved using SLEPc shift-and-invert spectral transformation. This is usually what gives the best performance for GenEO, cf.
[SDH+11] [JHNPrudhomme13]. As
stated above, SLEPc options are available through -pc_hpddm_levels_%d_, e.g., -pc_hpddm_levels_1_eps_type arpack -pc_hpddm_levels_1_eps_nev 10
-pc_hpddm_levels_1_st_type sinvert. There are furthermore three options related to the (subdomain-wise local) eigensolver that are not described in
SLEPc documentation since they are specific to PCHPDDM
.
-pc_hpddm_levels_1_st_share_sub_ksp
-pc_hpddm_levels_%d_eps_threshold
-pc_hpddm_levels_1_eps_use_inertia
The first option from the list only applies to the fine-level eigensolver, see PCHPDDMSetSTShareSubKSP()
. The second option from the list is
used to filter eigenmodes retrieved after convergence of EPSSolve()
at “level N” such that eigenvectors used to define a “level N+1” coarse
correction are associated to eigenvalues whose magnitude are lower or equal than -pc_hpddm_levels_N_eps_threshold. When using an EPS
which cannot
determine a priori the proper -pc_hpddm_levels_N_eps_nev such that all wanted eigenmodes are retrieved, it is possible to get an estimation of the
correct value using the third option from the list, -pc_hpddm_levels_1_eps_use_inertia, see MatGetInertia()
. In that case, there is no need
to supply -pc_hpddm_levels_1_eps_nev. This last option also only applies to the fine-level (N = 1) eigensolver.
References#
Hussam Al Daas, Laura Grigori, Pierre Jolivet, and Pierre-Henri Tournier. A multilevel Schwarz preconditioner based on a hierarchy of robust coarse spaces. SIAM Journal on Scientific Computing, 43(3):A1907–A1928, 2021.
Hussam Al Daas and Pierre Jolivet. A robust algebraic multilevel domain decomposition preconditioner for sparse symmetric positive definite matrices. SIAM Journal on Scientific Computing, 44(4):A2582–A2598, 2022.
Hussam Al Daas, Pierre Jolivet, and Jennifer A Scott. A robust algebraic domain decomposition preconditioner for sparse normal equations. SIAM Journal on Scientific Computing, 44(3):A1047–A1068, 2022.
Victorita Dolean, Pierre Jolivet, and Frédéric Nataf. An introduction to domain decomposition methods: algorithms, theory, and parallel implementation. SIAM, 2015.
Pierre Jolivet, Frédéric Hecht, Frédéric Nataf, and Christophe Prud'homme. Scalable domain decomposition preconditioners for heterogeneous elliptic problems. In Proceedings of SC13: International Conference for High Performance Computing, Networking, Storage and Analysis, 80. ACM, 2013.
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.
Frédéric Nataf and Pierre-Henri Tournier. Recent advances in domain decomposition methods for large-scale saddle point problems. Comptes Rendus. Mécanique, 350(S1):1–15, 2022.
See Also#
KSP: Linear System Solvers, PCCreate()
, PCSetType()
, PCType
, PC
, PCHPDDMSetAuxiliaryMat()
, MATIS
, PCBDDC
, PCDEFLATION
, PCTELESCOPE
, PCASM
,
PCHPDDMSetCoarseCorrectionType()
, PCHPDDMHasNeumannMat()
, PCHPDDMSetRHSMat()
, PCHPDDMSetDeflationMat()
, PCHPDDMSetSTShareSubKSP()
,
PCHPDDMGetSTShareSubKSP()
, PCHPDDMGetCoarseCorrectionType()
, PCHPDDMGetComplexities()
Level#
intermediate
Location#
Examples#
src/ksp/ksp/tutorials/ex82.c
src/ksp/ksp/tutorials/ex76f.F90
src/ksp/ksp/tutorials/ex27.c
src/ksp/ksp/tutorials/ex76.c
src/snes/tutorials/ex12.c
src/ksp/ksp/tutorials/ex87.c
Index of all PC routines
Table of Contents for all manual pages
Index of all manual pages