[petsc-users] Petsc option and scalablity problems on ML

Jed Brown jed at 59A2.org
Wed May 25 18:59:06 CDT 2011


On Thu, May 26, 2011 at 01:42, Li, Zhisong (lizs) <lizs at mail.uc.edu> wrote:

> I have been working for a while on applying ML precoditioner to improve
> performance of a KSP solving a Poisson-style problem.  So far the only
> speedup gain attributed to resetting the number of multigrid levels.  I
> wonder if all the available Petsc options to control ML are limited to those
> listed on the Petsc manual (about 10 totally). Do we have any other controls
> more than that?
>

Run with "-help |grep pc_ml_" to see the full list of options that can be
passed to ML. ML constructs the hierarchy and then you can control it
through PCMG. Relevant options are under the prefixes -pc_mg_*,
-mg_coarse_*, and -mg_levels_*.


> In ML user's guide, they state another tip to improve ML's parallel speed:
> "Instead of doing a direct solve on the coarsest level, try a few smoothing
> sweeps instead".  I don't think the available options "-pc_mg_smoothup" and
> "-pc_mg_smoothdown" correspond to this control.  Can any PETSC interface
> option do that?
>

You can control the coarse level with -mg_coarse_ksp_*, for example,

-mg_coarse_ksp_type cg -mg_coarse_ksp_max_it 3 -mg_coarse_pc_type bjacobi
-mg_coarse_sub_pc_type lu

to do 3 cycles of CG preconditioned by block Jacobi with direct subdomain
solves (trivially cheap on a coarse level) as the coarse level solver. Note
that Krylov iterations here will make the preconditioner nonlinear so you
should use -ksp_type fgmres on the outside. If you know something about the
spectrum, you could use something like -mg_coarse_ksp_type chebychev which
is linear for a fixed number of iterations.
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20110526/fd48b613/attachment.htm>


More information about the petsc-users mailing list