PCHPDDMSetAuxiliaryMat#

Sets the auxiliary matrix used by PCHPDDM for the concurrent GenEO problems at the finest level. As an example, in a finite element context with nonoverlapping subdomains plus (overlapping) ghost elements, this could be the unassembled (Neumann) local overlapping operator. As opposed to the assembled (Dirichlet) local overlapping operator obtained by summing neighborhood contributions at the interface of ghost elements.

Synopsis#

#include "petscpc.h" 
PetscErrorCode PCHPDDMSetAuxiliaryMat(PC pc, IS is, Mat A, PetscErrorCode (*setup)(Mat, PetscReal, Vec, Vec, PetscReal, IS, void *), void *setup_ctx)

Input Parameters#

  • pc - preconditioner context

  • is - index set of the local auxiliary, e.g., Neumann, matrix

  • A - auxiliary sequential matrix

  • setup - function for generating the auxiliary matrix

  • setup_ctx - context for setup

See Also#

PCHPDDM, PCCreate(), PCSetType(), PCType, PC, PCHPDDMSetRHSMat(), MATIS

Level#

intermediate

Location#

src/ksp/pc/impls/hpddm/pchpddm.cxx

Examples#

src/ksp/ksp/tutorials/ex76.c
src/ksp/ksp/tutorials/ex76f.F90
src/ksp/ksp/tutorials/ex82.c

Implementations#

PCHPDDMSetAuxiliaryMat_HPDDM(PC pc, IS is, Mat A, PetscErrorCode (*setup) in src/ksp/pc/impls/hpddm/pchpddm.cxx


Edit on GitLab

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