petsc-3.6.1 2015-08-06
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 (or NULL to use the default of 1 subdomain per process)

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.

Keywords

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

See Also

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

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/ex8.c.html