petsc-3.13.6 2020-09-29
DMPlexGetSupport
Return the points on the out-edges for this point in the DAG
Synopsis
#include "petscdmplex.h"
PetscErrorCode DMPlexGetSupport(DM dm, PetscInt p, const PetscInt *support[])
Not collective
Input Parameters
| mesh | - The DMPlex
|
| p | - The point, which must lie in the chart set with DMPlexSetChart()
|
Output Parameter
support -An array of points which are on the out-edges for point p
Fortran Notes
Since it returns an array, this routine is only available in Fortran 90, and you must
include petsc.h90 in your code.
You must also call DMPlexRestoreSupport() after you finish using the returned array.
DMPlexRestoreSupport() is not needed/available in C.
See Also
DMPlexCreate(), DMPlexSetCone(), DMPlexSetChart(), DMPlexGetCone()
Level
beginner
Location
src/dm/impls/plex/plex.c
Examples
src/ts/tutorials/ex11.c.html
Index of all DMPLEX routines
Table of Contents for all manual pages
Index of all manual pages