petsc-3.9.4 2018-09-11
Report Typos and Errors

PCHYPRE

Allows you to use the matrix element based preconditioners in the LLNL package hypre

Options Database Keys

-pc_hypre_type - One of pilut, parasails, boomeramg, ams, ads
Too many others to list, run with - pc_type hypre -pc_hypre_type XXX -help to see options for the XXX preconditioner

Notes: Apart from pc_hypre_type (for which there is PCHYPRESetType()), the many hypre options can ONLY be set via the options database (e.g. the command line or with PetscOptionsSetValue(), there are no functions to set them)

The options -pc_hypre_boomeramg_max_iter and -pc_hypre_boomeramg_tol refer to the number of iterations (V-cycles) and tolerance that boomeramg does EACH time it is called. So for example, if -pc_hypre_boomeramg_max_iter is set to 2 then 2-V-cycles are being used to define the preconditioner (-pc_hypre_boomeramg_tol should be set to 0.0 - the default - to strictly use a fixed number of iterations per hypre call). -ksp_max_it and -ksp_rtol STILL determine the total number of iterations and tolerance for the Krylov solver. For example, if -pc_hypre_boomeramg_max_iter is 2 and -ksp_max_it is 10 then AT MOST twenty V-cycles of boomeramg will be called.

Note that the option -pc_hypre_boomeramg_relax_type_all defaults to symmetric relaxation (symmetric-SOR/Jacobi), which is required for Krylov solvers like CG that expect symmetry. Otherwise, you may want to use -pc_hypre_boomeramg_relax_type_all SOR/Jacobi. If you wish to use BoomerAMG WITHOUT a Krylov method use -ksp_type richardson NOT -ksp_type preonly and use -ksp_max_it to control the number of V-cycles. (see the PETSc FAQ.html at the PETSc website under the Documentation tab).

2007-02-03 Using HYPRE-1.11.1b, the routine HYPRE_BoomerAMGSolveT and the option -pc_hypre_parasails_reuse were failing with SIGSEGV. Dalcin L.

MatSetNearNullSpace() - if you provide a near null space to your matrix it is ignored by hypre UNLESS you also use

the two options

-pc_hypre_boomeramg_nodal_coarsen <n> - where n is from 1 to 6 (see HYPRE_BOOMERAMGSetNodal())
-pc_hypre_boomeramg_vec_interp_variant <v> where v is from 1 to 3 (see HYPRE_BoomerAMGSetInterpVecVariant())- Depending on the linear system you may see the same or different convergence depending on the values you use.

See PCPFMG for access to the hypre Struct PFMG solver

See Also

PCCreate(), PCSetType(), PCType (for list of available types), PC,
PCHYPRESetType(), PCPFMG

Level

intermediate

Location

src/ksp/pc/impls/hypre/hypre.c
Index of all PC routines
Table of Contents for all manual pages
Index of all manual pages