DMPlexComputeProjection3Dto1D#
Rewrite coordinates to be the 1D projection of the 3D coordinates
Synopsis#
#include "petscdmplex.h"
#include "petscfe.h"
PetscErrorCode DMPlexComputeProjection3Dto1D(PetscScalar coords[], PetscReal R[])
Not Collective
Input/Output Parameter#
coords - The coordinates of a segment; on output, the new y-coordinate, and 0 for x and z, an array of size 6, the other entries are unchanged
Output Parameter#
R - The rotation which accomplishes the projection, an array of size 9
Note#
This uses the basis completion described by Frisvad [Fri12]
References#
[Fri12]
Jeppe Revall Frisvad. Building an orthonormal basis from a 3D unit vector without normalization. Journal of Graphics Tools, 16(3):151–159, 2012.
See Also#
DMPLEX
, DMPlexComputeProjection2Dto1D()
, DMPlexComputeProjection3Dto2D()
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