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