PCFieldSplitSchurGetS#

extract the MATSCHURCOMPLEMENT object used by this PCFIELDSPLIT in case it needs to be configured separately

Synopsis#

#include "petscpc.h" 
PetscErrorCode PCFieldSplitSchurGetS(PC pc, Mat *S)

Not Collective

Input Parameter#

  • pc - the preconditioner context

Output Parameter#

  • S - the Schur complement matrix

Note#

This matrix should not be destroyed using MatDestroy(); rather, use PCFieldSplitSchurRestoreS().

See Also#

Solving Block Matrices with PCFIELDSPLIT, PC, PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSchurPreType, PCFieldSplitSetSchurPre(), MATSCHURCOMPLEMENT, PCFieldSplitSchurRestoreS(), MatCreateSchurComplement(), MatSchurComplementGetKSP(), MatSchurComplementComputeExplicitOperator(), MatGetSchurComplement()

Level#

advanced

Location#

src/ksp/pc/impls/fieldsplit/fieldsplit.c


Index of all PC routines
Table of Contents for all manual pages
Index of all manual pages