petsc-3.9.4 2018-09-11
PCGASMGetSubKSP
Gets the local KSP contexts for all blocks on this processor.
Synopsis
#include "petscpc.h"
PetscErrorCode PCGASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
Collective on PC iff first_local is requested
Input Parameter
pc -the preconditioner context
Output Parameters
| n_local | - the number of blocks on this processor or NULL
|
| first_local | - the global number of the first block on this processor or NULL,
all processors 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 processor
is supported.
You must call KSPSetUp() before calling PCGASMGetSubKSP().
Keywords
PC, GASM, additive Schwarz, get, sub, KSP, context
See Also
PCGASMSetSubdomains(), PCGASMSetOverlap(),
PCGASMCreateSubdomains2D(),
Level
advanced
Location
src/ksp/pc/impls/gasm/gasm.c
Examples
src/ksp/ksp/examples/tutorials/ex62.c.html
Implementations
PCGASMGetSubKSP_GASM in src/ksp/pc/impls/gasm/gasm.c
Index of all PC routines
Table of Contents for all manual pages
Index of all manual pages