#include "petscpc.h" PetscErrorCode PCASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[])Collective on pc
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) |
-pc_asm_local_blocks <blks> | - Sets local blocks |
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.