1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108
|
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML3.2 EN">
<HTML>
<HEAD> <link rel="canonical" href="http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCDEFLATION.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<TITLE>PCDEFLATION</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
<div id="version" align=right><b>petsc-3.14.5 2021-03-03</b></div>
<div id="bugreport" align=right><a href="mailto:petsc-maint@mcs.anl.gov?subject=Typo or Error in Documentation &body=Please describe the typo or error in the documentation: petsc-3.14.5 v3.14.5 docs/manualpages/PC/PCDEFLATION.html "><small>Report Typos and Errors</small></a></div>
<A NAME="PCDEFLATION"><H1>PCDEFLATION</H1></A>
Deflation preconditioner shifts (deflates) part of the spectrum to zero or to a predefined value.
<H3><FONT COLOR="#CC3333">Options Database Keys</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-pc_deflation_init_only <false> </B></TD><TD>- if true computes only the special guess
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-pc_deflation_max_lvl <0> </B></TD><TD>- maximum number of levels for multilevel deflation
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-pc_deflation_reduction_factor <-1> </B></TD><TD>- reduction factor on bottom level coarse problem for <A HREF="../PC/PCTELESCOPE.html#PCTELESCOPE">PCTELESCOPE</A> (default based on the size of the coarse problem)
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-pc_deflation_correction <false> </B></TD><TD>- if true apply coarse problem correction
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-pc_deflation_correction_factor <1.0> </B></TD><TD>- sets coarse problem correction factor
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-pc_deflation_compute_space <haar> </B></TD><TD>- compute <A HREF="../PC/PCDeflationSpaceType.html#PCDeflationSpaceType">PCDeflationSpaceType</A> deflation space
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-pc_deflation_compute_space_size <1> </B></TD><TD>- size of the deflation space (corresponds to number of levels for wavelet-based deflation)
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#CC3333">Notes</FONT></H3>
Given a (complex - transpose is always Hermitian) full rank deflation matrix W, the deflation (introduced in [1,2])
preconditioner uses projections Q = W*(W'*A*W)^{-1}*W' and P = I - Q*A, where A is pmat.
<P>
The deflation computes initial guess x0 = x_{-1} - Q*r_{-1}, which is the solution on the deflation space.
If <A HREF="../PC/PCDeflationSetInitOnly.html#PCDeflationSetInitOnly">PCDeflationSetInitOnly</A>() or -pc_deflation_init_only is set to <A HREF="../Sys/PETSC_TRUE.html#PETSC_TRUE">PETSC_TRUE</A> (InitDef scheme), the application of the
preconditioner consists only of application of the additional preconditioner M^{-1}. Otherwise, the preconditioner
application consists of P*M^{-1} + factor*Q. The first part of the preconditioner (PM^{-1}) shifts some eigenvalues
to zero while the addition of the coarse problem correction (factor*Q) 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.
<P>
The deflation matrix is by default automatically computed. The type of deflation matrix and its size to compute can
be controlled by <A HREF="../PC/PCDeflationSetSpaceToCompute.html#PCDeflationSetSpaceToCompute">PCDeflationSetSpaceToCompute</A>() or -pc_deflation_compute_space and -pc_deflation_compute_space_size.
User can set an arbitrary deflation space matrix with <A HREF="../PC/PCDeflationSetSpace.html#PCDeflationSetSpace">PCDeflationSetSpace</A>(). If the deflation matrix
is a multiplicative <A HREF="../Mat/MATCOMPOSITE.html#MATCOMPOSITE">MATCOMPOSITE</A>, a multilevel deflation [3] is used. The first matrix in the composite is used as the
deflation matrix, and the coarse problem (W'*A*W)^{-1} is solved by <A HREF="../KSP/KSPFCG.html#KSPFCG">KSPFCG</A> (if A is <A HREF="../Mat/MatOption.html#MatOption">MAT_SPD</A>) or <A HREF="../KSP/KSPFGMRES.html#KSPFGMRES">KSPFGMRES</A> preconditioned
by deflation with deflation matrix being the next matrix in the <A HREF="../Mat/MATCOMPOSITE.html#MATCOMPOSITE">MATCOMPOSITE</A>. 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
<A HREF="../PC/PCDeflationSetLevels.html#PCDeflationSetLevels">PCDeflationSetLevels</A>() or by -pc_deflation_levels.
<P>
The coarse problem <A HREF="../KSP/KSP.html#KSP">KSP</A> 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
<A HREF="../PC/PCDeflationGetCoarseKSP.html#PCDeflationGetCoarseKSP">PCDeflationGetCoarseKSP</A>() to control it from code. The bottom level <A HREF="../KSP/KSP.html#KSP">KSP</A> defaults to
<A HREF="../KSP/KSPPREONLY.html#KSPPREONLY">KSPPREONLY</A> with <A HREF="../PC/PCLU.html#PCLU">PCLU</A> direct solver (<A HREF="../Mat/MATSOLVERSUPERLU.html#MATSOLVERSUPERLU">MATSOLVERSUPERLU</A>/<A HREF="../Mat/MATSOLVERSUPERLU_DIST.html#MATSOLVERSUPERLU_DIST">MATSOLVERSUPERLU_DIST</A> if available) wrapped into <A HREF="../PC/PCTELESCOPE.html#PCTELESCOPE">PCTELESCOPE</A>.
For convenience, the reduction factor can be set by <A HREF="../PC/PCDeflationSetReductionFactor.html#PCDeflationSetReductionFactor">PCDeflationSetReductionFactor</A>()
or -pc_deflation_recduction_factor. The default is chosen heuristically based on the coarse problem size.
<P>
The additional preconditioner can be controlled from command line with prefix -deflation_[lvl]_pc (same rules used for
coarse problem <A HREF="../KSP/KSP.html#KSP">KSP</A> apply for [lvl]_ part of prefix), e.g., -deflation_1_pc_pc_type bjacobi. You can also use
<A HREF="../PC/PCDeflationGetPC.html#PCDeflationGetPC">PCDeflationGetPC</A>() to control the additional preconditioner from code. It defaults to <A HREF="../PC/PCNONE.html#PCNONE">PCNONE</A>.
<P>
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 <A HREF="../PC/PCDeflationSetCorrectionFactor.html#PCDeflationSetCorrectionFactor">PCDeflationSetCorrectionFactor</A>(). 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.
<P>
The options are automatically inherited from the previous deflation level.
<P>
The preconditioner supports <A HREF="../KSP/KSPMonitorDynamicTolerance.html#KSPMonitorDynamicTolerance">KSPMonitorDynamicTolerance</A>(). This is useful for the multilevel scheme for which we also
recommend limiting the number of iterations for the coarse problems.
<P>
See section 3 of [4] for additional references and decription of the algorithm when used for conjugate gradients.
Section 4 describes some possible choices for the deflation space.
<P>
<H3><FONT COLOR="#CC3333">Developer Notes</FONT></H3>
Contributed by Jakub Kruzik (PERMON), Institute of Geonics of the Czech
Academy of Sciences and VSB - TU Ostrava.
<P>
Developed from PERMON code used in [4] while on a research stay with
Prof. Reinhard Nabben at the Institute of Mathematics, TU Berlin.
<P>
<H3><FONT COLOR="#CC3333">References</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>[1] </B></TD><TD>- A. Nicolaides. "Deflation of conjugate gradients with applications to boundary value problems", SIAM J. Numer. Anal. 24.2, 1987.
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>[2] </B></TD><TD>- Z. Dostal. "Conjugate gradient method with preconditioning by projector", Int J. Comput. Math. 23.3-4, 1988.
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>[3] </B></TD><TD>- Y. A. Erlangga and R. Nabben. "Multilevel Projection-Based Nested Krylov Iteration for Boundary Value Problems", SIAM J. Sci. Comput. 30.3, 2008.
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>[4] </B></TD><TD>- 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
</TD></TR></TABLE>
<P>
<P>
<H3><FONT COLOR="#CC3333">See Also</FONT></H3>
<A HREF="../PC/PCCreate.html#PCCreate">PCCreate</A>(), <A HREF="../PC/PCSetType.html#PCSetType">PCSetType</A>(), <A HREF="../PC/PCType.html#PCType">PCType</A> (for list of available types), <A HREF="../PC/PC.html#PC">PC</A>,
<BR><A HREF="../PC/PCDeflationSetInitOnly.html#PCDeflationSetInitOnly">PCDeflationSetInitOnly</A>(), <A HREF="../PC/PCDeflationSetLevels.html#PCDeflationSetLevels">PCDeflationSetLevels</A>(), <A HREF="../PC/PCDeflationSetReductionFactor.html#PCDeflationSetReductionFactor">PCDeflationSetReductionFactor</A>(),
<A HREF="../PC/PCDeflationSetCorrectionFactor.html#PCDeflationSetCorrectionFactor">PCDeflationSetCorrectionFactor</A>(), <A HREF="../PC/PCDeflationSetSpaceToCompute.html#PCDeflationSetSpaceToCompute">PCDeflationSetSpaceToCompute</A>(),
<A HREF="../PC/PCDeflationSetSpace.html#PCDeflationSetSpace">PCDeflationSetSpace</A>(), <A HREF="../PC/PCDeflationSpaceType.html#PCDeflationSpaceType">PCDeflationSpaceType</A>, <A HREF="../PC/PCDeflationSetProjectionNullSpaceMat.html#PCDeflationSetProjectionNullSpaceMat">PCDeflationSetProjectionNullSpaceMat</A>(),
<A HREF="../PC/PCDeflationSetCoarseMat.html#PCDeflationSetCoarseMat">PCDeflationSetCoarseMat</A>(), <A HREF="../PC/PCDeflationGetCoarseKSP.html#PCDeflationGetCoarseKSP">PCDeflationGetCoarseKSP</A>(), <A HREF="../PC/PCDeflationGetPC.html#PCDeflationGetPC">PCDeflationGetPC</A>()
<P><B></B><H3><FONT COLOR="#CC3333">Level</FONT></H3>intermediate<BR>
<H3><FONT COLOR="#CC3333">Location</FONT></H3>
</B><A HREF="../../../src/ksp/pc/impls/deflation/deflation.c.html#PCDEFLATION">src/ksp/pc/impls/deflation/deflation.c</A>
<BR><A HREF="./index.html">Index of all PC routines</A>
<BR><A HREF="../../index.html">Table of Contents for all manual pages</A>
<BR><A HREF="../singleindex.html">Index of all manual pages</A>
</BODY></HTML>
|