DMPlexComputeProjection3Dto2D#

Rewrite coordinates of 3 or more coplanar 3D points to a common 2D basis for the plane. The normal is defined by positive orientation of the first 3 points.

Synopsis#

#include "petscdmplex.h"   
#include "petscfe.h"       
PetscErrorCode DMPlexComputeProjection3Dto2D(PetscInt coordSize, PetscScalar coords[], PetscReal R[])

Not Collective

Input Parameter#

  • coordSize - Length of coordinate array (3x number of points); must be at least 9 (3 points)

Input/Output Parameter#

  • coords - The interlaced coordinates of each coplanar 3D point; on output the first 2*coordSize/3 entries contain interlaced 2D points, with the rest undefined

Output Parameter#

  • R - 3x3 row-major rotation matrix whose columns are the tangent basis [t1, t2, n]. Multiplying by R^T transforms from original frame to tangent frame.

See Also#

DMPLEX, DMPlexComputeProjection2Dto1D(), DMPlexComputeProjection3Dto1D()

Level#

developer

Location#

src/dm/impls/plex/plexgeometry.c


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