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
- number of levels including finest-pc_mg_cycle_type <v,w> - provide the cycle desired
-pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
-pc_mg_log - log information about time spent on each level of the solver
-pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes)
-pc_mg_galerkin <both,pmat,mat,none> - 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 a
PETSCVIEWERSOCKET
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#
The Krylov solver (if any) and preconditioner (smoother) and their parameters are controlled from the options database with the standard
options database keywords prefixed with -mg_levels_
to affect all the levels but the coarsest, which is controlled with -mg_coarse_
,
and the finest where -mg_fine_
can override -mg_levels_
. One can set different preconditioners etc on specific levels with the prefix
-mg_levels_n_
where n
is the level number (zero being the coarse level. For example
-mg_levels_ksp_type gmres -mg_levels_pc_type bjacobi -mg_coarse_pc_type svd -mg_levels_2_pc_type sor
These options also work for controlling the smoothers etc inside PCGAMG
If one uses a Krylov method such KSPGMRES
or KSPCG
as the smoother then one must use KSPFGMRES
, KSPGCR
, or KSPRICHARDSON
as the outer Krylov method
When run with a single level the smoother options are used on that level NOT the coarse grid solver options
When run with KSPRICHARDSON
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.
See Also#
Multigrid Preconditioners, PCCreate()
, PCSetType()
, PCType
, PC
, PCMGType
, PCEXOTIC
, PCGAMG
, PCML
, PCHYPRE
PCMGSetLevels()
, PCMGGetLevels()
, PCMGSetType()
, PCMGSetCycleType()
,
PCMGSetDistinctSmoothUp()
, PCMGGetCoarseSolve()
, PCMGSetResidual()
, PCMGSetInterpolation()
,
PCMGSetRestriction()
, PCMGGetSmoother()
, PCMGGetSmootherUp()
, PCMGGetSmootherDown()
,
PCMGSetCycleTypeOnLevel()
, PCMGSetRhs()
, PCMGSetX()
, PCMGSetR()
,
PCMGSetAdaptCR()
, PCMGGetAdaptInterpolation()
, PCMGSetGalerkin()
, PCMGGetAdaptCoarseSpaceType()
, PCMGSetAdaptCoarseSpaceType()
Level#
intermediate
Location#
Examples#
src/ksp/ksp/tutorials/ex36.cxx
src/snes/tutorials/ex20.c
src/dm/impls/stag/tutorials/ex4.c
src/ksp/ksp/tutorials/ex42.c
src/ksp/ksp/tutorials/ex65.c
src/ksp/ksp/tutorials/ex35.cxx
Index of all PC routines
Table of Contents for all manual pages
Index of all manual pages