[petsc-users] Seting advanced smoothers for BoomerAMG through PETSc

Barry Smith bsmith at mcs.anl.gov
Wed May 6 07:27:23 CDT 2015


  Eike,

   You'll need to edit src/ksp/pc/impls/hypre/hypre.c and add an option to set this (where you see all the other BoomerAMG options). Then run make in
the root PETSc directory.

   If you do this in the development copy of PETSc and make a pull request (see https://bitbucket.org/petsc/petsc/wiki/Home Contributing to PETSc) then all future users can benefit from your additions.

  Barry

> On May 6, 2015, at 3:14 AM, Eike Mueller <E.Mueller at bath.ac.uk> wrote:
> 
> Dear PETSc,
> 
> I’m trying to use hypre’s BoomerAMG preconditioner through PETSc. However, I would like to use a non-standard smoother, namely an incomplete LU factorisation, which is implemented in hypre’s Euclid (that’s because the problem is highly anisotropic, and GAMG worked very weill with ILU smoother). In native hypre I can set this with HYPRE_BoomerAMGSetSmoothType(). However, I can’t see how I can pass this through the PETSc interface, since the only boomeramg options I can set are listed below. As far as I can see they allow me to change the pointwise smoother via HYPRE_BoomerAMGSetRelaxType().
> 
> Does anyone have any advice on what I can do to use Euclid as a smoother via setting the correct PETSc options?
> 
> Thank you very much,
> 
> Eike
> 
> eikemueller at 138-38-173-198 $ ./ex6 -f -f helmholtz-sphere.dat  -ksp_type cg -ksp_convergence_test skip -ksp_max_it 2  -ksp_monitor -table -pc_hypre_type boomeramg -pc_type hypre -log_summary -ksp_view -help | grep boomeramg
>   -pc_hypre_type <boomeramg> (choose one of) pilut parasails boomeramg ams (PCHYPRESetType)
>   -pc_hypre_boomeramg_cycle_type <V> (choose one of) V W (None)
>   -pc_hypre_boomeramg_max_levels <25>: Number of levels (of grids) allowed (None)
>   -pc_hypre_boomeramg_max_iter <1>: Maximum iterations used PER hypre call (None)
>   -pc_hypre_boomeramg_tol <0>: Convergence tolerance PER hypre call (0.0 = use a fixed number of iterations) (None)
>   -pc_hypre_boomeramg_truncfactor <0>: Truncation factor for interpolation (0=no truncation) (None)
>   -pc_hypre_boomeramg_P_max <0>: Max elements per row for interpolation operator (0=unlimited) (None)
>   -pc_hypre_boomeramg_agg_nl <0>: Number of levels of aggressive coarsening (None)
>   -pc_hypre_boomeramg_agg_num_paths <1>: Number of paths for aggressive coarsening (None)
>   -pc_hypre_boomeramg_strong_threshold <0.25>: Threshold for being strongly connected (None)
>   -pc_hypre_boomeramg_max_row_sum <0.9>: Maximum row sum (None)
>   -pc_hypre_boomeramg_grid_sweeps_all <1>: Number of sweeps for the up and down grid levels (None)
>   -pc_hypre_boomeramg_grid_sweeps_down <1>: Number of sweeps for the down cycles (None)
>   -pc_hypre_boomeramg_grid_sweeps_up <1>: Number of sweeps for the up cycles (None)
>   -pc_hypre_boomeramg_grid_sweeps_coarse <1>: Number of sweeps for the coarse level (None)
>   -pc_hypre_boomeramg_relax_type_all <symmetric-SOR/Jacobi> (choose one of) Jacobi sequential-Gauss-Seidel seqboundary-Gauss-Seidel SOR/Jacobi backward-SOR/Jacobi  symmetric-SOR/Jacobi  l1scaled-SOR/Jacobi Gaussian-elimination      CG Chebyshev FCF-Jacobi l1scaled-Jacobi (None)
>   -pc_hypre_boomeramg_relax_type_down <symmetric-SOR/Jacobi> (choose one of) Jacobi sequential-Gauss-Seidel seqboundary-Gauss-Seidel SOR/Jacobi backward-SOR/Jacobi  symmetric-SOR/Jacobi  l1scaled-SOR/Jacobi Gaussian-elimination      CG Chebyshev FCF-Jacobi l1scaled-Jacobi (None)
>   -pc_hypre_boomeramg_relax_type_up <symmetric-SOR/Jacobi> (choose one of) Jacobi sequential-Gauss-Seidel seqboundary-Gauss-Seidel SOR/Jacobi backward-SOR/Jacobi  symmetric-SOR/Jacobi  l1scaled-SOR/Jacobi Gaussian-elimination      CG Chebyshev FCF-Jacobi l1scaled-Jacobi (None)
>   -pc_hypre_boomeramg_relax_type_coarse <Gaussian-elimination> (choose one of) Jacobi sequential-Gauss-Seidel seqboundary-Gauss-Seidel SOR/Jacobi backward-SOR/Jacobi  symmetric-SOR/Jacobi  l1scaled-SOR/Jacobi Gaussian-elimination      CG Chebyshev FCF-Jacobi l1scaled-Jacobi (None)
>   -pc_hypre_boomeramg_relax_weight_all <1>: Relaxation weight for all levels (0 = hypre estimates, -k = determined with k CG steps) (None)
>   -pc_hypre_boomeramg_relax_weight_level <1>: Set the relaxation weight for a particular level (weight,level) (None)
>   -pc_hypre_boomeramg_outer_relax_weight_all <1>: Outer relaxation weight for all levels (-k = determined with k CG steps) (None)
>   -pc_hypre_boomeramg_outer_relax_weight_level <1>: Set the outer relaxation weight for a particular level (weight,level) (None)
>   -pc_hypre_boomeramg_no_CF: <FALSE> Do not use CF-relaxation (None)
>   -pc_hypre_boomeramg_measure_type <local> (choose one of) local global (None)
>   -pc_hypre_boomeramg_coarsen_type <Falgout> (choose one of) CLJP Ruge-Stueben  modifiedRuge-Stueben   Falgout  PMIS  HMIS (None)
>   -pc_hypre_boomeramg_interp_type <classical> (choose one of) classical   direct multipass multipass-wts ext+i ext+i-cc standard standard-wts   FF FF1 (None)
>   -pc_hypre_boomeramg_print_statistics: Print statistics (None)
>   -pc_hypre_boomeramg_print_debug: Print debug information (None)
>   -pc_hypre_boomeramg_nodal_coarsen: <FALSE> HYPRE_BoomerAMGSetNodal() (None)
>   -pc_hypre_boomeramg_nodal_relaxation: <FALSE> Nodal relaxation via Schwarz (None)
> -pc_hypre_type boomeramg
> 
> --
> 
> Dr Eike Hermann Mueller
> Lecturer in Scientific Computing
> 
> Department of Mathematical Sciences
> University of Bath
> Bath BA2 7AY, United Kingdom
> 
> +44 1225 38 6241
> e.mueller at bath.ac.uk
> http://people.bath.ac.uk/em459/
> 
> 
> 
> 
> 



More information about the petsc-users mailing list