petsc-3.8.4 2018-03-24
Report Typos and Errors

PCASMSetLocalSubdomains

Sets the local subdomains (for this processor only) for the additive Schwarz preconditioner.

Synopsis

#include "petscpc.h" 
PetscErrorCode  PCASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[])
Collective on PC

Input Parameters

pc - the preconditioner context
n - the number of subdomains for this processor (default value = 1)
is - the index set that defines the subdomains for this processor (or NULL for PETSc to determine subdomains)
is_local - the index sets that define the local part of the subdomains for this processor, not used unless PCASMType is PC_ASM_RESTRICT (or NULL to not provide these)

Notes

The IS numbering is in the parallel, global numbering of the vector for both is and is_local

By default the ASM preconditioner uses 1 block per processor.

Use PCASMSetTotalSubdomains() to set the subdomains for all processors.

If is_local is provided and PCASMType is PC_ASM_RESTRICT then the solution only over the is_local region is interpolated back to form the global solution (this is the standard restricted additive Schwarz method)

If the is_local is provided and PCASMType is PC_ASM_INTERPOLATE or PC_ASM_NONE then an error is generated since there is no code to handle that case.

Keywords

PC, ASM, set, local, subdomains, additive Schwarz

See Also

PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), PCASMType, PCASMSetType()

Level:advanced
Location:
src/ksp/pc/impls/asm/asm.c
Index of all PC routines
Table of Contents for all manual pages
Index of all manual pages

Examples

src/ksp/ksp/examples/tutorials/ex6.c.html
src/ksp/ksp/examples/tutorials/ex8.c.html