DMSnapToGeomModel#
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 "petscdm.h"
PetscErrorCode DMSnapToGeomModel(DM dm, PetscInt p, PetscInt dE, const PetscScalar mcoords[], PetscScalar gcoords[])
Not Collective
Input Parameters#
dm - The
DMPLEX
objectp - 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.
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#
Implementations#
DMSnapToGeomModel_EGADSLite() in src/dm/impls/plex/plexegads.c
DMSnapToGeomModel_EGADS() in src/dm/impls/plex/plexegads.c
Index of all DM routines
Table of Contents for all manual pages
Index of all manual pages