DMPlexSnapToGeomModel#

Given a coordinate point ‘mcoords’ on the mesh point ‘p’, return the closest coordinate point ‘gcoords’ on the geometry model associated with that point.

Synopsis#

#include "petscdmplex.h"   
PetscErrorCode DMPlexSnapToGeomModel(DM dm, PetscInt p, PetscInt dE, const PetscScalar mcoords[], PetscScalar gcoords[])

Not Collective

Input Parameters#

  • dm - The DMPLEX object

  • p - The mesh point

  • dE - The coordinate dimension

  • mcoords - A coordinate point lying on the mesh point

Output Parameter#

  • gcoords - The closest coordinate point on the geometry model associated with ‘p’ to the given point

Note#

Returns the original coordinates if no geometry model is found. Right now the only supported geometry model is EGADS.

The coordinate dimension may be different from the coordinate dimension of the dm, for example if the transformation is extrusion.

See Also#

DMPlex: Unstructured Grids, DM, DMPLEX, DMRefine(), DMPlexCreate(), DMPlexSetRefinementUniform()

Level#

intermediate

Location#

src/dm/impls/plex/plexegads.c


Edit on GitLab

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