petsc-3.3-p7 2013-05-11

PCMG

Use multigrid preconditioning. This preconditioner requires you provide additional information about the coarser grid matrices and restriction/interpolation operators.

Options Database Keys

-pc_mg_levels <nlevels> - number of levels including finest
-pc_mg_cycles v or w- . -pc_mg_smoothup <n> - number of smoothing steps after interpolation
-pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
-pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default
-pc_mg_log - log information about time spent on each level of the solver
-pc_mg_monitor - print information on the multigrid convergence
-pc_mg_galerkin - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
-pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
-pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices to the Socket viewer for reading from MATLAB.
-pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices to the binary output file called binaryoutput

Notes: By default this uses GMRES on the fine grid smoother so this should be used with KSPFGMRES or the smoother changed to not use GMRES

When run with a single level the smoother options are used on that level NOT the coarse grid solver options

See Also

PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()

Level:intermediate
Location:
src/ksp/pc/impls/mg/mg.c
Index of all PC routines
Table of Contents for all manual pages
Index of all manual pages

Examples

src/ksp/ksp/examples/tutorials/ex42.c.html