:orphan: # DMSwarmSetPointCoordinatesCellwise Insert point coordinates (defined over the reference cell) within each cell ## Synopsis ``` #include "petscdmswarm.h" PETSC_EXTERN PetscErrorCode DMSwarmSetPointCoordinatesCellwise(DM dm, PetscInt npoints, PetscReal xi[]) ``` Not Collective ## Input parameters - ***dm -*** the `DMSWARM` - ***celldm -*** the cell `DM` - ***npoints -*** the number of points to insert in each cell - ***xi -*** the coordinates (defined in the local coordinate system for each cell) to insert ## Notes The method will reset any previous defined points within the `DMSWARM`. Only supported for `DMPLEX`. If you are using a `DMDA` it is recommended to either use `DMSwarmInsertPointsUsingCellDM()`, or extract and set the coordinates yourself the following code ```none PetscReal *coor; DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&coor); // user code to define the coordinates here DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&coor); ``` ## See Also `DMSWARM`, `DMSwarmSetCellDM()`, `DMSwarmInsertPointsUsingCellDM()` ## Level beginner ## Location src/dm/impls/swarm/swarmpic.c --- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/dm/impls/swarm/swarmpic.c) [Index of all DMSwarm routines](index.md) [Table of Contents for all manual pages](/manualpages/index.md) [Index of all manual pages](/manualpages/singleindex.md)