PCHPDDMSetAuxiliaryMat#

Sets the auxiliary matrix used by PCHPDDM for the concurrent GenEO problems at the finest level.

Synopsis#

#include "petscpc.h" 
PetscErrorCode PCHPDDMSetAuxiliaryMat(PC pc, IS is, Mat A, PetscErrorCode (*setup)(Mat J, PetscReal t, Vec X, Vec X_t, PetscReal s, IS ovl, void *ctx), void *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 entries, may be NULL

  • ctx - context for setup, may be NULL

Calling sequence of setup#

  • J - matrix whose values are to be set

  • t - time

  • X - linearization point

  • X_t - time-derivative of the linearization point

  • s - step

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

  • ctx - context for setup, may be NULL

Note#

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.

Fortran Notes#

Only PETSC_NULL_FUNCTION is supported for setup and ctx is never accessed

See Also#

KSP: Linear System Solvers, PCHPDDM, PCCreate(), PCSetType(), PCType, PC, PCHPDDMSetRHSMat(), MATIS

Level#

intermediate

Location#

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

Examples#

src/ksp/ksp/tutorials/ex82.c
src/ksp/ksp/tutorials/ex87.c
src/ksp/ksp/tutorials/ex76.c
src/ksp/ksp/tutorials/ex27.c

Implementations#

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


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