PCGASMGetSubKSP#

Gets the local KSP contexts for all subdomains on this MPI rank.

Synopsis#

#include "petscpc.h" 
PetscErrorCode PCGASMGetSubKSP(PC pc, PetscInt *n_local, PetscInt *first_local, KSP *ksp[])

Collective iff first_local is requested

Input Parameter#

  • pc - the preconditioner context

Output Parameters#

  • n_local - the number of blocks on this MPI rank or NULL

  • first_local - the global number of the first block on this rank or NULL, all ranks must request or all must pass NULL

  • ksp - the array of KSP contexts

Note#

After PCGASMGetSubKSP() the array of KSPes is not to be freed

Currently for some matrix implementations only 1 block per MPI process is supported.

You must call KSPSetUp() before calling PCGASMGetSubKSP().

See Also#

PCGASM, PCGASMSetSubdomains(), PCGASMSetOverlap(), PCGASMCreateSubdomains2D(),

Level#

advanced

Location#

src/ksp/pc/impls/gasm/gasm.c

Examples#

src/ksp/ksp/tutorials/ex62.c

Implementations#

PCGASMGetSubKSP_GASM in src/ksp/pc/impls/gasm/gasm.c


Edit on GitLab

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