[petsc-users] Memory and Speed Issue of using MG as preconditioner

Dave May dave.mayhem23 at gmail.com
Wed Nov 6 12:40:16 CST 2013


To speed up the smoother you could:
1) do less smoothing, but this might incur needing to do more KSP iterations
2) run on more cores
3) if you are happy with cheby/jacobi, you could write a matrix free
implementation for your discretisation which gets used on the fine level.
Depending on your discretisation, this might improve the speed and
scalability.

As I said, the prefix for the solver is mentioned in the result of -ksp_view
e.g.
PC Object:    (*mg_coarse_*)     32 MPI processes
      type: redundant

The string in brackets (bolded) is the solver prefix.
So if you want to change the coarse grid solver, do this

-mg_coarse_ksp_type XXXX
-mg_coarse_pc_type YYYY


Cheers,
  Dave



On 6 November 2013 19:35, Alan Z. Wei <zhenglun.wei at gmail.com> wrote:

>  Thanks Dave,
>     I further simulated the problem with -pc_mg_log and output these files
> in the attachments.
>     I found that the smoothing process of the last level always consumes
> the most time, i.e. 'MGSmooth Level 5' in out-level5 and "MGSmooth Level 2'
> in out-level2. However, as I tested several other -mg_levels_pc_type such
> as 'bjacobi', 'asm' etc. The default one, which is 'jacobi', actually works
> the best. Therefore, I decide to keep using it. However, do you have any
> suggestions to speed up this smoothing process other than changing
> -mg_levels_pc_type?
>      Also, as you suggested to change -mg_levels_ksp_type, it does not
> influence much if replacing 'chebyshev' by 'cg'. However, this part never
> change while I modify '-mg_levels_ksp_type':
> PC Object:    (mg_coarse_)     32 MPI processes
>       type: redundant
>         Redundant preconditioner: First (color=0) of 32 PCs follows
>       KSP Object:      (mg_coarse_redundant_)       1 MPI processes
>         type: preonly
>         maximum iterations=10000, initial guess is zero
>         tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
>         left preconditioning
>         using NONE norm type for convergence test
>       PC Object:      (mg_coarse_redundant_)       1 MPI processes
>         type: lu
>           LU: out-of-place factorization
>       As you mentioned that the redundant LU for the coarse grid solver
> primarily cause the large memory request for the 2-level case. How could I
> change the coarse grid solver to reduce the memory requirement or speed up
> the solver.
>
> thanks again,
> Alan
>
>      Hey Alan,
>
>  1/ One difference in the memory footprint is likely coming from your
> coarse grid solver which is redundant LU.
>  The 2 level case has a coarse grid problem with 70785 unknowns whilst the
> 5 level case has a coarse grid problem with 225 unknowns.
>
> 2/ The solve time difference will be affected by your coarse grid size.
> Add the command line argument
>   -pc_mg_log
> to profile the setup time spent on the coarse grid and all other levels.
>  See
>   http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCMG.html
>
> 3/ You can change the smoother on all levels by using the command line
> argument with the appropriate prefix, eg
>    -mg_levels_ksp_type cg
>  Note the prefix is displayed in the result of -ksp_view
>
>  Also, your mesh size can be altered at run time using arguments like
>  -da_grid_x 5
>  You shouldn't have to modify the source code each time
> See
>
> http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMDACreate3d.html
>
>
>  Cheers,
>    Dave
>
>
> On 6 November 2013 04:21, Alan <zhenglun.wei at gmail.com> wrote:
>
>> Dear all,
>> I hope you're having a nice day.
>> Recently, I came across a problem on using MG as preconditioner.
>> Basically, to achieve the same finest mesh with pc_type = mg, the memory
>> usage for -da_refine 2 is much more than that for -da_refine 5. To my
>> limited knowledge, more refinement should consume more memory, which is
>> contradict to the behavior of pc_type = mg in PETsc.
>> Here, I provide two output files. They are all from
>> /src/ksp/ksp/example/tutorial/ex45.c with 32 processes.
>> The execute file for out-level2 is
>> mpiexec -np 32 ./ex45 -pc_type mg -ksp_type cg -da_refine 2
>> -pc_mg_galerkin -ksp_rtol 1.0e-7 -mg_levels_pc_type jacobi
>> -mg_levels_ksp_type chebyshev -dm_view -log_summary -pc_mg_log
>> -pc_mg_monitor -ksp_view -ksp_monitor > out &
>> and in ex45.c, KSPCreate is changed as:
>> ierr =
>>
>> DMDACreate3d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE,DMDA_BOUNDARY_NONE,DMDA_BOUNDARY_NONE,DMDA_STENCIL_STAR,-65,-33,-33,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,1,1,0,0,0,&da);CHKERRQ(ierr);
>> On the other hand, the execute file for out-level5 is
>> mpiexec -np 32 ./ex45 -pc_type mg -ksp_type cg -da_refine 5
>> -pc_mg_galerkin -ksp_rtol 1.0e-7 -mg_levels_pc_type jacobi
>> -mg_levels_ksp_type chebyshev -dm_view -log_summary -pc_mg_log
>> -pc_mg_monitor -ksp_view -ksp_monitor > out &
>> and in ex45.c, KSPCreate is changed as:
>> ierr =
>>
>> DMDACreate3d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE,DMDA_BOUNDARY_NONE,DMDA_BOUNDARY_NONE,DMDA_STENCIL_STAR,-9,-5,-5,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,1,1,0,0,0,&da);CHKERRQ(ierr);
>> In summary, the final finest meshes obtained for both cases are
>> 257*129*129 as documented in both files. However, the out-level2 shows
>> that the Matrix requested 822871308 memory while out-level5 only need
>> 36052932.
>> Furthermore, although the total iterations for KSP solver are shown as 5
>> times in both files. the wall time elapsed for out-level2 is around
>> 150s, while out-level5 is about 4.7s.
>> At last, there is a minor question. In both files, under 'Down solver
>> (pre-smoother) on level 1' and 'Down solver (pre-smoother) on level 2',
>> the type of "KSP Object: (mg_levels_1_est_)" and "KSP Object:
>> (mg_levels_2_est_)" are all 'gmres'. Since I'm using uniformly Cartesian
>> mesh, would it be helpful to speed up the solver if the 'gmres' is
>> replaced by 'cg' here? If so, which PETSc option can change the type of
>> KSP object.
>>
>> sincerely appreciate,
>> Alan
>>
>>
>>
>>
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20131106/c21870d7/attachment.html>


More information about the petsc-users mailing list