PCASMSetOverlap#
Sets the overlap between a pair of subdomains for the additive Schwarz preconditioner, PCASM
.
Synopsis#
#include "petscpc.h"
PetscErrorCode PCASMSetOverlap(PC pc, PetscInt ovl)
Logically Collective
Input Parameters#
pc - the preconditioner context
ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1)
Options Database Key#
-pc_asm_overlap
- Sets overlap
Notes#
By default the PCASM
preconditioner uses 1 block per processor. To use
multiple blocks per perocessor, see PCASMSetTotalSubdomains()
and
PCASMSetLocalSubdomains()
(and the option -pc_asm_blocks
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 PCASM
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.
One can define initial index sets with any overlap via
PCASMSetLocalSubdomains()
; the routine
PCASMSetOverlap()
merely allows PETSc to extend that overlap further
if desired.
See Also#
PCASM
, PCASMSetTotalSubdomains()
, PCASMSetLocalSubdomains()
, PCASMGetSubKSP()
,
PCASMCreateSubdomains2D()
, PCASMGetLocalSubdomains()
, MatIncreaseOverlap()
, PCGASM
Level#
intermediate
Location#
Examples#
Implementations#
PCASMSetOverlap_ASM in src/ksp/pc/impls/asm/asm.c
Index of all PC routines
Table of Contents for all manual pages
Index of all manual pages