petsc-3.14.6 2021-03-30
PCHPDDMHasNeumannMat
Informs PCHPDDM that the Mat passed to PCHPDDMSetAuxiliaryMat() is the local Neumann matrix. This may be used to bypass a call to MatCreateSubMatrices() and to MatConvert() for MATMPISBAIJ matrices. If a DMCreateNeumannOverlap() implementation is available in the DM attached to the Pmat, or the Amat, or the PC, the flag is internally set to PETSC_TRUE. Its default value is otherwise PETSC_FALSE.
Synopsis
#include "petscpc.h"
PetscErrorCode PCHPDDMHasNeumannMat(PC pc, PetscBool has)
Input Parameters
| pc | - preconditioner context
|
| has | - Boolean value
|
See Also
PCHPDDM, PCHPDDMSetAuxiliaryMat()
Level
intermediate
Location
src/ksp/pc/impls/hpddm/hpddm.cxx
Examples
src/ksp/ksp/tutorials/ex76.c.html
src/ksp/ksp/tutorials/ex76f.F90.html
Implementations
PCHPDDMHasNeumannMat_HPDDM in src/ksp/pc/impls/hpddm/hpddm.cxx
Index of all PC routines
Table of Contents for all manual pages
Index of all manual pages