:orphan: # DMMoabPToRMapping Compute the mapping from the physical coordinate system for a given element to the canonical reference element. In addition to finding the inverse mapping evaluation through Newton iteration, the basis function at the parametric point is also evaluated optionally. ## Synopsis ``` #include "petscdt.h" #include "petscdmmoab.h" PetscErrorCode DMMoabPToRMapping(const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscReal *xphy, PetscReal *natparam, PetscReal *phi) ``` ## Input Parameters - ***PetscInt dim -*** the element dimension (1=EDGE, 2=QUAD/TRI, 3=HEX/TET) - ***PetscInt nverts -*** the number of vertices in the physical element - ***PetscReal coordinates -*** the coordinates of vertices in the physical element - ***PetscReal[3] xphy -*** the coordinates of physical point for which natural coordinates (in reference frame) are sought ## Output Parameters - ***PetscReal[3] natparam -*** the natural coordinates (in reference frame) corresponding to xphy - ***PetscReal[nverts] phi -*** the basis functions evaluated at the natural coordinates (natparam) ## Level advanced ## Location src/dm/impls/moab/dmmbfem.cxx --- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/dm/impls/moab/dmmbfem.cxx) [Index of all DMMOAB routines](index.md) [Table of Contents for all manual pages](/manualpages/index.md) [Index of all manual pages](/manualpages/singleindex.md)