PCGASMGetSubdomains#
Gets the subdomains supported on this MPI process for the PCGASM
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 MPI process (default value = 1)
iis - the index sets that define the inner subdomains (without overlap) supported on this process (can be
NULL
)ois - the index sets that define the outer subdomains (with overlap) supported on this process (can be
NULL
)
Notes#
The user is responsible for destroying the IS
s and freeing the returned arrays, this can be done with
PCGASMDestroySubdomains()
The IS
numbering is in the parallel, global numbering of the vector.
See Also#
KSP: Linear System Solvers, PCGASM
, PCGASMSetOverlap()
, PCGASMGetSubKSP()
, PCGASMCreateSubdomains2D()
,
PCGASMSetSubdomains()
, PCGASMGetSubmatrices()
, PCGASMDestroySubdomains()
Level#
advanced
Location#
Index of all PC routines
Table of Contents for all manual pages
Index of all manual pages