PCASMGetLocalSubmatrices#
Gets the local submatrices (for this processor only) for the additive Schwarz preconditioner, PCASM
.
Synopsis#
#include "petscpc.h"
PetscErrorCode PCASMGetLocalSubmatrices(PC pc, PetscInt *n, Mat *mat[])
Not Collective
Input Parameter#
pc - the preconditioner context
Output Parameters#
n - if requested, the number of matrices for this processor (default value = 1)
mat - if requested, the matrices
Notes#
Call after PCSetUp()
(or KSPSetUp()
) but before PCApply()
and before PCSetUpOnBlocks()
)
Usually one would use PCSetModifySubMatrices()
to change the submatrices in building the preconditioner.
See Also#
KSP: Linear System Solvers, PCASM
, PCASMSetTotalSubdomains()
, PCASMSetOverlap()
, PCASMGetSubKSP()
,
PCASMCreateSubdomains2D()
, PCASMSetLocalSubdomains()
, PCASMGetLocalSubdomains()
, PCSetModifySubMatrices()
Level#
advanced
Location#
Index of all PC routines
Table of Contents for all manual pages
Index of all manual pages