petsc-3.11.4 2019-09-28
PCGASMGetSubmatrices
Gets the local submatrices (for this processor only) for the additive Schwarz preconditioner.
Synopsis
#include "petscpc.h"
PetscErrorCode PCGASMGetSubmatrices(PC pc,PetscInt *n,Mat *mat[])
Not Collective
Input Parameter
pc -the preconditioner context
Output Parameters
| n | - the number of matrices for this processor (default value = 1)
|
| mat | - the matrices
|
Notes
matrices returned by this routine have the same communicators as the index sets (IS)
used to define subdomains in PCGASMSetSubdomains()
Keywords
PC, GASM, set, local, subdomains, additive Schwarz, block Jacobi
See Also
PCGASMSetOverlap(), PCGASMGetSubKSP(),
PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubdomains()
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