PCASMSetLocalSubdomains#
Sets the local subdomains (for this processor only) for the additive Schwarz preconditioner PCASM
.
Synopsis#
#include "petscpc.h"
PetscErrorCode PCASMSetLocalSubdomains(PC pc, PetscInt n, IS is[], IS is_local[])
Collective
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) the values of theis
array are copied so you can free the array (not theIS
in the array) after this callis_local - the index sets that define the local part of the subdomains for this processor, not used unless
PCASMType
isPC_ASM_RESTRICT
(orNULL
to not provide these). The values of theis_local
array are copied so you can free the array (not theIS
in the array) after this call
Options Database Key#
-pc_asm_local_blocks
- Sets number of local blocks
Notes#
The IS
numbering is in the parallel, global numbering of the vector for both is
and is_local
By default the PCASM
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, RASM)
If 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.
See Also#
KSP: Linear System Solvers, PCASM
, PCASMSetTotalSubdomains()
, PCASMSetOverlap()
, PCASMGetSubKSP()
,
PCASMCreateSubdomains2D()
, PCASMGetLocalSubdomains()
, PCASMType
, PCASMSetType()
, PCGASM
Level#
advanced
Location#
Examples#
src/ksp/ksp/tutorials/ex76.c
src/ksp/ksp/tutorials/ex8.c
src/ksp/ksp/tutorials/ex19.c
Implementations#
PCASMSetLocalSubdomains_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