Poor performance with BoomerAMG?

Barry Smith bsmith at mcs.anl.gov
Tue Feb 19 15:56:18 CST 2008


    BoomerAMG works like a charm. Likely you forgot the -pc_hypre_type  
boomeramg

    Hmm, I think I'll change the default solver to boomeramg

    Barry

barry-smiths-macbook-pro-17:ksp/examples/tutorials] bsmith% ./ex1f - 
ksp_monitor -pc_type hypre -pc_hypre_type boomeramg -m 513 -n 513 - 
ksp_type richardson -ksp_view
  p=           1
   0 KSP Residual norm 4.213878296084e+03
   1 KSP Residual norm 2.135189837330e+02
   2 KSP Residual norm 1.225934028865e+01
   3 KSP Residual norm 7.255859884400e-01
   4 KSP Residual norm 4.353504737395e-02
   5 KSP Residual norm 2.643035146258e-03
   6 KSP Residual norm 1.628271972668e-04
KSP Object:
   type: richardson
     Richardson: damping factor=1
   maximum iterations=10000, initial guess is zero
   tolerances:  relative=1e-07, 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=seqaij, rows=263169, cols=263169
     total: nonzeros=1313793, allocated nonzeros=1315845
       not using I-node routines
  Iterations:           7

[barry-smiths-macbook-pro-17:ksp/examples/tutorials] bsmith% ./ex1f - 
ksp_monitor -pc_type hypre -pc_hypre_type boomeramg -m 513 -n 513 - 
ksp_type gmres -ksp_view
  p=           1
   0 KSP Residual norm 4.213878296084e+03
   1 KSP Residual norm 5.272381634094e+01
   2 KSP Residual norm 8.107668116258e-01
   3 KSP Residual norm 1.807380875232e-02
   4 KSP Residual norm 4.068259191532e-04
KSP Object:
   type: gmres
     GMRES: restart=30, using Classical (unmodified) Gram-Schmidt  
Orthogonalization with no iterative refinement
     GMRES: happy breakdown tolerance 1e-30
   maximum iterations=10000, initial guess is zero
   tolerances:  relative=1e-07, 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=seqaij, rows=263169, cols=263169
     total: nonzeros=1313793, allocated nonzeros=1315845
       not using I-node routines
  Iterations:           4


On Feb 18, 2008, at 1:57 AM, knutert at stud.ntnu.no wrote:

> Thank you for the reply, Barry.
>
> The same thing happens if I use hypre with the DMMG solver.
> As you say, with hypre, the convergence is extremely slow, requiring
> a lot of iterations, 1413 iterations (1820 if I use richardson) for  
> a 257x257
> problem, while the default only needs 5.
>
> I use the same way of handling boundary conditions in the two codes.
> I've also compared the coeff matrix and rhs, and they are equal.
>
> -Knut Erik-
>
> Siterer Barry Smith <bsmith at mcs.anl.gov>:
>
>>
>>  Run with the DMMG solver with the option -pc_type hypre
>> What happens? Then run again with the additional option -ksp_type  
>> richardson
>>
>> Is hypre taking many, many iterations which is causing the slow  
>> speed?
>>
>> I expect there is something wrong with your code that does not use  
>> DMMG.
>> Be careful how you handle boundary conditions; you need to make sure
>> they have the same scaling as the other equations.
>>
>>   Barry
>>
>>
>>
>> On Feb 15, 2008, at 8:36 AM, knutert at stud.ntnu.no wrote:
>>
>>> Hi Ben,
>>>
>>> Thank you for answering. With gmres and boomeramg I get a run time  
>>> of
>>> 2s, so that is much better. However, if I increase the grid size to
>>> 513x513, I get a run time of one minute. With richardson, it  
>>> fails  to converge.
>>> LU gives 6 seconds, CG and ICC gives 7s, and the DMMG solver 3s  
>>> for  the 513x513 problem.
>>>
>>> When using the DMMG framework, I just used the default solvers.
>>> I use the Galerkin process to generate the coarse matrices for
>>> the multigrid cycle.
>>>
>>> Best,
>>> Knut
>>>
>>> Siterer Ben Tay <zonexo at gmail.com>:
>>>
>>>> Hi Knut,
>>>>
>>>> I'm currently using boomeramg to solve my poisson eqn too. I'm  
>>>> using it
>>>> on my structured C-grid. I found it to be faster than LU,  
>>>> especially as
>>>> the grid size increases. However I use it as a preconditioner with
>>>> GMRES as the solver. Have you tried this option? Although it's  
>>>> faster,
>>>> the speed increase is usually less than double. It seems to be  
>>>> worse if
>>>> there is a lot of stretching in the grid.
>>>>
>>>> Btw, your mention using the DMMG framework and it takes less than a
>>>> sec. What solver or preconditioner did you use? It's 4 times faster
>>>> than GMRES...
>>>>
>>>> thanks!
>>>>
>>>> knutert at stud.ntnu.no wrote:
>>>>> Hello,
>>>>>
>>>>> I am trying to use the hypre multigrid solver to solve a Poisson  
>>>>> equation.
>>>>> However, on a test case with grid size 257x257 it takes 40   
>>>>> seconds  to converge
>>>>> on one processor when I run with
>>>>> ./run -ksp_type richardson -pc_type hypre -pc_type_hypre boomeramg
>>>>>
>>>>> Using the DMMG framework, the same problem takes less than a  
>>>>> second,
>>>>> and the default gmres solver uses only four seconds.
>>>>>
>>>>> Am I somehow using the solver the wrong way, or is this   
>>>>> performance  expected?
>>>>>
>>>>> Regards
>>>>> Knut Erik Teigen
>>>>>
>>>>>
>>>
>>>
>>>
>
>
>




More information about the petsc-users mailing list