My settings are overriden when I call the solver more than once?

Barry Smith bsmith at mcs.anl.gov
Tue Apr 18 12:16:47 CDT 2006


    What controls the iteration max and tolerance for Richardson 
solver is -inner_ksp_max_its <its> and -inner_ksp_rtol <rtol>
This is true for any preconditioner including hypre boomeramg.
The defaults are 10000 and 1.e-5.

The problem is that the options
-inner_pc_hypre_boomeramg_max_iter 4
-inner_pc_hypre_boomeramg_tol 1.0e-06
sure look like THEY should control this. Hence the confusion.

Here is the issue, pc_hypre_boomeramg_max_iter and pc_hypre_boomeramg_tol are
suppose to control (with PETSc) how many iterations boomeramg does WITHIN a 
single call to ITS solve. Say we used a ksp of gmres and a pc_hypre_boomeramg_max_iter 4
this means we are using preconditioned gmres with a preconditioner of FOUR V cycles
of boomeramg (this is different then using preconditioned gmres with a preconditioner of
ONE cycle). Completely possible and maybe a reasonable thing to do.

The error in PETSc is that paragraph three above does NOT hold for Richardsons method
(but it does hold for all other KSP methods). This is because Richardson's method
has a "special method" PCApplyRichardson() that avoids the overhead of explicitly
applying Richardson's method (see src/ksp/ksp/impls/rich/rich.c). I did not handle
this properly when I wrote the PCApplyRichardson_BoomerAMG. I will fix this error.

Anyways, the short answer is when using Richardson with hypre 
use -inner_ksp_max_its <its> and -inner_ksp_rtol <rtol> to control 
the iterations and tolerance, NOT 
-inner_pc_hypre_boomeramg_max_iter
and -inner_pc_hypre_boomeramg_tol
I will attempt to make the docs clearer.


    Barry





On Tue, 18 Apr 2006, Sh.M wrote:

> Hi all,
>
> I am solving a matrix with boomerAMG, and I set the boomerAMG settings
> thru the console. The matrix(it is a preconditioner matrix) is solved
> several times... I set the boomerAMG preconditioner/solver to have a
> convergence tolerance of 1.0e-06 and max number of iterations 4. I print out boomerAMG status to
> verify that it works as it should.. But for some reason it doesnt..
>
> Here is a print from boomerAMG when it starts, I guess when it constructs
> the boomerAMG solver/preconditioner:
>
> .
> .
> .
> .
> .
> .
> BoomerAMG SOLVER PARAMETERS:
>
>  Maximum number of cycles:         4
>  Stopping Tolerance:               1.000000e-06
>  Cycle type (1 = V, 2 = W, etc.):  1
>
>  Relaxation Parameters:
>   Visiting Grid:                     fine  down   up  coarse
>            Number of partial sweeps:   1     1    1     1
>   Type 0=Jac, 1=GS, 3=Hybrid 9=GE:     3     3    3     9
>   Point types, partial sweeps (1=C, -1=F):
>                               Finest grid:   1  -1
>                  Pre-CG relaxation (down):   1  -1
>                   Post-CG relaxation (up):  -1   1
>                             Coarsest grid:   0
>
>
>
> Immediately after the above message, wich is when it starts to
> solve, it changes its parameters to this:
>
>
> BoomerAMG SOLVER PARAMETERS:
>
>  Maximum number of cycles:         10000
>  Stopping Tolerance:               1.000000e-05
>  Cycle type (1 = V, 2 = W, etc.):  1
>
>  Relaxation Parameters:
>   Visiting Grid:                     fine  down   up  coarse
>            Number of partial sweeps:   1     1    1     1
>   Type 0=Jac, 1=GS, 3=Hybrid 9=GE:     3     3    3     9
>   Point types, partial sweeps (1=C, -1=F):
>                               Finest grid:   1  -1
>                  Pre-CG relaxation (down):   1  -1
>                   Post-CG relaxation (up):  -1   1
>                             Coarsest grid:   0
> .
> .
> .
> .
> .
> .
>
>
> As you see my settings have been changed.
>
> I have used boomerAMG in the past, not exactly this way.,... but my
> settings were preserved before and after solve.
>
> here is what I run from the command line:
>
> mprun -np 1 petscSolver
> -a ../hypre/data/nr_hvl40.csr
> -b ../hypre/data/nr_hvl40.csr.blockMatrix -inner_ksp_type richardson
> -inner_pc_type hypre -inner_pc_hypre_type boomeramg
> -inner_pc_hypre_boomeramg_max_iter 4
> -inner_pc_hypre_boomeramg_tol 1.0e-06
> -ksp_monitor
> -ksp_type gmres
> -inner_pc_hypre_boomeramg_print_statistics
>
>
> And a piece of code:
>
> .
> .
> .
> .
> .
> .
> BLOCK_PC  user_pc;
> ,
> ,
> ,
> ,
> ,
> .
> KSPSetOptionsPrefix(user_pc.block_solver,"inner_");
> KSPSetOperators(user_pc.block_solver,user_pc.M,user_pc.M,DIFFERENT_NONZERO_PATTERN);
> KSPSetFromOptions(user_pc.block_solver);
>
>
> I am running on solaris64.
>
> With best regards, Shaman Mahmoudi
>
>




More information about the petsc-users mailing list