petsc-3.4.5 2014-06-29

PCGASMGetSubdomains

Gets the subdomains supported on this processor for the additive Schwarz preconditioner.

Synopsis

#include "petscpc.h" 
PetscErrorCode  PCGASMGetSubdomains(PC pc,PetscInt *n,IS *iis[],IS *ois[])
Not Collective

Input Parameter

pc -the preconditioner context

Output Parameters

n - the number of subdomains for this processor (default value = 1)
iis - the index sets that define the inner subdomains (without overlap) supported on this processor (can be NULL)
ois - the index sets that define the outer subdomains (with overlap) supported on this processor (can be NULL)

Notes

The user is responsible for destroying the ISs and freeing the returned arrays. The IS numbering is in the parallel, global numbering of the vector.

Keywords

PC, GASM, get, subdomains, additive Schwarz

See Also

PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubmatrices()

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