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
|
<!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/PCML.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<TITLE>PCML</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
<div id="version" align=right><b>petsc-3.7.5 2017-01-01</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.7.5 v3.7.5 docs/manualpages/PC/PCML.html "><small>Report Typos and Errors</small></a></div>
<A NAME="PCML"><H1>PCML</H1></A>
Use algebraic multigrid preconditioning. This preconditioner requires you provide fine grid discretization matrix. The coarser grid matrices and restriction/interpolation operators are computed by ML, with the matrices coverted to PETSc matrices in aij format and the restriction/interpolation operators wrapped as PETSc shell matrices.
<H3><FONT COLOR="#CC3333">Options Database Key</FONT></H3>
<H3><FONT COLOR="#CC3333">Multigrid options(inherited)</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-pc_mg_cycles <1>: 1 for V cycle, 2 for W</B></TD><TD>- cycle (MGSetCycles)
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-pc_mg_smoothup <1>: Number of post</B></TD><TD>- smoothing steps (MGSetNumberSmoothUp)
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-pc_mg_smoothdown <1>: Number of pre</B></TD><TD>- smoothing steps (MGSetNumberSmoothDown)
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade</B></TD><TD>- ML options:
+ -pc_ml_PrintLevel <0>: Print level (ML_Set_PrintLevel)
</TD></TR></TABLE>
<DT><B>-pc_ml_maxNlevels <10>: Maximum number of levels (None)</B> -. -pc_ml_maxCoarseSize <1>: Maximum coarsest mesh size (ML_Aggregate_Set_MaxCoarseSize)
<br>
<DT><B>-pc_ml_CoarsenScheme <Uncoupled>: (one of) Uncoupled Coupled MIS METIS</B> -. -pc_ml_DampingFactor <1.33333>: P damping factor (ML_Aggregate_Set_DampingFactor)
<br>
<DT><B>-pc_ml_Threshold <0>: Smoother drop tol (ML_Aggregate_Set_Threshold)</B> -. -pc_ml_SpectralNormScheme_Anorm <false>: Method used for estimating spectral radius (ML_Set_SpectralNormScheme_Anorm)
<br>
<DT><B>-pc_ml_repartition <false>: Allow ML to repartition levels of the heirarchy (ML_Repartition_Activate)</B> -. -pc_ml_repartitionMaxMinRatio <1.3>: Acceptable ratio of repartitioned sizes (ML_Repartition_Set_LargestMinMaxRatio)
<br>
<DT><B>-pc_ml_repartitionMinPerProc <512>: Smallest repartitioned size (ML_Repartition_Set_MinPerProc)</B> -. -pc_ml_repartitionPutOnSingleProc <5000>: Problem size automatically repartitioned to one processor (ML_Repartition_Set_PutOnSingleProc)
<br>
<DT><B>-pc_ml_repartitionType <Zoltan>: Repartitioning library to use (ML_Repartition_Set_Partitioner)</B> -. -pc_ml_repartitionZoltanScheme <RCB>: Repartitioning scheme to use (None)
<br>
<DT><B>-pc_ml_Aux <false>: Aggregate using auxiliary coordinate</B> -based laplacian (None)
<br>
<DT><B>-pc_ml_AuxThreshold <0.0>: Auxiliary smoother drop tol (None)</B> -
<br>
</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>, <A HREF="../PC/PCMGType.html#PCMGType">PCMGType</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>(), MPSetCycles(), <A HREF="../PC/PCMGSetNumberSmoothDown.html#PCMGSetNumberSmoothDown">PCMGSetNumberSmoothDown</A>(),
<A HREF="../PC/PCMGSetNumberSmoothUp.html#PCMGSetNumberSmoothUp">PCMGSetNumberSmoothUp</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>(),
PCMGSetCycleTypeOnLevel(), <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><P><B><FONT COLOR="#CC3333">Level:</FONT></B>intermediate
<BR><FONT COLOR="#CC3333">Location:</FONT></B><A HREF="../../../src/ksp/pc/impls/ml/ml.c.html#PCML">src/ksp/pc/impls/ml/ml.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>
|