xref: /petsc/src/dm/impls/swarm/swarmpic.c (revision 3d96e9964ff330fd2a9eee374bcd4376da7efe60)
1 
2 #define PETSCDM_DLL
3 #include <petsc/private/dmswarmimpl.h>    /*I   "petscdmswarm.h"   I*/
4 #include <petscsf.h>
5 #include <petscdmda.h>
6 #include <petscdmplex.h>
7 #include "data_bucket.h"
8 
9 /*
10  Error chceking macto to ensure the swarm type is correct and that a cell DM has been set
11 */
12 #define DMSWARMPICVALID(dm) \
13 { \
14   DM_Swarm *_swarm = (DM_Swarm*)(dm)->data; \
15   if (_swarm->swarm_type != DMSWARM_PIC) SETERRQ(PetscObjectComm((PetscObject)(dm)),PETSC_ERR_SUP,"Only valid for DMSwarm-PIC. You must call DMSwarmSetType(dm,DMSWARM_PIC)"); \
16   else \
17     if (!_swarm->dmcell) SETERRQ(PetscObjectComm((PetscObject)(dm)),PETSC_ERR_SUP,"Only valid for DMSwarmPIC if the cell DM is set. You must call DMSwarmSetCellDM(dm,celldm)"); \
18 }
19 
20 /* Coordinate insertition/addition API */
21 /*@C
22    DMSwarmSetPointsUniformCoordinates - Set point coordinates in a DMSwarm on a regular (ijk) grid
23 
24    Collective on DM
25 
26    Input parameters:
27 +  dm - the DMSwarm
28 .  min - minimum coordinate values in the x, y, z directions (array of length dim)
29 .  max - maximum coordinate values in the x, y, z directions (array of length dim)
30 .  npoints - number of points in each spatial direction (array of length dim)
31 -  mode - indicates whether to append points to the swarm (ADD_VALUES), or over-ride existing points (INSERT_VALUES)
32 
33    Level: beginner
34 
35    Notes:
36    When using mode = INSERT_VALUES, this method will reset the number of particles in the DMSwarm
37    to be npoints[0]*npoints[1] (2D) or npoints[0]*npoints[1]*npoints[2] (3D). When using mode = ADD_VALUES,
38    new points will be appended to any already existing in the DMSwarm
39 
40 .seealso: DMSwarmSetType(), DMSwarmSetCellDM(), DMSwarmType
41 @*/
42 PETSC_EXTERN PetscErrorCode DMSwarmSetPointsUniformCoordinates(DM dm,PetscReal min[],PetscReal max[],PetscInt npoints[],InsertMode mode)
43 {
44   PetscErrorCode ierr;
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   ierr = DMSwarmGetCellDM(dm,&celldm);CHKERRQ(ierr);
63   ierr = DMGetCoordinatesLocal(celldm,&coorlocal);CHKERRQ(ierr);
64   ierr = VecGetSize(coorlocal,&N);CHKERRQ(ierr);
65   ierr = VecGetBlockSize(coorlocal,&bs);CHKERRQ(ierr);
66   N = N / bs;
67   ierr = VecGetArrayRead(coorlocal,&_coor);CHKERRQ(ierr);
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   ierr = VecRestoreArrayRead(coorlocal,&_coor);CHKERRQ(ierr);
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 
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   ierr = VecCreate(PETSC_COMM_SELF,&pos);CHKERRQ(ierr);
112   ierr = VecSetSizes(pos,bs*n_estimate,PETSC_DECIDE);CHKERRQ(ierr);
113   ierr = VecSetBlockSize(pos,bs);CHKERRQ(ierr);
114   ierr = VecSetFromOptions(pos);CHKERRQ(ierr);
115   ierr = VecGetArray(pos,&_pos);CHKERRQ(ierr);
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   ierr = VecRestoreArray(pos,&_pos);CHKERRQ(ierr);
145 
146   /* locate points */
147   ierr = DMLocatePoints(celldm,pos,DM_POINTLOCATION_NONE,&sfcell);CHKERRQ(ierr);
148 
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 
190   PetscFunctionReturn(0);
191 }
192 
193 /*@C
194    DMSwarmSetPointCoordinates - Set point coordinates in a DMSwarm from a user defined list
195 
196    Collective on DM
197 
198    Input parameters:
199 +  dm - the DMSwarm
200 .  npoints - the number of points to insert
201 .  coor - the coordinate values
202 .  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
203 -  mode - indicates whether to append points to the swarm (ADD_VALUES), or over-ride existing points (INSERT_VALUES)
204 
205    Level: beginner
206 
207    Notes:
208    If the user has specified redundant = PETSC_FALSE, the cell DM will attempt to locate the coordinates provided by coor[] within
209    its sub-domain. If they any values within coor[] are not located in the sub-domain, they will be ignored and will not get
210    added to the DMSwarm.
211 
212 .seealso: DMSwarmSetType(), DMSwarmSetCellDM(), DMSwarmType, DMSwarmSetPointsUniformCoordinates()
213 @*/
214 PETSC_EXTERN PetscErrorCode DMSwarmSetPointCoordinates(DM dm,PetscInt npoints,PetscReal coor[],PetscBool redundant,InsertMode mode)
215 {
216   PetscErrorCode ierr;
217   PetscReal gmin[] = {PETSC_MAX_REAL ,PETSC_MAX_REAL, PETSC_MAX_REAL};
218   PetscReal gmax[] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL};
219   PetscInt i,N,bs,b,n_estimate,n_curr,n_new_est,p,n_found;
220   Vec coorlocal;
221   const PetscScalar *_coor;
222   DM celldm;
223   Vec pos;
224   PetscScalar *_pos;
225   PetscReal *swarm_coor;
226   PetscInt *swarm_cellid;
227   PetscSF sfcell = NULL;
228   const PetscSFNode *LA_sfcell;
229   PetscReal *my_coor;
230   PetscInt my_npoints;
231   PetscMPIInt rank;
232   MPI_Comm comm;
233 
234   PetscFunctionBegin;
235   DMSWARMPICVALID(dm);
236   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
237   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
238 
239   ierr = DMSwarmGetCellDM(dm,&celldm);CHKERRQ(ierr);
240   ierr = DMGetCoordinatesLocal(celldm,&coorlocal);CHKERRQ(ierr);
241   ierr = VecGetSize(coorlocal,&N);CHKERRQ(ierr);
242   ierr = VecGetBlockSize(coorlocal,&bs);CHKERRQ(ierr);
243   N = N / bs;
244   ierr = VecGetArrayRead(coorlocal,&_coor);CHKERRQ(ierr);
245   for (i=0; i<N; i++) {
246     for (b=0; b<bs; b++) {
247       gmin[b] = PetscMin(gmin[b],PetscRealPart(_coor[bs*i+b]));
248       gmax[b] = PetscMax(gmax[b],PetscRealPart(_coor[bs*i+b]));
249     }
250   }
251   ierr = VecRestoreArrayRead(coorlocal,&_coor);CHKERRQ(ierr);
252 
253   /* broadcast points from rank 0 if requested */
254   if (redundant) {
255     my_npoints = npoints;
256     ierr = MPI_Bcast(&my_npoints,1,MPIU_INT,0,comm);CHKERRQ(ierr);
257 
258     if (rank > 0) { /* allocate space */
259       ierr = PetscMalloc1(bs*my_npoints,&my_coor);CHKERRQ(ierr);
260     } else {
261       my_coor = coor;
262     }
263     ierr = MPI_Bcast(my_coor,bs*my_npoints,MPIU_REAL,0,comm);CHKERRQ(ierr);
264   } else {
265     my_npoints = npoints;
266     my_coor = coor;
267   }
268 
269   /* determine the number of points living in the bounding box */
270   n_estimate = 0;
271   for (i=0; i<my_npoints; i++) {
272     PetscBool point_inside = PETSC_TRUE;
273 
274     for (b=0; b<bs; b++) {
275       if (my_coor[bs*i+b] < gmin[b]) { point_inside = PETSC_FALSE; }
276       if (my_coor[bs*i+b] > gmax[b]) { point_inside = PETSC_FALSE; }
277     }
278     if (point_inside) { n_estimate++; }
279   }
280 
281   /* create candidate list */
282   ierr = VecCreate(PETSC_COMM_SELF,&pos);CHKERRQ(ierr);
283   ierr = VecSetSizes(pos,bs*n_estimate,PETSC_DECIDE);CHKERRQ(ierr);
284   ierr = VecSetBlockSize(pos,bs);CHKERRQ(ierr);
285   ierr = VecSetFromOptions(pos);CHKERRQ(ierr);
286   ierr = VecGetArray(pos,&_pos);CHKERRQ(ierr);
287 
288   n_estimate = 0;
289   for (i=0; i<my_npoints; i++) {
290     PetscBool point_inside = PETSC_TRUE;
291 
292     for (b=0; b<bs; b++) {
293       if (my_coor[bs*i+b] < gmin[b]) { point_inside = PETSC_FALSE; }
294       if (my_coor[bs*i+b] > gmax[b]) { point_inside = PETSC_FALSE; }
295     }
296     if (point_inside) {
297       for (b=0; b<bs; b++) {
298         _pos[bs*n_estimate+b] = my_coor[bs*i+b];
299       }
300       n_estimate++;
301     }
302   }
303   ierr = VecRestoreArray(pos,&_pos);CHKERRQ(ierr);
304 
305   /* locate points */
306   ierr = DMLocatePoints(celldm,pos,DM_POINTLOCATION_NONE,&sfcell);CHKERRQ(ierr);
307 
308   ierr = PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell);CHKERRQ(ierr);
309   n_found = 0;
310   for (p=0; p<n_estimate; p++) {
311     if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) {
312       n_found++;
313     }
314   }
315 
316   /* adjust size */
317   if (mode == ADD_VALUES) {
318     ierr = DMSwarmGetLocalSize(dm,&n_curr);CHKERRQ(ierr);
319     n_new_est = n_curr + n_found;
320     ierr = DMSwarmSetLocalSizes(dm,n_new_est,-1);CHKERRQ(ierr);
321   }
322   if (mode == INSERT_VALUES) {
323     n_curr = 0;
324     n_new_est = n_found;
325     ierr = DMSwarmSetLocalSizes(dm,n_new_est,-1);CHKERRQ(ierr);
326   }
327 
328   /* initialize new coords, cell owners, pid */
329   ierr = VecGetArrayRead(pos,&_coor);CHKERRQ(ierr);
330   ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr);
331   ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr);
332   n_found = 0;
333   for (p=0; p<n_estimate; p++) {
334     if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) {
335       for (b=0; b<bs; b++) {
336         swarm_coor[bs*(n_curr + n_found) + b] = PetscRealPart(_coor[bs*p+b]);
337       }
338       swarm_cellid[n_curr + n_found] = LA_sfcell[p].index;
339       n_found++;
340     }
341   }
342   ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr);
343   ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr);
344   ierr = VecRestoreArrayRead(pos,&_coor);CHKERRQ(ierr);
345 
346   if (redundant) {
347     if (rank > 0) {
348       ierr = PetscFree(my_coor);CHKERRQ(ierr);
349     }
350   }
351   ierr = PetscSFDestroy(&sfcell);CHKERRQ(ierr);
352   ierr = VecDestroy(&pos);CHKERRQ(ierr);
353 
354   PetscFunctionReturn(0);
355 }
356 
357 extern PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_DA(DM,DM,DMSwarmPICLayoutType,PetscInt);
358 extern PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX(DM,DM,DMSwarmPICLayoutType,PetscInt);
359 
360 /*@C
361    DMSwarmInsertPointsUsingCellDM - Insert point coordinates within each cell
362 
363    Not collective
364 
365    Input parameters:
366 +  dm - the DMSwarm
367 .  layout_type - method used to fill each cell with the cell DM
368 -  fill_param - parameter controlling how many points per cell are added (the meaning of this parameter is dependent on the layout type)
369 
370  Level: beginner
371 
372  Notes:
373  The insert method will reset any previous defined points within the DMSwarm
374 
375 .seealso: DMSwarmPICLayoutType, DMSwarmSetType(), DMSwarmSetCellDM(), DMSwarmType
376 @*/
377 PETSC_EXTERN PetscErrorCode DMSwarmInsertPointsUsingCellDM(DM dm,DMSwarmPICLayoutType layout_type,PetscInt fill_param)
378 {
379   PetscErrorCode ierr;
380   DM celldm;
381   PetscBool isDA,isPLEX;
382 
383   PetscFunctionBegin;
384   DMSWARMPICVALID(dm);
385   ierr = DMSwarmGetCellDM(dm,&celldm);CHKERRQ(ierr);
386   ierr = PetscObjectTypeCompare((PetscObject)celldm,DMDA,&isDA);CHKERRQ(ierr);
387   ierr = PetscObjectTypeCompare((PetscObject)celldm,DMPLEX,&isPLEX);CHKERRQ(ierr);
388   if (isDA) {
389     ierr = private_DMSwarmInsertPointsUsingCellDM_DA(dm,celldm,layout_type,fill_param);CHKERRQ(ierr);
390   } else if (isPLEX) {
391     ierr = private_DMSwarmInsertPointsUsingCellDM_PLEX(dm,celldm,layout_type,fill_param);CHKERRQ(ierr);
392   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only supported for cell DMs of type DMDA and DMPLEX");
393 
394   PetscFunctionReturn(0);
395 }
396 
397 /*
398 PETSC_EXTERN PetscErrorCode DMSwarmAddPointCoordinatesCellWise(DM dm,PetscInt cell,PetscInt npoints,PetscReal xi[],PetscBool proximity_initialization)
399 {
400   PetscFunctionBegin;
401   PetscFunctionReturn(0);
402 }
403 */
404 
405 /* Field projection API */
406 extern PetscErrorCode private_DMSwarmProjectFields_DA(DM swarm,DM celldm,PetscInt project_type,PetscInt nfields,DataField dfield[],Vec vecs[]);
407 extern PetscErrorCode private_DMSwarmProjectFields_PLEX(DM swarm,DM celldm,PetscInt project_type,PetscInt nfields,DataField dfield[],Vec vecs[]);
408 
409 /*@C
410    DMSwarmProjectFields - Project a set of swarm fields onto the cell DM
411 
412    Collective on Vec
413 
414    Input parameters:
415 +  dm - the DMSwarm
416 .  nfields - the number of swarm fields to project
417 .  fieldnames - the textual names of the swarm fields to project
418 .  fields - an array of Vec's of length nfields
419 -  reuse - flag indicating whether the array and contents of fields should be re-used or internally allocated
420 
421    Currently, the only available projection method consists of
422      phi_i = \sum_{p=0}^{np} N_i(x_p) phi_p dJ / \sum_{p=0}^{np} N_i(x_p) dJ
423    where phi_p is the swarm field at point p,
424      N_i() is the cell DM basis function at vertex i,
425      dJ is the determinant of the cell Jacobian and
426      phi_i is the projected vertex value of the field phi.
427 
428    Level: beginner
429 
430    Notes:
431    - If reuse = PETSC_FALSE, this function will allocate the array of Vec's, and each individual Vec.
432      The user is responsible for destroying both the array and the individual Vec objects.
433    - Only swarm fields registered with data type = PETSC_REAL can be projected onto the cell DM.
434    - Only swarm fields of block size = 1 can currently be projected.
435    - The only projection methods currently only support the DA (2D) and PLEX (triangles 2D).
436 
437 .seealso: DMSwarmSetType(), DMSwarmSetCellDM(), DMSwarmType
438 @*/
439 PETSC_EXTERN PetscErrorCode DMSwarmProjectFields(DM dm,PetscInt nfields,const char *fieldnames[],Vec **fields,PetscBool reuse)
440 {
441   DM_Swarm *swarm = (DM_Swarm*)dm->data;
442   DataField *gfield;
443   DM celldm;
444   PetscBool isDA,isPLEX;
445   Vec *vecs;
446   PetscInt f,nvecs;
447   PetscInt project_type = 0;
448   PetscErrorCode ierr;
449 
450   PetscFunctionBegin;
451   DMSWARMPICVALID(dm);
452   ierr = DMSwarmGetCellDM(dm,&celldm);CHKERRQ(ierr);
453   ierr = PetscMalloc1(nfields,&gfield);CHKERRQ(ierr);
454   nvecs = 0;
455   for (f=0; f<nfields; f++) {
456     ierr = DataBucketGetDataFieldByName(swarm->db,fieldnames[f],&gfield[f]);CHKERRQ(ierr);
457     if (gfield[f]->petsc_type != PETSC_REAL) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Projection only valid for fields using a data type = PETSC_REAL");
458     if (gfield[f]->bs != 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Projection only valid for fields with block size = 1");
459     nvecs += gfield[f]->bs;
460   }
461   if (!reuse) {
462     ierr = PetscMalloc1(nvecs,&vecs);CHKERRQ(ierr);
463     for (f=0; f<nvecs; f++) {
464       ierr = DMCreateGlobalVector(celldm,&vecs[f]);CHKERRQ(ierr);
465       ierr = PetscObjectSetName((PetscObject)vecs[f],gfield[f]->name);CHKERRQ(ierr);
466     }
467   } else {
468     vecs = *fields;
469   }
470 
471   ierr = PetscObjectTypeCompare((PetscObject)celldm,DMDA,&isDA);CHKERRQ(ierr);
472   ierr = PetscObjectTypeCompare((PetscObject)celldm,DMPLEX,&isPLEX);CHKERRQ(ierr);
473   if (isDA) {
474     ierr = private_DMSwarmProjectFields_DA(dm,celldm,project_type,nfields,gfield,vecs);CHKERRQ(ierr);
475   } else if (isPLEX) {
476     ierr = private_DMSwarmProjectFields_PLEX(dm,celldm,project_type,nfields,gfield,vecs);CHKERRQ(ierr);
477   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only supported for cell DMs of type DMDA and DMPLEX");
478 
479   ierr = PetscFree(gfield);CHKERRQ(ierr);
480   if (!reuse) {
481     *fields = vecs;
482   }
483 
484   PetscFunctionReturn(0);
485 }
486 
487 /*@C
488    DMSwarmCreatePointPerCellCount - Count the number of points within all cells in the cell DM
489 
490    Not collective
491 
492    Input parameter:
493 .  dm - the DMSwarm
494 
495    Output parameters:
496 +  ncells - the number of cells in the cell DM (optional argument, pass NULL to ignore)
497 -  count - array of length ncells containing the number of points per cell
498 
499    Level: beginner
500 
501    Notes:
502    The array count is allocated internally and must be free'd by the user.
503 
504 .seealso: DMSwarmSetType(), DMSwarmSetCellDM(), DMSwarmType
505 @*/
506 PETSC_EXTERN PetscErrorCode DMSwarmCreatePointPerCellCount(DM dm,PetscInt *ncells,PetscInt **count)
507 {
508   PetscErrorCode ierr;
509   PetscBool      isvalid;
510   PetscInt       nel;
511   PetscInt       *sum;
512 
513   PetscFunctionBegin;
514   ierr = DMSwarmSortGetIsValid(dm,&isvalid);CHKERRQ(ierr);
515   nel = 0;
516   if (isvalid) {
517     PetscInt e;
518 
519     ierr = DMSwarmSortGetSizes(dm,&nel,NULL);CHKERRQ(ierr);
520 
521     ierr = PetscMalloc1(nel,&sum);CHKERRQ(ierr);
522     for (e=0; e<nel; e++) {
523       ierr = DMSwarmSortGetNumberOfPointsPerCell(dm,e,&sum[e]);CHKERRQ(ierr);
524     }
525   } else {
526     DM        celldm;
527     PetscBool isda,isplex,isshell;
528     PetscInt  p,npoints;
529     PetscInt *swarm_cellid;
530 
531     /* get the number of cells */
532     ierr = DMSwarmGetCellDM(dm,&celldm);CHKERRQ(ierr);
533     ierr = PetscObjectTypeCompare((PetscObject)celldm,DMDA,&isda);CHKERRQ(ierr);
534     ierr = PetscObjectTypeCompare((PetscObject)celldm,DMPLEX,&isplex);CHKERRQ(ierr);
535     ierr = PetscObjectTypeCompare((PetscObject)celldm,DMSHELL,&isshell);CHKERRQ(ierr);
536     if (isda) {
537       PetscInt _nel,_npe;
538       const PetscInt *_element;
539 
540       ierr = DMDAGetElements(celldm,&_nel,&_npe,&_element);CHKERRQ(ierr);
541       nel = _nel;
542       ierr = DMDARestoreElements(celldm,&_nel,&_npe,&_element);CHKERRQ(ierr);
543     } else if (isplex) {
544       PetscInt ps,pe;
545 
546       ierr = DMPlexGetHeightStratum(celldm,0,&ps,&pe);CHKERRQ(ierr);
547       nel = pe - ps;
548     } else if (isshell) {
549       PetscErrorCode (*method_DMShellGetNumberOfCells)(DM,PetscInt*);
550 
551       ierr = PetscObjectQueryFunction((PetscObject)celldm,"DMGetNumberOfCells_C",&method_DMShellGetNumberOfCells);CHKERRQ(ierr);
552       if (method_DMShellGetNumberOfCells) {
553         ierr = method_DMShellGetNumberOfCells(celldm,&nel);CHKERRQ(ierr);
554       } 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 );");
555     } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot determine the number of cells for a DM not of type DA, PLEX or SHELL");
556 
557     ierr = PetscMalloc1(nel,&sum);CHKERRQ(ierr);
558     ierr = PetscMemzero(sum,sizeof(PetscInt)*nel);CHKERRQ(ierr);
559     ierr = DMSwarmGetLocalSize(dm,&npoints);CHKERRQ(ierr);
560     ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr);
561     for (p=0; p<npoints; p++) {
562       if (swarm_cellid[p] != DMLOCATEPOINT_POINT_NOT_FOUND) {
563         sum[ swarm_cellid[p] ]++;
564       }
565     }
566     ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr);
567   }
568   if (ncells) { *ncells = nel; }
569   *count  = sum;
570   PetscFunctionReturn(0);
571 }
572