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