PCHPDDM#

Interface with the HPDDM library. This PC may be used to build multilevel spectral domain decomposition methods based on the GenEO framework [2011, 2019]. It may be viewed as an alternative to spectral AMGe or PCBDDC with adaptive selection of constraints. The interface is explained in details in [2021] (see references below)

The matrix to be preconditioned (Pmat) may be unassembled (MATIS), assembled (MATAIJ, MATBAIJ, or MATSBAIJ), hierarchical (MATHTOOL), or MATNORMAL.

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#

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_chop

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 or -pc_hpddm_levels_N_eps_threshold . The default -pc_hpddm_coarse_p value is 1, meaning that the coarse operator is aggregated on a single process.

This preconditioner requires that you build PETSc 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. [2011, 2013]. 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#

  • 2011 - A robust two-level domain decomposition preconditioner for systems of PDEs. Spillane, Dolean, Hauret, Nataf, Pechstein, and Scheichl. Comptes Rendus Mathematique.

  • 2013 - Scalable domain decomposition preconditioners for heterogeneous elliptic problems. Jolivet, Hecht, Nataf, and Prud’homme. SC13.

  • 2015 - An introduction to domain decomposition methods: algorithms, theory, and parallel implementation. Dolean, Jolivet, and Nataf. SIAM.

  • 2019 - A multilevel Schwarz preconditioner based on a hierarchy of robust coarse spaces. Al Daas, Grigori, Jolivet, and Tournier. SIAM Journal on Scientific Computing.

  • 2021 - KSPHPDDM and PCHPDDM: extending PETSc with advanced Krylov methods and robust multilevel overlapping Schwarz preconditioners. Jolivet, Roman, and Zampini. Computer & Mathematics with Applications.

  • 2022a - A robust algebraic domain decomposition preconditioner for sparse normal equations. Al Daas, Jolivet, and Scott. SIAM Journal on Scientific Computing.

  • 2022b - A robust algebraic multilevel domain decomposition preconditioner for sparse symmetric positive definite matrices. Al Daas and Jolivet.

See Also#

PCCreate(), PCSetType(), PCType, PC, PCHPDDMSetAuxiliaryMat(), MATIS, PCBDDC, PCDEFLATION, PCTELESCOPE, PCASM, PCHPDDMSetCoarseCorrectionType(), PCHPDDMHasNeumannMat(), PCHPDDMSetRHSMat(), PCHPDDMSetDeflationMat(), PCHPDDMSetSTShareSubKSP(), PCHPDDMGetSTShareSubKSP(), PCHPDDMGetCoarseCorrectionType(), PCHPDDMGetComplexities()

Level#

intermediate

Location#

src/ksp/pc/impls/hpddm/pchpddm.cxx

Examples#

src/ksp/ksp/tutorials/ex76.c
src/ksp/ksp/tutorials/ex76f.F90
src/ksp/ksp/tutorials/ex82.c


Edit on GitLab

Index of all PC routines
Table of Contents for all manual pages
Index of all manual pages