PCBDDCCreateFETIDPOperators#

Create FETI-DP operators

Synopsis#

#include "petscpc.h" 
PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, PetscBool fully_redundant, const char *prefix, Mat *fetidp_mat, PC *fetidp_pc)

Collective

Input Parameters#

  • pc - the PCBDDC preconditioning context (setup should have been called before)

  • fully_redundant - true for a fully redundant set of Lagrange multipliers

  • prefix - optional options database prefix for the objects to be created (can be NULL)

Output Parameters#

  • fetidp_mat - shell FETI-DP matrix object

  • fetidp_pc - shell Dirichlet preconditioner for FETI-DP matrix

Note#

Currently the only operations provided for FETI-DP matrix are MatMult() and MatMultTranspose()

See Also#

PCBDDC, PCBDDCMatFETIDPGetRHS(), PCBDDCMatFETIDPGetSolution()

Level#

developer

Location#

src/ksp/pc/impls/bddc/bddc.c

Examples#

src/ksp/ksp/tutorials/ex59.c

Implementations#

PCBDDCCreateFETIDPOperators_BDDC in src/ksp/pc/impls/bddc/bddc.c


Edit on GitLab

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