petsc-3.11.4 2019-09-28
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
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
DMSwarmSetCellDM(), DMSwarmInsertPointsUsingCellDM()
Level
beginner
Location
src/dm/impls/swarm/swarmpic.c
Index of all DMSWARM routines
Table of Contents for all manual pages
Index of all manual pages