PETSc version 3.15.5
Fix/Edit manual page

PCGASMCreateSubdomains

Creates n index sets defining n nonoverlapping subdomains for the 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 processor
iis - the array of index sets defining the local inner subdomains (on which the correction is applied)

Note: When N >= A's communicator size, each subdomain is local -- contained within a single processor. When N < size, the subdomains are 'straddling' (processor 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.

See Also

PCGASMSetSubdomains(), PCGASMDestroySubdomains()

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