ksp_monitor and hypre/boomeramg

Barry Smith bsmith at mcs.anl.gov
Tue Nov 11 21:26:27 CST 2008


   You are not doing anything wrong. If you run with the additional
option -pc_hypre_boomeramg_tol <rtol> where you use the same rtol
here as you use for -ksp_rtol <rtol> (default is 1.e-5) then you will
get the behavior you expect.

   This is a "feature"/"design bug" in our 2.3.3 interface to HYPRE's  
BoomerAMG.
It uses a default -pc_hypre_boomeramg_tol of 0.0 which causes it to
run the BoomerAMG cycle for the full default 10,000 cycles when run
with Richardson.

   In src/ksp/pc/impls/hypre/hypre.c line 608 you will find
   ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,PetscMin(rtol,jac- 
 >tol));CHKERRQ(ierr);
it has been changed in petsc-dev to
   ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,rtol);CHKERRQ(ierr);
this causes it to only run BoomerAMG cycles until it gets below the  
KSP rtol.

    Barry

Note: running Richardson with and without monitoring is different code
for many preconditioners. This is because monitoring reguires the || 
residual||
at each iteration. Many Richardson application of preconditioners like  
SOR,
and BoomerAMG never compute the residual hence cannot compute || 
residual||.
These preconditioners have a separate PCApplyRichardson() method rather
then just calling PCApply() and computing the residual at each  
iteration.


On Nov 11, 2008, at 6:16 PM, Shao-Ching Huang wrote:

> Hi,
>
> I am using PETSc-2.3.3p15, compiled with Hypre 2.4.0b, to solve a
> Poisson equation (100x100 mesh) using Hypre's boomeramg. The same
> matrix equation has been solved by PETSc's CG or GMRES successfully.
>
> When I run with "-ksp_monitor" (see below), the iterations show
> convergence and everything seems fine.
>
> When I remove the "ksp_monitor" flag (keeping all other flags the
> same), it takes relatively long time to finish. It appears to me (just
> a guess) that the KSP runs through the default ksp_max_it=10000
> iterations without looking at the residual level (even when
> ksp_rtol=1.e-5).
>
> My question is: why does "-ksp_monitor" change the convergence
> criteria/behavior? Or am I doing something wrong?
>
> Some more details are attached below. Thanks,
>
> Shao-Ching
>
> ---------------------------------
>
> In the following, "solve time" is in seconds.
>
> #### with "-ksp_monitor"
>
> $ mpiexec -np 1 -machinefile mach ./ibm  -ksp_type richardson \
> -pc_type hypre -pc_hypre_type boomeramg -ksp_monitor |grep  
> poisson_solve
>
> [poisson_solve] solve time = 0.043164
>
> The convergence looks like:
>
>  0 KSP Residual norm 1.263975733203e-01
>  1 KSP Residual norm 4.825832685316e-03
>  2 KSP Residual norm 1.980816153418e-04
>  3 KSP Residual norm 9.187303727632e-06
>  4 KSP Residual norm 5.022156966290e-07
>
> #### without "-ksp_monitor"
>
> $ mpiexec -np 1 -machinefile mach ./ibm  -ksp_type richardson \
> -pc_type hypre -pc_hypre_type boomeramg  |grep poisson_solve
>
> [poisson_solve] solve time = 66.463
>
>
> #### without "-ksp_monitor" but set "-ksp_max_it 10"
>
> $ mpiexec -np 1 -machinefile mach ./ibm  -ksp_type richardson \
> -pc_type hypre -pc_hypre_type boomeramg  -ksp_max_it 10 |grep  
> poisson_solve
>
> [poisson_solve] solve time = 0.0682499
>
> --------------------------------------
>
> ksp_view:
>
> KSP Object:
>  type: richardson
>    Richardson: damping factor=1
>  maximum iterations=10000, initial guess is zero
>  tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
>  left preconditioning
> PC Object:
>  type: hypre
>    HYPRE BoomerAMG preconditioning
>    HYPRE BoomerAMG: Cycle type V
>    HYPRE BoomerAMG: Maximum number of levels 25
>    HYPRE BoomerAMG: Maximum number of iterations PER hypre call 1
>    HYPRE BoomerAMG: Convergence tolerance PER hypre call 0
>    HYPRE BoomerAMG: Threshold for strong coupling 0.25
>    HYPRE BoomerAMG: Interpolation truncation factor 0
>    HYPRE BoomerAMG: Interpolation: max elements per row 0
>    HYPRE BoomerAMG: Number of levels of aggressive coarsening 0
>    HYPRE BoomerAMG: Number of paths for aggressive coarsening 1
>    HYPRE BoomerAMG: Maximum row sums 0.9
>    HYPRE BoomerAMG: Sweeps down         1
>    HYPRE BoomerAMG: Sweeps up           1
>    HYPRE BoomerAMG: Sweeps on coarse    1
>    HYPRE BoomerAMG: Relax down          symmetric-SOR/Jacobi
>    HYPRE BoomerAMG: Relax up            symmetric-SOR/Jacobi
>    HYPRE BoomerAMG: Relax on coarse     Gaussian-elimination
>    HYPRE BoomerAMG: Relax weight  (all)      1
>    HYPRE BoomerAMG: Outer relax weight (all) 1
>    HYPRE BoomerAMG: Using CF-relaxation
>    HYPRE BoomerAMG: Measure type        local
>    HYPRE BoomerAMG: Coarsen type        Falgout
>    HYPRE BoomerAMG: Interpolation type  classical
>  linear system matrix = precond matrix:
>  Matrix Object:
>    type=mpiaij, rows=10000, cols=10000
>    total: nonzeros=90000, allocated nonzeros=90000
>      not using I-node (on process 0) routines
>
> --------------------------------------
>
> $ uname -a
> Linux n125 2.6.18-92.1.6.el5 #1 SMP Wed Jun 25 13:45:47 EDT 2008
> x86_64 x86_64 x86_64 GNU/Linux
>
>
>
>
>




More information about the petsc-users mailing list