# DMPlexCoordinatesToReference
Pull coordinates back from the mesh to the reference element using a single element map.  This inversion will be accurate inside the reference element, but may be inaccurate for mappings that do not extend uniquely outside the reference cell (e.g, most non-affine maps) 
## Synopsis
```
#include "petscdmplex.h"   
#include "petscfe.h"       
PetscErrorCode DMPlexCoordinatesToReference(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal realCoords[], PetscReal refCoords[])
```
Not collective


## Input Parameters

- ***dm         -*** The mesh, with coordinate maps defined either by a PetscDS for the coordinate DM (see DMGetCoordinateDM()) or
implicitly by the coordinates of the corner vertices of the cell: as an affine map for simplicial elements, or
as a multilinear map for tensor-product elements
- ***cell       -*** the cell whose map is used.
- ***numPoints  -*** the number of points to locate
- ***realCoords -*** (numPoints x coordinate dimension) array of coordinates (see DMGetCoordinateDim())



## Output Parameters

- ***refCoords  -*** (numPoints x dimension) array of reference coordinates (see DMGetDimension())





## See Also
 `DMPlexReferenceToCoordinates()`

## Level
intermediate

## Location
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/dm/impls/plex/plexgeometry.c.html#DMPlexCoordinatesToReference">src/dm/impls/plex/plexgeometry.c</A>

## Implementations

<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/dm/impls/plex/plexgeometry.c.html#DMPlexCoordinatesToReference_NewtonUpdate">DMPlexCoordinatesToReference_NewtonUpdate in src/dm/impls/plex/plexgeometry.c</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/dm/impls/plex/plexgeometry.c.html#DMPlexCoordinatesToReference_Tensor">DMPlexCoordinatesToReference_Tensor in src/dm/impls/plex/plexgeometry.c</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/dm/impls/plex/plexgeometry.c.html#DMPlexCoordinatesToReference_FE">DMPlexCoordinatesToReference_FE in src/dm/impls/plex/plexgeometry.c</A><BR>


---
[Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/dm/impls/plex/plexgeometry.c)


[Index of all DMPlex routines](index.md)  
[Table of Contents for all manual pages](/docs/manualpages/index.md)  
[Index of all manual pages](/docs/manualpages/singleindex.md)  
