DMPlexRemapGeometry#
This function maps the original DM coordinates to new coordinates.
Synopsis#
#include "petscdmplex.h"   
#include "petscfe.h"       
PetscErrorCode DMPlexRemapGeometry(DM dm, PetscReal time, void (*func)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[]))
Not Collective
Input Parameters#
- dm - The - DM
- time - The time 
- func - The function transforming current coordinates to new coordinates 
Calling sequence of func#
- dim - The spatial dimension 
- Nf - The number of input fields (here 1) 
- NfAux - The number of input auxiliary fields 
- uOff - The offset of the coordinates in u[] (here 0) 
- uOff_x - The offset of the coordinates in u_x[] (here 0) 
- u - The coordinate values at this point in space 
- u_t - The coordinate time derivative at this point in space (here - NULL)
- u_x - The coordinate derivatives at this point in space 
- aOff - The offset of each auxiliary field in u[] 
- aOff_x - The offset of each auxiliary field in u_x[] 
- a - The auxiliary field values at this point in space 
- a_t - The auxiliary field time derivative at this point in space (or - NULL)
- a_x - The auxiliary field derivatives at this point in space 
- t - The current time 
- x - The coordinates of this point (here not used) 
- numConstants - The number of constants 
- constants - The value of each constant 
- f - The new coordinates at this point in space 
See Also#
DMPLEX, DMGetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM(), DMProjectFieldLocal(), DMProjectFieldLabelLocal()
Level#
intermediate
Location#
src/dm/impls/plex/plexgeometry.c
Index of all DMPlex routines
Table of Contents for all manual pages
Index of all manual pages