petsc-3.12.5 2020-03-29
Report Typos and Errors

DMPlexFindVertices

Try to find DAG points based on their coordinates.

Synopsis

#include "petscdmplex.h"   
#include "petscfe.h"       
PetscErrorCode DMPlexFindVertices(DM dm, PetscInt npoints, const PetscReal coord[], PetscReal eps, PetscInt dagPoints[])
Not Collective (provided DMGetCoordinatesLocalSetUp() has been called already)

Input Parameters

dm - The DMPlex object
npoints - The number of sought points
coords - The array of coordinates of the sought points
eps - The tolerance or PETSC_DEFAULT

Output Parameters

dagPoints -The array of found DAG points, or -1 if not found

Notes

The length of the array coords must be npoints * dim where dim is the spatial dimension returned by DMGetDimension().

The output array dagPoints is NOT newly allocated; the user must pass an array of length npoints.

Each rank does the search independently; a nonnegative value is returned only if this rank's local DMPlex portion contains the point.

The tolerance is interpreted as the maximum Euclidean (L2) distance of the sought point from the specified coordinates.

See Also

DMPlexCreate(), DMGetCoordinates()

Level

intermediate

Location

src/dm/impls/plex/plexgeometry.c
Index of all DMPLEX routines
Table of Contents for all manual pages
Index of all manual pages