xref: /petsc/src/dm/impls/swarm/swarmpic.c (revision b75c6efc21bfcba5897c8ca359bc3d0e82c122c1)
1 #define PETSCDM_DLL
2 #include <petsc/private/dmswarmimpl.h> /*I   "petscdmswarm.h"   I*/
3 #include <petscsf.h>
4 #include <petscdmda.h>
5 #include <petscdmplex.h>
6 #include <petscdt.h>
7 #include "../src/dm/impls/swarm/data_bucket.h"
8 
9 #include <petsc/private/petscfeimpl.h> /* For CoordinatesRefToReal() */
10 
11 /*
12  Error checking to ensure the swarm type is correct and that a cell DM has been set
13 */
14 #define DMSWARMPICVALID(dm) \
15   do { \
16     DM_Swarm *_swarm = (DM_Swarm *)(dm)->data; \
17     PetscCheck(_swarm->swarm_type == DMSWARM_PIC, PetscObjectComm((PetscObject)(dm)), PETSC_ERR_SUP, "Valid only for DMSwarm-PIC. You must call DMSwarmSetType(dm,DMSWARM_PIC)"); \
18     PetscCheck(_swarm->dmcell, PetscObjectComm((PetscObject)(dm)), PETSC_ERR_SUP, "Valid only for DMSwarmPIC if the cell DM is set. You must call DMSwarmSetCellDM(dm,celldm)"); \
19   } while (0)
20 
21 /* Coordinate insertition/addition API */
22 /*@C
23    DMSwarmSetPointsUniformCoordinates - Set point coordinates in a DMSwarm on a regular (ijk) grid
24 
25    Collective on dm
26 
27    Input parameters:
28 +  dm - the DMSwarm
29 .  min - minimum coordinate values in the x, y, z directions (array of length dim)
30 .  max - maximum coordinate values in the x, y, z directions (array of length dim)
31 .  npoints - number of points in each spatial direction (array of length dim)
32 -  mode - indicates whether to append points to the swarm (ADD_VALUES), or over-ride existing points (INSERT_VALUES)
33 
34    Level: beginner
35 
36    Notes:
37    When using mode = INSERT_VALUES, this method will reset the number of particles in the DMSwarm
38    to be npoints[0]*npoints[1] (2D) or npoints[0]*npoints[1]*npoints[2] (3D). When using mode = ADD_VALUES,
39    new points will be appended to any already existing in the DMSwarm
40 
41 .seealso: `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType`
42 @*/
43 PETSC_EXTERN PetscErrorCode DMSwarmSetPointsUniformCoordinates(DM dm, PetscReal min[], PetscReal max[], PetscInt npoints[], InsertMode mode)
44 {
45   PetscReal          gmin[] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL};
46   PetscReal          gmax[] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL};
47   PetscInt           i, j, k, N, bs, b, n_estimate, n_curr, n_new_est, p, n_found;
48   Vec                coorlocal;
49   const PetscScalar *_coor;
50   DM                 celldm;
51   PetscReal          dx[3];
52   PetscInt           _npoints[] = {0, 0, 1};
53   Vec                pos;
54   PetscScalar       *_pos;
55   PetscReal         *swarm_coor;
56   PetscInt          *swarm_cellid;
57   PetscSF            sfcell = NULL;
58   const PetscSFNode *LA_sfcell;
59 
60   PetscFunctionBegin;
61   DMSWARMPICVALID(dm);
62   PetscCall(DMSwarmGetCellDM(dm, &celldm));
63   PetscCall(DMGetCoordinatesLocal(celldm, &coorlocal));
64   PetscCall(VecGetSize(coorlocal, &N));
65   PetscCall(VecGetBlockSize(coorlocal, &bs));
66   N = N / bs;
67   PetscCall(VecGetArrayRead(coorlocal, &_coor));
68   for (i = 0; i < N; i++) {
69     for (b = 0; b < bs; b++) {
70       gmin[b] = PetscMin(gmin[b], PetscRealPart(_coor[bs * i + b]));
71       gmax[b] = PetscMax(gmax[b], PetscRealPart(_coor[bs * i + b]));
72     }
73   }
74   PetscCall(VecRestoreArrayRead(coorlocal, &_coor));
75 
76   for (b = 0; b < bs; b++) {
77     if (npoints[b] > 1) {
78       dx[b] = (max[b] - min[b]) / ((PetscReal)(npoints[b] - 1));
79     } else {
80       dx[b] = 0.0;
81     }
82     _npoints[b] = npoints[b];
83   }
84 
85   /* determine number of points living in the bounding box */
86   n_estimate = 0;
87   for (k = 0; k < _npoints[2]; k++) {
88     for (j = 0; j < _npoints[1]; j++) {
89       for (i = 0; i < _npoints[0]; i++) {
90         PetscReal xp[] = {0.0, 0.0, 0.0};
91         PetscInt  ijk[3];
92         PetscBool point_inside = PETSC_TRUE;
93 
94         ijk[0] = i;
95         ijk[1] = j;
96         ijk[2] = k;
97         for (b = 0; b < bs; b++) xp[b] = min[b] + ijk[b] * dx[b];
98         for (b = 0; b < bs; b++) {
99           if (xp[b] < gmin[b]) point_inside = PETSC_FALSE;
100           if (xp[b] > gmax[b]) point_inside = PETSC_FALSE;
101         }
102         if (point_inside) n_estimate++;
103       }
104     }
105   }
106 
107   /* create candidate list */
108   PetscCall(VecCreate(PETSC_COMM_SELF, &pos));
109   PetscCall(VecSetSizes(pos, bs * n_estimate, PETSC_DECIDE));
110   PetscCall(VecSetBlockSize(pos, bs));
111   PetscCall(VecSetFromOptions(pos));
112   PetscCall(VecGetArray(pos, &_pos));
113 
114   n_estimate = 0;
115   for (k = 0; k < _npoints[2]; k++) {
116     for (j = 0; j < _npoints[1]; j++) {
117       for (i = 0; i < _npoints[0]; i++) {
118         PetscReal xp[] = {0.0, 0.0, 0.0};
119         PetscInt  ijk[3];
120         PetscBool point_inside = PETSC_TRUE;
121 
122         ijk[0] = i;
123         ijk[1] = j;
124         ijk[2] = k;
125         for (b = 0; b < bs; b++) xp[b] = min[b] + ijk[b] * dx[b];
126         for (b = 0; b < bs; b++) {
127           if (xp[b] < gmin[b]) point_inside = PETSC_FALSE;
128           if (xp[b] > gmax[b]) point_inside = PETSC_FALSE;
129         }
130         if (point_inside) {
131           for (b = 0; b < bs; b++) _pos[bs * n_estimate + b] = xp[b];
132           n_estimate++;
133         }
134       }
135     }
136   }
137   PetscCall(VecRestoreArray(pos, &_pos));
138 
139   /* locate points */
140   PetscCall(DMLocatePoints(celldm, pos, DM_POINTLOCATION_NONE, &sfcell));
141   PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell));
142   n_found = 0;
143   for (p = 0; p < n_estimate; p++) {
144     if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) n_found++;
145   }
146 
147   /* adjust size */
148   if (mode == ADD_VALUES) {
149     PetscCall(DMSwarmGetLocalSize(dm, &n_curr));
150     n_new_est = n_curr + n_found;
151     PetscCall(DMSwarmSetLocalSizes(dm, n_new_est, -1));
152   }
153   if (mode == INSERT_VALUES) {
154     n_curr    = 0;
155     n_new_est = n_found;
156     PetscCall(DMSwarmSetLocalSizes(dm, n_new_est, -1));
157   }
158 
159   /* initialize new coords, cell owners, pid */
160   PetscCall(VecGetArrayRead(pos, &_coor));
161   PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor));
162   PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid));
163   n_found = 0;
164   for (p = 0; p < n_estimate; p++) {
165     if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) {
166       for (b = 0; b < bs; b++) swarm_coor[bs * (n_curr + n_found) + b] = PetscRealPart(_coor[bs * p + b]);
167       swarm_cellid[n_curr + n_found] = LA_sfcell[p].index;
168       n_found++;
169     }
170   }
171   PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid));
172   PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor));
173   PetscCall(VecRestoreArrayRead(pos, &_coor));
174 
175   PetscCall(PetscSFDestroy(&sfcell));
176   PetscCall(VecDestroy(&pos));
177   PetscFunctionReturn(PETSC_SUCCESS);
178 }
179 
180 /*@C
181    DMSwarmSetPointCoordinates - Set point coordinates in a DMSwarm from a user defined list
182 
183    Collective on dm
184 
185    Input parameters:
186 +  dm - the DMSwarm
187 .  npoints - the number of points to insert
188 .  coor - the coordinate values
189 .  redundant - if set to PETSC_TRUE, it is assumed that npoints and coor[] are only valid on rank 0 and should be broadcast to other ranks
190 -  mode - indicates whether to append points to the swarm (ADD_VALUES), or over-ride existing points (INSERT_VALUES)
191 
192    Level: beginner
193 
194    Notes:
195    If the user has specified redundant = PETSC_FALSE, the cell DM will attempt to locate the coordinates provided by coor[] within
196    its sub-domain. If they any values within coor[] are not located in the sub-domain, they will be ignored and will not get
197    added to the DMSwarm.
198 
199 .seealso: `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType`, `DMSwarmSetPointsUniformCoordinates()`
200 @*/
201 PETSC_EXTERN PetscErrorCode DMSwarmSetPointCoordinates(DM dm, PetscInt npoints, PetscReal coor[], PetscBool redundant, InsertMode mode)
202 {
203   PetscReal          gmin[] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL};
204   PetscReal          gmax[] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL};
205   PetscInt           i, N, bs, b, n_estimate, n_curr, n_new_est, p, n_found;
206   Vec                coorlocal;
207   const PetscScalar *_coor;
208   DM                 celldm;
209   Vec                pos;
210   PetscScalar       *_pos;
211   PetscReal         *swarm_coor;
212   PetscInt          *swarm_cellid;
213   PetscSF            sfcell = NULL;
214   const PetscSFNode *LA_sfcell;
215   PetscReal         *my_coor;
216   PetscInt           my_npoints;
217   PetscMPIInt        rank;
218   MPI_Comm           comm;
219 
220   PetscFunctionBegin;
221   DMSWARMPICVALID(dm);
222   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
223   PetscCallMPI(MPI_Comm_rank(comm, &rank));
224 
225   PetscCall(DMSwarmGetCellDM(dm, &celldm));
226   PetscCall(DMGetCoordinatesLocal(celldm, &coorlocal));
227   PetscCall(VecGetSize(coorlocal, &N));
228   PetscCall(VecGetBlockSize(coorlocal, &bs));
229   N = N / bs;
230   PetscCall(VecGetArrayRead(coorlocal, &_coor));
231   for (i = 0; i < N; i++) {
232     for (b = 0; b < bs; b++) {
233       gmin[b] = PetscMin(gmin[b], PetscRealPart(_coor[bs * i + b]));
234       gmax[b] = PetscMax(gmax[b], PetscRealPart(_coor[bs * i + b]));
235     }
236   }
237   PetscCall(VecRestoreArrayRead(coorlocal, &_coor));
238 
239   /* broadcast points from rank 0 if requested */
240   if (redundant) {
241     my_npoints = npoints;
242     PetscCallMPI(MPI_Bcast(&my_npoints, 1, MPIU_INT, 0, comm));
243 
244     if (rank > 0) { /* allocate space */
245       PetscCall(PetscMalloc1(bs * my_npoints, &my_coor));
246     } else {
247       my_coor = coor;
248     }
249     PetscCallMPI(MPI_Bcast(my_coor, bs * my_npoints, MPIU_REAL, 0, comm));
250   } else {
251     my_npoints = npoints;
252     my_coor    = coor;
253   }
254 
255   /* determine the number of points living in the bounding box */
256   n_estimate = 0;
257   for (i = 0; i < my_npoints; i++) {
258     PetscBool point_inside = PETSC_TRUE;
259 
260     for (b = 0; b < bs; b++) {
261       if (my_coor[bs * i + b] < gmin[b]) point_inside = PETSC_FALSE;
262       if (my_coor[bs * i + b] > gmax[b]) point_inside = PETSC_FALSE;
263     }
264     if (point_inside) n_estimate++;
265   }
266 
267   /* create candidate list */
268   PetscCall(VecCreate(PETSC_COMM_SELF, &pos));
269   PetscCall(VecSetSizes(pos, bs * n_estimate, PETSC_DECIDE));
270   PetscCall(VecSetBlockSize(pos, bs));
271   PetscCall(VecSetFromOptions(pos));
272   PetscCall(VecGetArray(pos, &_pos));
273 
274   n_estimate = 0;
275   for (i = 0; i < my_npoints; i++) {
276     PetscBool point_inside = PETSC_TRUE;
277 
278     for (b = 0; b < bs; b++) {
279       if (my_coor[bs * i + b] < gmin[b]) point_inside = PETSC_FALSE;
280       if (my_coor[bs * i + b] > gmax[b]) point_inside = PETSC_FALSE;
281     }
282     if (point_inside) {
283       for (b = 0; b < bs; b++) _pos[bs * n_estimate + b] = my_coor[bs * i + b];
284       n_estimate++;
285     }
286   }
287   PetscCall(VecRestoreArray(pos, &_pos));
288 
289   /* locate points */
290   PetscCall(DMLocatePoints(celldm, pos, DM_POINTLOCATION_NONE, &sfcell));
291 
292   PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell));
293   n_found = 0;
294   for (p = 0; p < n_estimate; p++) {
295     if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) n_found++;
296   }
297 
298   /* adjust size */
299   if (mode == ADD_VALUES) {
300     PetscCall(DMSwarmGetLocalSize(dm, &n_curr));
301     n_new_est = n_curr + n_found;
302     PetscCall(DMSwarmSetLocalSizes(dm, n_new_est, -1));
303   }
304   if (mode == INSERT_VALUES) {
305     n_curr    = 0;
306     n_new_est = n_found;
307     PetscCall(DMSwarmSetLocalSizes(dm, n_new_est, -1));
308   }
309 
310   /* initialize new coords, cell owners, pid */
311   PetscCall(VecGetArrayRead(pos, &_coor));
312   PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor));
313   PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid));
314   n_found = 0;
315   for (p = 0; p < n_estimate; p++) {
316     if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) {
317       for (b = 0; b < bs; b++) swarm_coor[bs * (n_curr + n_found) + b] = PetscRealPart(_coor[bs * p + b]);
318       swarm_cellid[n_curr + n_found] = LA_sfcell[p].index;
319       n_found++;
320     }
321   }
322   PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid));
323   PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor));
324   PetscCall(VecRestoreArrayRead(pos, &_coor));
325 
326   if (redundant) {
327     if (rank > 0) PetscCall(PetscFree(my_coor));
328   }
329   PetscCall(PetscSFDestroy(&sfcell));
330   PetscCall(VecDestroy(&pos));
331   PetscFunctionReturn(PETSC_SUCCESS);
332 }
333 
334 extern PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_DA(DM, DM, DMSwarmPICLayoutType, PetscInt);
335 extern PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX(DM, DM, DMSwarmPICLayoutType, PetscInt);
336 
337 /*@C
338    DMSwarmInsertPointsUsingCellDM - Insert point coordinates within each cell
339 
340    Not collective
341 
342    Input parameters:
343 +  dm - the DMSwarm
344 .  layout_type - method used to fill each cell with the cell DM
345 -  fill_param - parameter controlling how many points per cell are added (the meaning of this parameter is dependent on the layout type)
346 
347    Level: beginner
348 
349    Notes:
350 
351    The insert method will reset any previous defined points within the DMSwarm.
352 
353    When using a DMDA both 2D and 3D are supported for all layout types provided you are using DMDA_ELEMENT_Q1.
354 
355    When using a DMPLEX the following case are supported:
356    (i) DMSWARMPIC_LAYOUT_REGULAR: 2D (triangle),
357    (ii) DMSWARMPIC_LAYOUT_GAUSS: 2D and 3D provided the cell is a tri/tet or a quad/hex,
358    (iii) DMSWARMPIC_LAYOUT_SUBDIVISION: 2D and 3D for quad/hex and 2D tri.
359 
360 .seealso: `DMSwarmPICLayoutType`, `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType`
361 @*/
362 PETSC_EXTERN PetscErrorCode DMSwarmInsertPointsUsingCellDM(DM dm, DMSwarmPICLayoutType layout_type, PetscInt fill_param)
363 {
364   DM        celldm;
365   PetscBool isDA, isPLEX;
366 
367   PetscFunctionBegin;
368   DMSWARMPICVALID(dm);
369   PetscCall(DMSwarmGetCellDM(dm, &celldm));
370   PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMDA, &isDA));
371   PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMPLEX, &isPLEX));
372   if (isDA) {
373     PetscCall(private_DMSwarmInsertPointsUsingCellDM_DA(dm, celldm, layout_type, fill_param));
374   } else if (isPLEX) {
375     PetscCall(private_DMSwarmInsertPointsUsingCellDM_PLEX(dm, celldm, layout_type, fill_param));
376   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only supported for cell DMs of type DMDA and DMPLEX");
377   PetscFunctionReturn(PETSC_SUCCESS);
378 }
379 
380 extern PetscErrorCode private_DMSwarmSetPointCoordinatesCellwise_PLEX(DM, DM, PetscInt, PetscReal *);
381 
382 /*@C
383    DMSwarmSetPointCoordinatesCellwise - Insert point coordinates (defined over the reference cell) within each cell
384 
385    Not collective
386 
387    Input parameters:
388 +  dm - the DMSwarm
389 .  celldm - the cell DM
390 .  npoints - the number of points to insert in each cell
391 -  xi - the coordinates (defined in the local coordinate system for each cell) to insert
392 
393  Level: beginner
394 
395  Notes:
396  The method will reset any previous defined points within the DMSwarm.
397  Only supported for DMPLEX. If you are using a DMDA it is recommended to either use
398  DMSwarmInsertPointsUsingCellDM(), or extract and set the coordinates yourself the following code
399 
400 $    PetscReal *coor;
401 $    DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&coor);
402 $    // user code to define the coordinates here
403 $    DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&coor);
404 
405 .seealso: `DMSwarmSetCellDM()`, `DMSwarmInsertPointsUsingCellDM()`
406 @*/
407 PETSC_EXTERN PetscErrorCode DMSwarmSetPointCoordinatesCellwise(DM dm, PetscInt npoints, PetscReal xi[])
408 {
409   DM        celldm;
410   PetscBool isDA, isPLEX;
411 
412   PetscFunctionBegin;
413   DMSWARMPICVALID(dm);
414   PetscCall(DMSwarmGetCellDM(dm, &celldm));
415   PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMDA, &isDA));
416   PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMPLEX, &isPLEX));
417   PetscCheck(!isDA, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only supported for cell DMs of type DMPLEX. Recommended you use DMSwarmInsertPointsUsingCellDM()");
418   if (isPLEX) {
419     PetscCall(private_DMSwarmSetPointCoordinatesCellwise_PLEX(dm, celldm, npoints, xi));
420   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only supported for cell DMs of type DMDA and DMPLEX");
421   PetscFunctionReturn(PETSC_SUCCESS);
422 }
423 
424 /* Field projection API */
425 extern PetscErrorCode private_DMSwarmProjectFields_DA(DM swarm, DM celldm, PetscInt project_type, PetscInt nfields, DMSwarmDataField dfield[], Vec vecs[]);
426 extern PetscErrorCode private_DMSwarmProjectFields_PLEX(DM swarm, DM celldm, PetscInt project_type, PetscInt nfields, DMSwarmDataField dfield[], Vec vecs[]);
427 
428 /*@C
429    DMSwarmProjectFields - Project a set of swarm fields onto the cell DM
430 
431    Collective on dm
432 
433    Input parameters:
434 +  dm - the DMSwarm
435 .  nfields - the number of swarm fields to project
436 .  fieldnames - the textual names of the swarm fields to project
437 .  fields - an array of Vec's of length nfields
438 -  reuse - flag indicating whether the array and contents of fields should be re-used or internally allocated
439 
440    Currently, the only available projection method consists of
441      phi_i = \sum_{p=0}^{np} N_i(x_p) phi_p dJ / \sum_{p=0}^{np} N_i(x_p) dJ
442    where phi_p is the swarm field at point p,
443      N_i() is the cell DM basis function at vertex i,
444      dJ is the determinant of the cell Jacobian and
445      phi_i is the projected vertex value of the field phi.
446 
447    Level: beginner
448 
449    Notes:
450 
451    If reuse = PETSC_FALSE, this function will allocate the array of Vec's, and each individual Vec.
452      The user is responsible for destroying both the array and the individual Vec objects.
453 
454    Only swarm fields registered with data type = PETSC_REAL can be projected onto the cell DM.
455 
456    Only swarm fields of block size = 1 can currently be projected.
457 
458    The only projection methods currently only support the DA (2D) and PLEX (triangles 2D).
459 
460 .seealso: `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType`
461 @*/
462 PETSC_EXTERN PetscErrorCode DMSwarmProjectFields(DM dm, PetscInt nfields, const char *fieldnames[], Vec **fields, PetscBool reuse)
463 {
464   DM_Swarm         *swarm = (DM_Swarm *)dm->data;
465   DMSwarmDataField *gfield;
466   DM                celldm;
467   PetscBool         isDA, isPLEX;
468   Vec              *vecs;
469   PetscInt          f, nvecs;
470   PetscInt          project_type = 0;
471 
472   PetscFunctionBegin;
473   DMSWARMPICVALID(dm);
474   PetscCall(DMSwarmGetCellDM(dm, &celldm));
475   PetscCall(PetscMalloc1(nfields, &gfield));
476   nvecs = 0;
477   for (f = 0; f < nfields; f++) {
478     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldnames[f], &gfield[f]));
479     PetscCheck(gfield[f]->petsc_type == PETSC_REAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Projection only valid for fields using a data type = PETSC_REAL");
480     PetscCheck(gfield[f]->bs == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Projection only valid for fields with block size = 1");
481     nvecs += gfield[f]->bs;
482   }
483   if (!reuse) {
484     PetscCall(PetscMalloc1(nvecs, &vecs));
485     for (f = 0; f < nvecs; f++) {
486       PetscCall(DMCreateGlobalVector(celldm, &vecs[f]));
487       PetscCall(PetscObjectSetName((PetscObject)vecs[f], gfield[f]->name));
488     }
489   } else {
490     vecs = *fields;
491   }
492 
493   PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMDA, &isDA));
494   PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMPLEX, &isPLEX));
495   if (isDA) {
496     PetscCall(private_DMSwarmProjectFields_DA(dm, celldm, project_type, nfields, gfield, vecs));
497   } else if (isPLEX) {
498     PetscCall(private_DMSwarmProjectFields_PLEX(dm, celldm, project_type, nfields, gfield, vecs));
499   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only supported for cell DMs of type DMDA and DMPLEX");
500 
501   PetscCall(PetscFree(gfield));
502   if (!reuse) *fields = vecs;
503   PetscFunctionReturn(PETSC_SUCCESS);
504 }
505 
506 /*@C
507    DMSwarmCreatePointPerCellCount - Count the number of points within all cells in the cell DM
508 
509    Not collective
510 
511    Input parameter:
512 .  dm - the DMSwarm
513 
514    Output parameters:
515 +  ncells - the number of cells in the cell DM (optional argument, pass NULL to ignore)
516 -  count - array of length ncells containing the number of points per cell
517 
518    Level: beginner
519 
520    Notes:
521    The array count is allocated internally and must be free'd by the user.
522 
523 .seealso: `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType`
524 @*/
525 PETSC_EXTERN PetscErrorCode DMSwarmCreatePointPerCellCount(DM dm, PetscInt *ncells, PetscInt **count)
526 {
527   PetscBool isvalid;
528   PetscInt  nel;
529   PetscInt *sum;
530 
531   PetscFunctionBegin;
532   PetscCall(DMSwarmSortGetIsValid(dm, &isvalid));
533   nel = 0;
534   if (isvalid) {
535     PetscInt e;
536 
537     PetscCall(DMSwarmSortGetSizes(dm, &nel, NULL));
538 
539     PetscCall(PetscMalloc1(nel, &sum));
540     for (e = 0; e < nel; e++) PetscCall(DMSwarmSortGetNumberOfPointsPerCell(dm, e, &sum[e]));
541   } else {
542     DM        celldm;
543     PetscBool isda, isplex, isshell;
544     PetscInt  p, npoints;
545     PetscInt *swarm_cellid;
546 
547     /* get the number of cells */
548     PetscCall(DMSwarmGetCellDM(dm, &celldm));
549     PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMDA, &isda));
550     PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMPLEX, &isplex));
551     PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMSHELL, &isshell));
552     if (isda) {
553       PetscInt        _nel, _npe;
554       const PetscInt *_element;
555 
556       PetscCall(DMDAGetElements(celldm, &_nel, &_npe, &_element));
557       nel = _nel;
558       PetscCall(DMDARestoreElements(celldm, &_nel, &_npe, &_element));
559     } else if (isplex) {
560       PetscInt ps, pe;
561 
562       PetscCall(DMPlexGetHeightStratum(celldm, 0, &ps, &pe));
563       nel = pe - ps;
564     } else if (isshell) {
565       PetscErrorCode (*method_DMShellGetNumberOfCells)(DM, PetscInt *);
566 
567       PetscCall(PetscObjectQueryFunction((PetscObject)celldm, "DMGetNumberOfCells_C", &method_DMShellGetNumberOfCells));
568       if (method_DMShellGetNumberOfCells) {
569         PetscCall(method_DMShellGetNumberOfCells(celldm, &nel));
570       } else
571         SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Cannot determine the number of cells for the DMSHELL object. User must provide a method via PetscObjectComposeFunction( (PetscObject)shelldm, \"DMGetNumberOfCells_C\", your_function_to_compute_number_of_cells);");
572     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Cannot determine the number of cells for a DM not of type DA, PLEX or SHELL");
573 
574     PetscCall(PetscMalloc1(nel, &sum));
575     PetscCall(PetscArrayzero(sum, nel));
576     PetscCall(DMSwarmGetLocalSize(dm, &npoints));
577     PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid));
578     for (p = 0; p < npoints; p++) {
579       if (swarm_cellid[p] != DMLOCATEPOINT_POINT_NOT_FOUND) sum[swarm_cellid[p]]++;
580     }
581     PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid));
582   }
583   if (ncells) *ncells = nel;
584   *count = sum;
585   PetscFunctionReturn(PETSC_SUCCESS);
586 }
587 
588 /*@
589   DMSwarmGetNumSpecies - Get the number of particle species
590 
591   Not collective
592 
593   Input parameter:
594 . dm - the DMSwarm
595 
596   Output parameters:
597 . Ns - the number of species
598 
599   Level: intermediate
600 
601 .seealso: `DMSwarmSetNumSpecies()`, `DMSwarmSetType()`, `DMSwarmType`
602 @*/
603 PetscErrorCode DMSwarmGetNumSpecies(DM sw, PetscInt *Ns)
604 {
605   DM_Swarm *swarm = (DM_Swarm *)sw->data;
606 
607   PetscFunctionBegin;
608   *Ns = swarm->Ns;
609   PetscFunctionReturn(PETSC_SUCCESS);
610 }
611 
612 /*@
613   DMSwarmSetNumSpecies - Set the number of particle species
614 
615   Not collective
616 
617   Input parameter:
618 + dm - the DMSwarm
619 - Ns - the number of species
620 
621   Level: intermediate
622 
623 .seealso: `DMSwarmGetNumSpecies()`, `DMSwarmSetType()`, `DMSwarmType`
624 @*/
625 PetscErrorCode DMSwarmSetNumSpecies(DM sw, PetscInt Ns)
626 {
627   DM_Swarm *swarm = (DM_Swarm *)sw->data;
628 
629   PetscFunctionBegin;
630   swarm->Ns = Ns;
631   PetscFunctionReturn(PETSC_SUCCESS);
632 }
633 
634 /*@C
635   DMSwarmGetCoordinateFunction - Get the function setting initial particle positions, if it exists
636 
637   Not collective
638 
639   Input parameter:
640 . dm - the DMSwarm
641 
642   Output Parameter:
643 . coordFunc - the function setting initial particle positions, or NULL
644 
645   Level: intermediate
646 
647 .seealso: `DMSwarmSetCoordinateFunction()`, `DMSwarmGetVelocityFunction()`, `DMSwarmInitializeCoordinates()`
648 @*/
649 PetscErrorCode DMSwarmGetCoordinateFunction(DM sw, PetscSimplePointFunc *coordFunc)
650 {
651   DM_Swarm *swarm = (DM_Swarm *)sw->data;
652 
653   PetscFunctionBegin;
654   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
655   PetscValidPointer(coordFunc, 2);
656   *coordFunc = swarm->coordFunc;
657   PetscFunctionReturn(PETSC_SUCCESS);
658 }
659 
660 /*@C
661   DMSwarmSetCoordinateFunction - Set the function setting initial particle positions
662 
663   Not collective
664 
665   Input parameters:
666 + dm - the DMSwarm
667 - coordFunc - the function setting initial particle positions
668 
669   Level: intermediate
670 
671 .seealso: `DMSwarmGetCoordinateFunction()`, `DMSwarmSetVelocityFunction()`, `DMSwarmInitializeCoordinates()`
672 @*/
673 PetscErrorCode DMSwarmSetCoordinateFunction(DM sw, PetscSimplePointFunc coordFunc)
674 {
675   DM_Swarm *swarm = (DM_Swarm *)sw->data;
676 
677   PetscFunctionBegin;
678   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
679   PetscValidFunction(coordFunc, 2);
680   swarm->coordFunc = coordFunc;
681   PetscFunctionReturn(PETSC_SUCCESS);
682 }
683 
684 /*@C
685   DMSwarmGetCoordinateFunction - Get the function setting initial particle velocities, if it exists
686 
687   Not collective
688 
689   Input parameter:
690 . dm - the DMSwarm
691 
692   Output Parameter:
693 . velFunc - the function setting initial particle velocities, or NULL
694 
695   Level: intermediate
696 
697 .seealso: `DMSwarmSetVelocityFunction()`, `DMSwarmGetCoordinateFunction()`, `DMSwarmInitializeVelocities()`
698 @*/
699 PetscErrorCode DMSwarmGetVelocityFunction(DM sw, PetscSimplePointFunc *velFunc)
700 {
701   DM_Swarm *swarm = (DM_Swarm *)sw->data;
702 
703   PetscFunctionBegin;
704   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
705   PetscValidPointer(velFunc, 2);
706   *velFunc = swarm->velFunc;
707   PetscFunctionReturn(PETSC_SUCCESS);
708 }
709 
710 /*@C
711   DMSwarmSetVelocityFunction - Set the function setting initial particle velocities
712 
713   Not collective
714 
715   Input parameters:
716 + dm - the DMSwarm
717 - coordFunc - the function setting initial particle velocities
718 
719   Level: intermediate
720 
721 .seealso: `DMSwarmGetVelocityFunction()`, `DMSwarmSetCoordinateFunction()`, `DMSwarmInitializeVelocities()`
722 @*/
723 PetscErrorCode DMSwarmSetVelocityFunction(DM sw, PetscSimplePointFunc velFunc)
724 {
725   DM_Swarm *swarm = (DM_Swarm *)sw->data;
726 
727   PetscFunctionBegin;
728   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
729   PetscValidFunction(velFunc, 2);
730   swarm->velFunc = velFunc;
731   PetscFunctionReturn(PETSC_SUCCESS);
732 }
733 
734 /*@C
735   DMSwarmComputeLocalSize - Compute the local number and distribution of particles based upon a density function
736 
737   Not collective
738 
739   Input Parameters:
740 + sw      - The DMSwarm
741 . N       - The target number of particles
742 - density - The density field for the particle layout, normalized to unity
743 
744   Note: One particle will be created for each species.
745 
746   Level: advanced
747 
748 .seealso: `DMSwarmComputeLocalSizeFromOptions()`
749 @*/
750 PetscErrorCode DMSwarmComputeLocalSize(DM sw, PetscInt N, PetscProbFunc density)
751 {
752   DM               dm;
753   PetscQuadrature  quad;
754   const PetscReal *xq, *wq;
755   PetscReal       *n_int;
756   PetscInt        *npc_s, *cellid, Ni;
757   PetscReal        gmin[3], gmax[3], xi0[3];
758   PetscInt         Ns, cStart, cEnd, c, dim, d, Nq, q, Np = 0, p, s;
759   PetscBool        simplex;
760 
761   PetscFunctionBegin;
762   PetscCall(DMSwarmGetNumSpecies(sw, &Ns));
763   PetscCall(DMSwarmGetCellDM(sw, &dm));
764   PetscCall(DMGetDimension(dm, &dim));
765   PetscCall(DMGetBoundingBox(dm, gmin, gmax));
766   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
767   PetscCall(DMPlexIsSimplex(dm, &simplex));
768   PetscCall(DMGetCoordinatesLocalSetUp(dm));
769   if (simplex) PetscCall(PetscDTStroudConicalQuadrature(dim, 1, 5, -1.0, 1.0, &quad));
770   else PetscCall(PetscDTGaussTensorQuadrature(dim, 1, 5, -1.0, 1.0, &quad));
771   PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, &xq, &wq));
772   PetscCall(PetscCalloc2(Ns, &n_int, (cEnd - cStart) * Ns, &npc_s));
773   /* Integrate the density function to get the number of particles in each cell */
774   for (d = 0; d < dim; ++d) xi0[d] = -1.0;
775   for (c = 0; c < cEnd - cStart; ++c) {
776     const PetscInt cell = c + cStart;
777     PetscReal      v0[3], J[9], invJ[9], detJ, detJp = 2. / (gmax[0] - gmin[0]), xr[3], den;
778 
779     /*Have to transform quadrature points/weights to cell domain*/
780     PetscCall(DMPlexComputeCellGeometryFEM(dm, cell, NULL, v0, J, invJ, &detJ));
781     PetscCall(PetscArrayzero(n_int, Ns));
782     for (q = 0; q < Nq; ++q) {
783       CoordinatesRefToReal(dim, dim, xi0, v0, J, &xq[q * dim], xr);
784       /* Have to transform mesh to domain of definition of PDF, [-1, 1], and weight PDF by |J|/2 */
785       xr[0] = detJp * (xr[0] - gmin[0]) - 1.;
786 
787       for (s = 0; s < Ns; ++s) {
788         PetscCall(density(xr, NULL, &den));
789         n_int[s] += (detJp * den) * (detJ * wq[q]) / (PetscReal)Ns;
790       }
791     }
792     for (s = 0; s < Ns; ++s) {
793       Ni = N;
794       npc_s[c * Ns + s] += (PetscInt)(Ni * n_int[s]);
795       Np += npc_s[c * Ns + s];
796     }
797   }
798   PetscCall(PetscQuadratureDestroy(&quad));
799   PetscCall(DMSwarmSetLocalSizes(sw, Np, 0));
800   PetscCall(DMSwarmGetField(sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid));
801   for (c = 0, p = 0; c < cEnd - cStart; ++c) {
802     for (s = 0; s < Ns; ++s) {
803       for (q = 0; q < npc_s[c * Ns + s]; ++q, ++p) cellid[p] = c;
804     }
805   }
806   PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid));
807   PetscCall(PetscFree2(n_int, npc_s));
808   PetscFunctionReturn(PETSC_SUCCESS);
809 }
810 
811 /*@
812   DMSwarmComputeLocalSizeFromOptions - Compute the local number and distribution of particles based upon a density function determined by options
813 
814   Not collective
815 
816   Input Parameters:
817 , sw - The DMSwarm
818 
819   Level: advanced
820 
821 .seealso: `DMSwarmComputeLocalSize()`
822 @*/
823 PetscErrorCode DMSwarmComputeLocalSizeFromOptions(DM sw)
824 {
825   PetscProbFunc pdf;
826   const char   *prefix;
827   char          funcname[PETSC_MAX_PATH_LEN];
828   PetscInt     *N, Ns, dim, n;
829   PetscBool     flg;
830   PetscMPIInt   size, rank;
831 
832   PetscFunctionBegin;
833   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)sw), &size));
834   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sw), &rank));
835   PetscCall(PetscCalloc1(size, &N));
836   PetscOptionsBegin(PetscObjectComm((PetscObject)sw), "", "DMSwarm Options", "DMSWARM");
837   n = size;
838   PetscCall(PetscOptionsIntArray("-dm_swarm_num_particles", "The target number of particles", "", N, &n, NULL));
839   PetscCall(DMSwarmGetNumSpecies(sw, &Ns));
840   PetscCall(PetscOptionsInt("-dm_swarm_num_species", "The number of species", "DMSwarmSetNumSpecies", Ns, &Ns, &flg));
841   if (flg) PetscCall(DMSwarmSetNumSpecies(sw, Ns));
842   PetscCall(PetscOptionsString("-dm_swarm_coordinate_function", "Function to determine particle coordinates", "DMSwarmSetCoordinateFunction", funcname, funcname, sizeof(funcname), &flg));
843   PetscOptionsEnd();
844   if (flg) {
845     PetscSimplePointFunc coordFunc;
846 
847     PetscCall(DMSwarmGetNumSpecies(sw, &Ns));
848     PetscCall(PetscDLSym(NULL, funcname, (void **)&coordFunc));
849     PetscCheck(coordFunc, PetscObjectComm((PetscObject)sw), PETSC_ERR_ARG_WRONG, "Could not locate function %s", funcname);
850     PetscCall(DMSwarmGetNumSpecies(sw, &Ns));
851     PetscCall(DMSwarmSetLocalSizes(sw, N[rank] * Ns, 0));
852     PetscCall(DMSwarmSetCoordinateFunction(sw, coordFunc));
853   } else {
854     PetscCall(DMGetDimension(sw, &dim));
855     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)sw, &prefix));
856     PetscCall(PetscProbCreateFromOptions(dim, prefix, "-dm_swarm_coordinate_density", &pdf, NULL, NULL));
857     PetscCall(DMSwarmComputeLocalSize(sw, N[rank], pdf));
858   }
859   PetscCall(PetscFree(N));
860   PetscFunctionReturn(PETSC_SUCCESS);
861 }
862 
863 /*@
864   DMSwarmInitializeCoordinates - Determine the initial coordinates of particles for a PIC method
865 
866   Not collective
867 
868   Input Parameters:
869 , sw - The DMSwarm
870 
871   Note: Currently, we randomly place particles in their assigned cell
872 
873   Level: advanced
874 
875 .seealso: `DMSwarmComputeLocalSize()`, `DMSwarmInitializeVelocities()`
876 @*/
877 PetscErrorCode DMSwarmInitializeCoordinates(DM sw)
878 {
879   PetscSimplePointFunc coordFunc;
880   PetscScalar         *weight;
881   PetscReal           *x;
882   PetscInt            *species;
883   void                *ctx;
884   PetscBool            removePoints = PETSC_TRUE;
885   PetscDataType        dtype;
886   PetscInt             Np, p, Ns, dim, d, bs;
887 
888   PetscFunctionBeginUser;
889   PetscCall(DMGetDimension(sw, &dim));
890   PetscCall(DMSwarmGetLocalSize(sw, &Np));
891   PetscCall(DMSwarmGetNumSpecies(sw, &Ns));
892   PetscCall(DMSwarmGetCoordinateFunction(sw, &coordFunc));
893 
894   PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, &bs, &dtype, (void **)&x));
895   PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&weight));
896   PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species));
897   if (coordFunc) {
898     PetscCall(DMGetApplicationContext(sw, &ctx));
899     for (p = 0; p < Np; ++p) {
900       PetscScalar X[3];
901 
902       PetscCall((*coordFunc)(dim, 0., NULL, p, X, ctx));
903       for (d = 0; d < dim; ++d) x[p * dim + d] = PetscRealPart(X[d]);
904       weight[p]  = 1.0;
905       species[p] = p % Ns;
906     }
907   } else {
908     DM          dm;
909     PetscRandom rnd;
910     PetscReal   xi0[3];
911     PetscInt    cStart, cEnd, c;
912 
913     PetscCall(DMSwarmGetCellDM(sw, &dm));
914     PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
915     PetscCall(DMGetApplicationContext(sw, &ctx));
916 
917     /* Set particle position randomly in cell, set weights to 1 */
918     PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rnd));
919     PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0));
920     PetscCall(PetscRandomSetFromOptions(rnd));
921     PetscCall(DMSwarmSortGetAccess(sw));
922     for (d = 0; d < dim; ++d) xi0[d] = -1.0;
923     for (c = cStart; c < cEnd; ++c) {
924       PetscReal v0[3], J[9], invJ[9], detJ;
925       PetscInt *pidx, Npc, q;
926 
927       PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &Npc, &pidx));
928       PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ));
929       for (q = 0; q < Npc; ++q) {
930         const PetscInt p = pidx[q];
931         PetscReal      xref[3];
932 
933         for (d = 0; d < dim; ++d) PetscCall(PetscRandomGetValueReal(rnd, &xref[d]));
934         CoordinatesRefToReal(dim, dim, xi0, v0, J, xref, &x[p * dim]);
935 
936         weight[p]  = 1.0 / Np;
937         species[p] = p % Ns;
938       }
939       PetscCall(PetscFree(pidx));
940     }
941     PetscCall(PetscRandomDestroy(&rnd));
942     PetscCall(DMSwarmSortRestoreAccess(sw));
943   }
944   PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&x));
945   PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&weight));
946   PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species));
947 
948   PetscCall(DMSwarmMigrate(sw, removePoints));
949   PetscCall(DMLocalizeCoordinates(sw));
950   PetscFunctionReturn(PETSC_SUCCESS);
951 }
952 
953 /*@C
954   DMSwarmInitializeVelocities - Set the initial velocities of particles using a distribution.
955 
956   Collective on dm
957 
958   Input Parameters:
959 + sw      - The DMSwarm object
960 . sampler - A function which uniformly samples the velocity PDF
961 - v0      - The velocity scale for nondimensionalization for each species
962 
963   Note: If v0 is zero for the first species, all velocities are set to zero. If it is zero for any other species, the effect will be to give that species zero velocity.
964 
965   Level: advanced
966 
967 .seealso: `DMSwarmComputeLocalSize()`, `DMSwarmInitializeCoordinates()`, `DMSwarmInitializeVelocitiesFromOptions()`
968 @*/
969 PetscErrorCode DMSwarmInitializeVelocities(DM sw, PetscProbFunc sampler, const PetscReal v0[])
970 {
971   PetscSimplePointFunc velFunc;
972   PetscReal           *v;
973   PetscInt            *species;
974   void                *ctx;
975   PetscInt             dim, Np, p;
976 
977   PetscFunctionBegin;
978   PetscCall(DMSwarmGetVelocityFunction(sw, &velFunc));
979 
980   PetscCall(DMGetDimension(sw, &dim));
981   PetscCall(DMSwarmGetLocalSize(sw, &Np));
982   PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&v));
983   PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species));
984   if (v0[0] == 0.) {
985     PetscCall(PetscArrayzero(v, Np * dim));
986   } else if (velFunc) {
987     PetscCall(DMGetApplicationContext(sw, &ctx));
988     for (p = 0; p < Np; ++p) {
989       PetscInt    s = species[p], d;
990       PetscScalar vel[3];
991 
992       PetscCall((*velFunc)(dim, 0., NULL, p, vel, ctx));
993       for (d = 0; d < dim; ++d) v[p * dim + d] = (v0[s] / v0[0]) * PetscRealPart(vel[d]);
994     }
995   } else {
996     PetscRandom rnd;
997 
998     PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)sw), &rnd));
999     PetscCall(PetscRandomSetInterval(rnd, 0, 1.));
1000     PetscCall(PetscRandomSetFromOptions(rnd));
1001 
1002     for (p = 0; p < Np; ++p) {
1003       PetscInt  s = species[p], d;
1004       PetscReal a[3], vel[3];
1005 
1006       for (d = 0; d < dim; ++d) PetscCall(PetscRandomGetValueReal(rnd, &a[d]));
1007       PetscCall(sampler(a, NULL, vel));
1008       for (d = 0; d < dim; ++d) v[p * dim + d] = (v0[s] / v0[0]) * vel[d];
1009     }
1010     PetscCall(PetscRandomDestroy(&rnd));
1011   }
1012   PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&v));
1013   PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species));
1014   PetscFunctionReturn(PETSC_SUCCESS);
1015 }
1016 
1017 /*@
1018   DMSwarmInitializeVelocitiesFromOptions - Set the initial velocities of particles using a distribution determined from options.
1019 
1020   Collective on dm
1021 
1022   Input Parameters:
1023 + sw      - The DMSwarm object
1024 - v0      - The velocity scale for nondimensionalization for each species
1025 
1026   Level: advanced
1027 
1028 .seealso: `DMSwarmComputeLocalSize()`, `DMSwarmInitializeCoordinates()`, `DMSwarmInitializeVelocities()`
1029 @*/
1030 PetscErrorCode DMSwarmInitializeVelocitiesFromOptions(DM sw, const PetscReal v0[])
1031 {
1032   PetscProbFunc sampler;
1033   PetscInt      dim;
1034   const char   *prefix;
1035   char          funcname[PETSC_MAX_PATH_LEN];
1036   PetscBool     flg;
1037 
1038   PetscFunctionBegin;
1039   PetscOptionsBegin(PetscObjectComm((PetscObject)sw), "", "DMSwarm Options", "DMSWARM");
1040   PetscCall(PetscOptionsString("-dm_swarm_velocity_function", "Function to determine particle velocities", "DMSwarmSetVelocityFunction", funcname, funcname, sizeof(funcname), &flg));
1041   PetscOptionsEnd();
1042   if (flg) {
1043     PetscSimplePointFunc velFunc;
1044 
1045     PetscCall(PetscDLSym(NULL, funcname, (void **)&velFunc));
1046     PetscCheck(velFunc, PetscObjectComm((PetscObject)sw), PETSC_ERR_ARG_WRONG, "Could not locate function %s", funcname);
1047     PetscCall(DMSwarmSetVelocityFunction(sw, velFunc));
1048   }
1049   PetscCall(DMGetDimension(sw, &dim));
1050   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)sw, &prefix));
1051   PetscCall(PetscProbCreateFromOptions(dim, prefix, "-dm_swarm_velocity_density", NULL, NULL, &sampler));
1052   PetscCall(DMSwarmInitializeVelocities(sw, sampler, v0));
1053   PetscFunctionReturn(PETSC_SUCCESS);
1054 }
1055