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
|
<!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/PCMG.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<TITLE>PCMG</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
<div id="version" align=right><b>petsc-3.10.3 2018-12-18</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.10.3 v3.10.3 docs/manualpages/PC/PCMG.html "><small>Report Typos and Errors</small></a></div>
<A NAME="PCMG"><H1>PCMG</H1></A>
Use multigrid preconditioning. This preconditioner requires you provide additional information about the coarser grid matrices and restriction/interpolation operators.
<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_mg_levels <nlevels> </B></TD><TD>- number of levels including finest
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-pc_mg_cycle_type <v,w> </B></TD><TD>- provide the cycle desired
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-pc_mg_type <additive,multiplicative,full,kaskade> </B></TD><TD>- multiplicative is the default
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-pc_mg_log </B></TD><TD>- log information about time spent on each level of the solver
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-pc_mg_distinct_smoothup </B></TD><TD>- configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes)
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-pc_mg_galerkin <both,pmat,mat,none> </B></TD><TD>- use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-pc_mg_multiplicative_cycles </B></TD><TD>- number of cycles to use as the preconditioner (defaults to 1)
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-pc_mg_dump_matlab </B></TD><TD>- dumps the matrices for each level and the restriction/interpolation matrices
to the Socket viewer for reading from MATLAB.
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-pc_mg_dump_binary </B></TD><TD>- dumps the matrices for each level and the restriction/interpolation matrices
to the binary output file called binaryoutput
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#CC3333">Notes</FONT></H3>
If one uses a Krylov method such GMRES or CG as the smoother then one must use <A HREF="../KSP/KSPFGMRES.html#KSPFGMRES">KSPFGMRES</A>, <A HREF="../KSP/KSPGCR.html#KSPGCR">KSPGCR</A>, or <A HREF="../KSP/KSPRICHARDSON.html#KSPRICHARDSON">KSPRICHARDSON</A> as the outer Krylov method
<P>
When run with a single level the smoother options are used on that level NOT the coarse grid solver options
<P>
When run with <A HREF="../KSP/KSPRICHARDSON.html#KSPRICHARDSON">KSPRICHARDSON</A> the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This
is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing
(because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the
residual is computed at the end of each cycle.
<P>
<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>, <A HREF="../PC/PCMGType.html#PCMGType">PCMGType</A>, <A HREF="../PC/PCEXOTIC.html#PCEXOTIC">PCEXOTIC</A>, <A HREF="../PC/PCGAMG.html#PCGAMG">PCGAMG</A>, <A HREF="../PC/PCML.html#PCML">PCML</A>, <A HREF="../PC/PCHYPRE.html#PCHYPRE">PCHYPRE</A>
<BR><A HREF="../PC/PCMGSetLevels.html#PCMGSetLevels">PCMGSetLevels</A>(), <A HREF="../PC/PCMGGetLevels.html#PCMGGetLevels">PCMGGetLevels</A>(), <A HREF="../PC/PCMGSetType.html#PCMGSetType">PCMGSetType</A>(), <A HREF="../PC/PCMGSetCycleType.html#PCMGSetCycleType">PCMGSetCycleType</A>(),
<A HREF="../PC/PCMGSetDistinctSmoothUp.html#PCMGSetDistinctSmoothUp">PCMGSetDistinctSmoothUp</A>(), <A HREF="../PC/PCMGGetCoarseSolve.html#PCMGGetCoarseSolve">PCMGGetCoarseSolve</A>(), <A HREF="../PC/PCMGSetResidual.html#PCMGSetResidual">PCMGSetResidual</A>(), <A HREF="../PC/PCMGSetInterpolation.html#PCMGSetInterpolation">PCMGSetInterpolation</A>(),
<A HREF="../PC/PCMGSetRestriction.html#PCMGSetRestriction">PCMGSetRestriction</A>(), <A HREF="../PC/PCMGGetSmoother.html#PCMGGetSmoother">PCMGGetSmoother</A>(), <A HREF="../PC/PCMGGetSmootherUp.html#PCMGGetSmootherUp">PCMGGetSmootherUp</A>(), <A HREF="../PC/PCMGGetSmootherDown.html#PCMGGetSmootherDown">PCMGGetSmootherDown</A>(),
<A HREF="../PC/PCMGSetCycleTypeOnLevel.html#PCMGSetCycleTypeOnLevel">PCMGSetCycleTypeOnLevel</A>(), <A HREF="../PC/PCMGSetRhs.html#PCMGSetRhs">PCMGSetRhs</A>(), <A HREF="../PC/PCMGSetX.html#PCMGSetX">PCMGSetX</A>(), <A HREF="../PC/PCMGSetR.html#PCMGSetR">PCMGSetR</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/mg/mg.c.html#PCMG">src/ksp/pc/impls/mg/mg.c</A>
<P><H3><FONT COLOR="#CC3333">Examples</FONT></H3>
<A HREF="../../../src/ksp/ksp/examples/tutorials/ex42.c.html">src/ksp/ksp/examples/tutorials/ex42.c.html</A><BR>
<A HREF="../../../src/snes/examples/tutorials/ex12.c.html">src/snes/examples/tutorials/ex12.c.html</A><BR>
<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>
|