petsc-3.9.4 2018-09-11
PCASMCreateSubdomains
Creates the index sets for the overlapping Schwarz preconditioner for a any problem on a general grid.
Synopsis
#include "petscpc.h"
PetscErrorCode PCASMCreateSubdomains(Mat A, PetscInt n, IS* outis[])
Collective
Input Parameters
| A | - The global matrix operator
|
| n | - the number of local blocks
|
Output Parameters
outis -the array of index sets defining the subdomains
Note: this generates nonoverlapping subdomains; the PCASM will generate the overlap
from these if you use PCASMSetLocalSubdomains()
In the Fortran version you must provide the array outis[] already allocated of length n.
Keywords
PC, ASM, additive Schwarz, create, subdomains, unstructured grid
See Also
PCASMSetLocalSubdomains(), PCASMDestroySubdomains()
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