petsc-3.13.6 2020-09-29
Report Typos and Errors

PCBDDCSetDivergenceMat

Sets the linear operator representing \int_\Omega \div {\bf u} \cdot p dx

Synopsis

#include "petscpc.h" 
PetscErrorCode PCBDDCSetDivergenceMat(PC pc, Mat divudotp, PetscBool trans, IS vl2l)
Collective on PC

Input Parameters

pc - the preconditioning context
divudotp - the matrix (must be of type MATIS)
trans - if trans if false (resp. true), then pressures are in the test (trial) space and velocities are in the trial (test) space.
vl2l - optional index set describing the local (wrt the local matrix in divudotp) to local (wrt the local matrix in the preconditioning matrix) map for the velocities

Notes

This auxiliary matrix is used to compute quadrature weights representing the net-flux across subdomain boundaries If vl2l is NULL, the local ordering for velocities in divudotp should match that of the preconditioning matrix

See Also

PCBDDC

Level

advanced

Location

src/ksp/pc/impls/bddc/bddc.c

Examples

src/ksp/ksp/tutorials/ex43.c.html

Implementations

PCBDDCSetDivergenceMat_BDDC in src/ksp/pc/impls/bddc/bddc.c

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