#include "petscpc.h" PetscErrorCode PCASMSetOverlap(PC pc,PetscInt ovl)Logically Collective on PC
pc | - the preconditioner context | |
ovl | - the amount of overlap between subdomains (ovl >= 0, default value = 1) |
The overlap defaults to 1, so if one desires that no additional overlap be computed beyond what may have been set with a call to PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(), then ovl must be set to be 0. In particular, if one does not explicitly set the subdomains an application code, then all overlap would be computed internally by PETSc, and using an overlap of 0 would result in an ASM variant that is equivalent to the block Jacobi preconditioner.
The default algorithm used by PETSc to increase overlap is fast, but not scalable, use the option -mat_increase_overlap_scalable when the problem and number of processes is large.
Note that one can define initial index sets with any overlap via PCASMSetLocalSubdomains(); the routine PCASMSetOverlap() merely allows PETSc to extend that overlap further if desired.