petsc-3.9.4 2018-09-11
PCASMGetLocalSubmatrices
Gets the local submatrices (for this processor only) for the additive Schwarz preconditioner.
Synopsis
#include "petscpc.h"
PetscErrorCode PCASMGetLocalSubmatrices(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: Call after PCSetUp() (or KSPSetUp()) but before PCApply() (or KSPApply()) and before PCSetUpOnBlocks())
Usually one would use PCSetModifySubmatrices() to change the submatrices in building the preconditioner.
Keywords
PC, ASM, set, local, subdomains, additive Schwarz, block Jacobi
See Also
PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains(), PCSetModifySubmatrices()
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