PCDEFLATION#

Deflation preconditioner shifts (deflates) part of the spectrum to zero or to a predefined value.

Options Database Keys#

  • -pc_deflation_init_only - if true computes only the special guess

  • -pc_deflation_max_lvl <0> - maximum number of levels for multilevel deflation

  • -pc_deflation_reduction_factor <-1> - reduction factor on bottom level coarse problem for PCTELESCOPE (default based on the size of the coarse problem)

  • -pc_deflation_correction - if true apply coarse problem correction

  • -pc_deflation_correction_factor <1.0> - sets coarse problem correction factor

  • -pc_deflation_compute_space - compute PCDeflationSpaceType deflation space

  • -pc_deflation_compute_space_size <1> - size of the deflation space (corresponds to number of levels for wavelet-based deflation)

Notes#

Given a (complex - transpose is always Hermitian) full rank deflation matrix W, the deflation (introduced in [1,2]) preconditioner uses projections Q = W*(W’AW)^{-1}W’ and P = I - QA, where A is pmat.

The deflation computes initial guess x0 = x_{-1} - Qr_{-1}, which is the solution on the deflation space. If PCDeflationSetInitOnly() or -pc_deflation_init_only is set to PETSC_TRUE (InitDef scheme), the application of the preconditioner consists only of application of the additional preconditioner M^{-1}. Otherwise, the preconditioner application consists of PM^{-1} + factorQ. The first part of the preconditioner (PM^{-1}) shifts some eigenvalues to zero while the addition of the coarse problem correction (factorQ) makes the preconditioner to shift some eigenvalues to the given factor. The InitDef scheme is recommended for deflation using high accuracy estimates of eigenvectors of A when it exhibits similar convergence to the full deflation but is cheaper.

The deflation matrix is by default automatically computed. The type of deflation matrix and its size to compute can be controlled by PCDeflationSetSpaceToCompute() or -pc_deflation_compute_space and -pc_deflation_compute_space_size. User can set an arbitrary deflation space matrix with PCDeflationSetSpace(). If the deflation matrix is a multiplicative MATCOMPOSITE, a multilevel deflation [3] is used. The first matrix in the composite is used as the deflation matrix, and the coarse problem (W’AW)^{-1} is solved by KSPFCG (if A is MAT_SPD) or KSPFGMRES preconditioned by deflation with deflation matrix being the next matrix in the MATCOMPOSITE. This scheme repeats until the maximum level is reached or there are no more matrices. If the maximum level is reached, the remaining matrices are merged (multiplied) to create the last deflation matrix. The maximum level defaults to 0 and can be set by PCDeflationSetLevels() or by -pc_deflation_levels.

The coarse problem KSP can be controlled from the command line with prefix -deflation_ for the first level and -deflation_[lvl-1] from the second level onward. You can also use PCDeflationGetCoarseKSP() to control it from code. The bottom level KSP defaults to KSPPREONLY with PCLU direct solver (MATSOLVERSUPERLU/MATSOLVERSUPERLU_DIST if available) wrapped into PCTELESCOPE. For convenience, the reduction factor can be set by PCDeflationSetReductionFactor() or -pc_deflation_recduction_factor. The default is chosen heuristically based on the coarse problem size.

The additional preconditioner can be controlled from command line with prefix -deflation_[lvl]pc (same rules used for coarse problem KSP apply for [lvl] part of prefix), e.g., -deflation_1_pc_pc_type bjacobi. You can also use PCDeflationGetPC() to control the additional preconditioner from code. It defaults to PCNONE.

The coarse problem correction term (factor*Q) can be turned on by -pc_deflation_correction and the factor value can be set by pc_deflation_correction_factor or by PCDeflationSetCorrectionFactor(). The coarse problem can significantly improve convergence when the deflation coarse problem is not solved with high enough accuracy. We recommend setting factor to some eigenvalue, e.g., the largest eigenvalue so that the preconditioner does not create an isolated eigenvalue.

The options are automatically inherited from the previous deflation level.

The preconditioner supports KSPMonitorDynamicTolerance(). This is useful for the multilevel scheme for which we also recommend limiting the number of iterations for the coarse problems.

See section 3 of [4] for additional references and description of the algorithm when used for conjugate gradients. Section 4 describes some possible choices for the deflation space.

Contributed by Jakub Kruzik (PERMON), Institute of Geonics of the Czech Academy of Sciences and VSB - TU Ostrava.

Developed from PERMON code used in [4] while on a research stay with Prof. Reinhard Nabben at the Institute of Mathematics, TU Berlin.

References#

  • **** -*** A. Nicolaides. “Deflation of conjugate gradients with applications to boundary value problems”, SIAM J. Numer. Anal. 24.2, 1987.

  • **** -*** Z. Dostal. “Conjugate gradient method with preconditioning by projector”, Int J. Comput. Math. 23.3-4, 1988.

  • **** -*** Y. A. Erlangga and R. Nabben. “Multilevel Projection-Based Nested Krylov Iteration for Boundary Value Problems”, SIAM J. Sci. Comput. 30.3, 2008.

  • **** -*** J. Kruzik “Implementation of the Deflated Variants of the Conjugate Gradient Method”, Master’s thesis, VSB-TUO, 2018 - http://dspace5.vsb.cz/bitstream/handle/10084/130303/KRU0097_USP_N2658_2612T078_2018.pdf

See Also#

PCCreate(), PCSetType(), PCType, PC, PCDeflationSetInitOnly(), PCDeflationSetLevels(), PCDeflationSetReductionFactor(), PCDeflationSetCorrectionFactor(), PCDeflationSetSpaceToCompute(), PCDeflationSetSpace(), PCDeflationSpaceType, PCDeflationSetProjectionNullSpaceMat(), PCDeflationSetCoarseMat(), PCDeflationGetCoarseKSP(), PCDeflationGetPC()

Level#

intermediate

Location#

src/ksp/pc/impls/deflation/deflation.c


Edit on GitLab

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