petsc-3.13.6 2020-09-29
MatSchurComplementGetAinvType
get the type of approximation for the inverse of the (0,0) block used in forming Sp in MatSchurComplementGetPmat()
Synopsis
#include "petscksp.h"
PetscErrorCode MatSchurComplementGetAinvType(Mat S,MatSchurComplementAinvType *ainvtype)
Not collective.
Input Parameter
S -matrix obtained with MatCreateSchurComplement() (or equivalent) and implementing the action of A11 - A10 ksp(A00,Ap00) A01
Output Parameter
ainvtype -type of approximation used to form A00inv from A00 when assembling Sp = A11 - A10 A00inv A01:
MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, or MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG
Note
Since the real Schur complement is usually dense, providing a good approximation to newpmat usually requires
Section 1.5 Writing Application Codes with PETSc-specific information. The default for assembled matrices is to use the inverse of the diagonal of
the (0,0) block A00 in place of A00^{-1}. This rarely produce a scalable algorithm. Optionally, A00 can be lumped
before forming inv(diag(A00)).
See Also
MatSchurComplementAinvType, MatCreateSchurComplement(), MatGetSchurComplement(), MatSchurComplementGetPmat(), MatSchurComplementSetAinvType()
Level
advanced
Location
src/ksp/ksp/utils/schurm/schurm.c
Index of all KSP routines
Table of Contents for all manual pages
Index of all manual pages