MatSchurComplementGetPmat#
Obtain a preconditioning matrix for the Schur complement by assembling \(Sp = A11 - A10 inv(DIAGFORM(A00)) A01\)
Synopsis#
#include "petscksp.h"
PetscErrorCode MatSchurComplementGetPmat(Mat S, MatReuse preuse, Mat *Sp)
Collective
Input Parameters#
S - matrix obtained with MatCreateSchurComplement() (or equivalent) that implements the action of \(A11 - A10 ksp(A00,Ap00) A01\)
preuse -
MAT_INITIAL_MATRIXfor a newSp, orMAT_REUSE_MATRIXto reuse an existingSp, orMAT_IGNORE_MATRIXto put nothing inSp
Output Parameter#
Sp - approximate Schur complement suitable for preconditioning the exact Schur complement \(S = A11 - A10 inv(A00) A01\)
Notes#
The approximation of Sp depends on the argument passed to MatSchurComplementSetAinvType()
MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG, or MAT_SCHUR_COMPLEMENT_AINV_FULL
-mat_schur_complement_ainv_type <diag,lump,blockdiag,full>
Sometimes users would like to provide problem-specific data in the Schur complement, usually only
for special row and column index sets. In that case, the user should call PetscObjectComposeFunction() to set
“MatSchurComplementGetPmat_C” to their function. If their function needs to fall back to the default implementation,
it should call MatSchurComplementGetPmat_Basic().
Developer Notes#
The API that includes MatGetSchurComplement(), MatCreateSchurComplement(), MatSchurComplementGetPmat() should be refactored to
remove redundancy and be clearer and simpler.
This routine should be called MatSchurComplementCreatePmat()
See Also#
KSP: Linear System Solvers, MatCreateSubMatrix(), PCFIELDSPLIT, MatGetSchurComplement(), MatCreateSchurComplement(), 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