PCGASMCreateSubdomains#
Creates n
index sets defining n
nonoverlapping subdomains on this MPI process for the PCGASM
additive Schwarz preconditioner for a any problem based on its matrix.
Synopsis#
#include "petscpc.h"
PetscErrorCode PCGASMCreateSubdomains(Mat A, PetscInt N, PetscInt *n, IS *iis[])
Collective
Input Parameters#
A - The global matrix operator
N - the number of global subdomains requested
Output Parameters#
n - the number of subdomains created on this MPI rank
iis - the array of index sets defining the local inner subdomains (on which the correction is applied)
Notes#
When N
>= A’s communicator size, each subdomain is local – contained within a single MPI process.
When N
< size, the subdomains are ‘straddling’ (rank boundaries) and are no longer local.
The resulting subdomains can be use in PCGASMSetSubdomains
(pc,n,iss,NULL
). The overlapping
outer subdomains will be automatically generated from these according to the requested amount of
overlap; this is currently supported only with local subdomains.
Use PCGASMDestroySubdomains()
to free the array and the list of index sets.
See Also#
Level#
advanced
Location#
Index of all PC routines
Table of Contents for all manual pages
Index of all manual pages