xref: /petsc/src/dm/impls/swarm/swarm_migrate.c (revision b16650c811c534e32f8711137bdce5226ab7d716)
1df21e3a8SDave May 
2df21e3a8SDave May #include <petsc/private/dmswarmimpl.h>    /*I   "petscdmswarm.h"   I*/
3df21e3a8SDave May #include "data_bucket.h"
4df21e3a8SDave May #include "data_ex.h"
5df21e3a8SDave May 
6df21e3a8SDave May 
7df21e3a8SDave May #undef __FUNCT__
8df21e3a8SDave May #define __FUNCT__ "DMSwarmMigrate_Push_Basic"
9df21e3a8SDave May PetscErrorCode DMSwarmMigrate_Push_Basic(DM dm,PetscBool remove_sent_points)
10df21e3a8SDave May {
11df21e3a8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
12df21e3a8SDave May   PetscErrorCode ierr;
13df21e3a8SDave May   DataEx de;
14df21e3a8SDave May   PetscInt p,npoints,*rankval,n_points_recv;
15df21e3a8SDave May   PetscMPIInt rank,nrank;
16df21e3a8SDave May   void *point_buffer,*recv_points;
17df21e3a8SDave May   size_t sizeof_dmswarm_point;
18df21e3a8SDave May 
19df21e3a8SDave May   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
20df21e3a8SDave May 
21df21e3a8SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
22df21e3a8SDave May   ierr = DMSwarmGetField(dm,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
23df21e3a8SDave May 
24df21e3a8SDave May   de = DataExCreate(PetscObjectComm((PetscObject)dm),0);CHKERRQ(ierr);
25df21e3a8SDave May 
26df21e3a8SDave May   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
27df21e3a8SDave May   for (p=0; p<npoints; p++) {
28df21e3a8SDave May     nrank = rankval[p];
29df21e3a8SDave May     if (nrank != rank) {
30df21e3a8SDave May       ierr = DataExTopologyAddNeighbour(de,nrank);CHKERRQ(ierr);
31df21e3a8SDave May     }
32df21e3a8SDave May   }
33df21e3a8SDave May   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
34df21e3a8SDave May 
35df21e3a8SDave May   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
36df21e3a8SDave May   for (p=0; p<npoints; p++) {
37df21e3a8SDave May     nrank = rankval[p];
38df21e3a8SDave May     if (nrank != rank) {
39df21e3a8SDave May       ierr = DataExAddToSendCount(de,nrank,1);CHKERRQ(ierr);
40df21e3a8SDave May     }
41df21e3a8SDave May   }
42df21e3a8SDave May   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
43df21e3a8SDave May 
44df21e3a8SDave May   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
45df21e3a8SDave May   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
46df21e3a8SDave May   for (p=0; p<npoints; p++) {
47df21e3a8SDave May     nrank = rankval[p];
48df21e3a8SDave May     if (nrank != rank) {
49df21e3a8SDave May       /* copy point into buffer */
50df21e3a8SDave May       ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr);
51df21e3a8SDave May       /* insert point buffer into DataExchanger */
52df21e3a8SDave May       ierr = DataExPackData(de,nrank,1,point_buffer);CHKERRQ(ierr);
53df21e3a8SDave May     }
54df21e3a8SDave May   }
55df21e3a8SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
56df21e3a8SDave May 
57df21e3a8SDave May 
58df21e3a8SDave May   if (remove_sent_points) {
59df21e3a8SDave May     /* remove points which left processor */
60df21e3a8SDave May     ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
61df21e3a8SDave May     for (p=0; p<npoints; p++) {
62df21e3a8SDave May       nrank = rankval[p];
63df21e3a8SDave May       if (nrank != rank) {
64df21e3a8SDave May         /* kill point */
65df21e3a8SDave May         ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr);
66df21e3a8SDave May         DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */
67df21e3a8SDave May         p--; /* check replacement point */
68df21e3a8SDave May       }
69df21e3a8SDave May     }
70df21e3a8SDave May   }
71df21e3a8SDave May   ierr = DMSwarmRestoreField(dm,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
72df21e3a8SDave May 
73df21e3a8SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
74df21e3a8SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
75df21e3a8SDave May 
76df21e3a8SDave May   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
77df21e3a8SDave May 
78df21e3a8SDave May 	ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
79df21e3a8SDave May 	ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr);
80df21e3a8SDave May 	for (p=0; p<n_points_recv; p++) {
81df21e3a8SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
82df21e3a8SDave May 
83df21e3a8SDave May     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
84df21e3a8SDave May   }
85df21e3a8SDave May   ierr = DataExView(de);CHKERRQ(ierr);
86df21e3a8SDave May   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
87df21e3a8SDave May   ierr = DataExDestroy(de);CHKERRQ(ierr);
88df21e3a8SDave May 
89df21e3a8SDave May   PetscFunctionReturn(0);
90df21e3a8SDave May }
912712d1f2SDave May 
922712d1f2SDave May #undef __FUNCT__
932712d1f2SDave May #define __FUNCT__ "DMSwarmMigrate_GlobalToLocal_Basic"
942712d1f2SDave May PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm,PetscInt *globalsize)
952712d1f2SDave May {
962712d1f2SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
972712d1f2SDave May   PetscErrorCode ierr;
982712d1f2SDave May   DataEx de;
992712d1f2SDave May   PetscInt p,npoints,*rankval,n_points_recv;
1002712d1f2SDave May   PetscMPIInt rank,nrank,negrank;
1012712d1f2SDave May   void *point_buffer,*recv_points;
1022712d1f2SDave May   size_t sizeof_dmswarm_point;
1032712d1f2SDave May 
1042712d1f2SDave May   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
1052712d1f2SDave May 
1062712d1f2SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
1072712d1f2SDave May   *globalsize = npoints;
1082712d1f2SDave May   ierr = DMSwarmGetField(dm,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
1092712d1f2SDave May 
1102712d1f2SDave May   de = DataExCreate(PetscObjectComm((PetscObject)dm),0);CHKERRQ(ierr);
1112712d1f2SDave May 
1122712d1f2SDave May   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
1132712d1f2SDave May   for (p=0; p<npoints; p++) {
1142712d1f2SDave May     negrank = rankval[p];
1152712d1f2SDave May     if (negrank < 0) {
1162712d1f2SDave May       nrank = -negrank - 1;
1172712d1f2SDave May       ierr = DataExTopologyAddNeighbour(de,nrank);CHKERRQ(ierr);
1182712d1f2SDave May     }
1192712d1f2SDave May   }
1202712d1f2SDave May   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
1212712d1f2SDave May 
1222712d1f2SDave May   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
1232712d1f2SDave May   for (p=0; p<npoints; p++) {
1242712d1f2SDave May     negrank = rankval[p];
1252712d1f2SDave May     if (negrank < 0) {
1262712d1f2SDave May       nrank = -negrank - 1;
1272712d1f2SDave May       ierr = DataExAddToSendCount(de,nrank,1);CHKERRQ(ierr);
1282712d1f2SDave May     }
1292712d1f2SDave May   }
1302712d1f2SDave May   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
1312712d1f2SDave May 
1322712d1f2SDave May   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
1332712d1f2SDave May   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
1342712d1f2SDave May   for (p=0; p<npoints; p++) {
1352712d1f2SDave May     negrank = rankval[p];
1362712d1f2SDave May     if (negrank < 0) {
1372712d1f2SDave May       nrank = -negrank - 1;
1382712d1f2SDave May       rankval[p] = nrank;
1392712d1f2SDave May       /* copy point into buffer */
1402712d1f2SDave May       ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr);
1412712d1f2SDave May       /* insert point buffer into DataExchanger */
1422712d1f2SDave May       ierr = DataExPackData(de,nrank,1,point_buffer);CHKERRQ(ierr);
1432712d1f2SDave May       rankval[p] = negrank;
1442712d1f2SDave May     }
1452712d1f2SDave May   }
1462712d1f2SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
1472712d1f2SDave May 
1482712d1f2SDave May   ierr = DMSwarmRestoreField(dm,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
1492712d1f2SDave May 
1502712d1f2SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
1512712d1f2SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
1522712d1f2SDave May 
1532712d1f2SDave May   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
1542712d1f2SDave May 
1552712d1f2SDave May 	ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
1562712d1f2SDave May 	ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr);
1572712d1f2SDave May 	for (p=0; p<n_points_recv; p++) {
1582712d1f2SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
1592712d1f2SDave May 
1602712d1f2SDave May     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
1612712d1f2SDave May   }
1622712d1f2SDave May   ierr = DataExView(de);CHKERRQ(ierr);
1632712d1f2SDave May   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
1642712d1f2SDave May   ierr = DataExDestroy(de);CHKERRQ(ierr);
1652712d1f2SDave May 
1662712d1f2SDave May   PetscFunctionReturn(0);
1672712d1f2SDave May }
168*b16650c8SDave May 
169*b16650c8SDave May typedef struct {
170*b16650c8SDave May   PetscMPIInt owner_rank;
171*b16650c8SDave May   PetscReal min[3],max[3];
172*b16650c8SDave May } CollectBBox;
173*b16650c8SDave May 
174*b16650c8SDave May #undef __FUNCT__
175*b16650c8SDave May #define __FUNCT__ "DMSwarmMigrate_GlobalToLocal_BoundingBox"
176*b16650c8SDave May PETSC_EXTERN PetscErrorCode DMSwarmMigrate_GlobalToLocal_BoundingBox(DM dm,PetscInt *globalsize)
177*b16650c8SDave May {
178*b16650c8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
179*b16650c8SDave May   PetscErrorCode ierr;
180*b16650c8SDave May   DataEx de;
181*b16650c8SDave May   PetscInt p,pk,npoints,*rankval,n_points_recv,n_bbox_recv,dim,neighbour_cells;
182*b16650c8SDave May   PetscMPIInt rank,nrank;
183*b16650c8SDave May   void *point_buffer,*recv_points;
184*b16650c8SDave May   size_t sizeof_dmswarm_point,sizeof_bbox_ctx;
185*b16650c8SDave May   PetscBool isdmda;
186*b16650c8SDave May   CollectBBox *bbox,*recv_bbox;
187*b16650c8SDave May   const PetscMPIInt *dmneighborranks;
188*b16650c8SDave May   DM dmcell;
189*b16650c8SDave May 
190*b16650c8SDave May   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
191*b16650c8SDave May 
192*b16650c8SDave May   //ierr = DMSwarmGetDMCell(swarm,&dmcell);CHKERRQ(ierr);
193*b16650c8SDave May   dmcell = swarm->dmcell;
194*b16650c8SDave May 
195*b16650c8SDave May   dim = 2;
196*b16650c8SDave May   if (dim == 2) {
197*b16650c8SDave May     neighbour_cells = 9;
198*b16650c8SDave May   }
199*b16650c8SDave May 
200*b16650c8SDave May   isdmda = PETSC_FALSE;
201*b16650c8SDave May   PetscObjectTypeCompare((PetscObject)dmcell,DMDA,&isdmda);
202*b16650c8SDave May   if (!isdmda) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only DMDA support for CollectBoundingBox");
203*b16650c8SDave May 
204*b16650c8SDave May   sizeof_bbox_ctx = sizeof(CollectBBox);
205*b16650c8SDave May   PetscMalloc1(1,&bbox);
206*b16650c8SDave May   bbox->owner_rank = rank;
207*b16650c8SDave May 
208*b16650c8SDave May   /* compute the bounding box based on the overlapping / stenctil size */
209*b16650c8SDave May   //ierr = DMDAGetLocalBoundingBox(dmcell,bbox->min,bbox->max);CHKERRQ(ierr);
210*b16650c8SDave May   {
211*b16650c8SDave May     Vec lcoor;
212*b16650c8SDave May 
213*b16650c8SDave May     ierr = DMGetCoordinatesLocal(dmcell,&lcoor);CHKERRQ(ierr);
214*b16650c8SDave May 
215*b16650c8SDave May     ierr = VecStrideMin(lcoor,0,NULL,&bbox->min[0]);CHKERRQ(ierr);
216*b16650c8SDave May     ierr = VecStrideMin(lcoor,1,NULL,&bbox->min[1]);CHKERRQ(ierr);
217*b16650c8SDave May     ierr = VecStrideMax(lcoor,0,NULL,&bbox->max[0]);CHKERRQ(ierr);
218*b16650c8SDave May     ierr = VecStrideMax(lcoor,1,NULL,&bbox->max[1]);CHKERRQ(ierr);
219*b16650c8SDave May   }
220*b16650c8SDave May 
221*b16650c8SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
222*b16650c8SDave May   *globalsize = npoints;
223*b16650c8SDave May   ierr = DMSwarmGetField(dm,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
224*b16650c8SDave May 
225*b16650c8SDave May   de = DataExCreate(PetscObjectComm((PetscObject)dm),0);CHKERRQ(ierr);
226*b16650c8SDave May 
227*b16650c8SDave May   /* use DMDA neighbours */
228*b16650c8SDave May 	ierr = DMDAGetNeighbors(dmcell,&dmneighborranks);CHKERRQ(ierr);
229*b16650c8SDave May 
230*b16650c8SDave May   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
231*b16650c8SDave May   for (p=0; p<neighbour_cells; p++) {
232*b16650c8SDave May 		if ( (dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank) ) {
233*b16650c8SDave May       ierr = DataExTopologyAddNeighbour(de,dmneighborranks[p]);CHKERRQ(ierr);
234*b16650c8SDave May     }
235*b16650c8SDave May   }
236*b16650c8SDave May   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
237*b16650c8SDave May 
238*b16650c8SDave May   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
239*b16650c8SDave May   for (p=0; p<neighbour_cells; p++) {
240*b16650c8SDave May 		if ( (dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank) ) {
241*b16650c8SDave May       ierr = DataExAddToSendCount(de,dmneighborranks[p],1);CHKERRQ(ierr);
242*b16650c8SDave May     }
243*b16650c8SDave May   }
244*b16650c8SDave May   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
245*b16650c8SDave May 
246*b16650c8SDave May 
247*b16650c8SDave May   /* send bounding boxes */
248*b16650c8SDave May   ierr = DataExPackInitialize(de,sizeof_bbox_ctx);CHKERRQ(ierr);
249*b16650c8SDave May   for (p=0; p<neighbour_cells; p++) {
250*b16650c8SDave May     nrank = dmneighborranks[p];
251*b16650c8SDave May 		if ( (nrank >= 0) && (nrank != rank) ) {
252*b16650c8SDave May       /* insert bbox buffer into DataExchanger */
253*b16650c8SDave May       ierr = DataExPackData(de,nrank,1,bbox);CHKERRQ(ierr);
254*b16650c8SDave May     }
255*b16650c8SDave May   }
256*b16650c8SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
257*b16650c8SDave May 
258*b16650c8SDave May   /* recv bounding boxes */
259*b16650c8SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
260*b16650c8SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
261*b16650c8SDave May 
262*b16650c8SDave May   ierr = DataExGetRecvData(de,&n_bbox_recv,(void**)&recv_bbox);CHKERRQ(ierr);
263*b16650c8SDave May 
264*b16650c8SDave May 
265*b16650c8SDave May   for (p=0; p<n_bbox_recv; p++) {
266*b16650c8SDave May     printf("[rank %d]: box from %d : range[%+1.4e,%+1.4e]x[%+1.4e,%+1.4e]\n",rank,recv_bbox[p].owner_rank,
267*b16650c8SDave May            recv_bbox[p].min[0],recv_bbox[p].max[0],recv_bbox[p].min[1],recv_bbox[p].max[1]);
268*b16650c8SDave May   }
269*b16650c8SDave May 
270*b16650c8SDave May   /* of course this is stupid as this "generic" function should have a better way to know what the coordinates are called */
271*b16650c8SDave May   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
272*b16650c8SDave May   for (pk=0; pk<n_bbox_recv; pk++) {
273*b16650c8SDave May     PetscReal *array_x,*array_y;
274*b16650c8SDave May 
275*b16650c8SDave May     ierr = DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
276*b16650c8SDave May     ierr = DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
277*b16650c8SDave May 
278*b16650c8SDave May     for (p=0; p<npoints; p++) {
279*b16650c8SDave May       if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0]) ) {
280*b16650c8SDave May         if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1]) ) {
281*b16650c8SDave May 
282*b16650c8SDave May           ierr = DataExAddToSendCount(de,recv_bbox[pk].owner_rank,1);CHKERRQ(ierr);
283*b16650c8SDave May 
284*b16650c8SDave May         }
285*b16650c8SDave May       }
286*b16650c8SDave May     }
287*b16650c8SDave May     ierr = DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
288*b16650c8SDave May     ierr = DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
289*b16650c8SDave May   }
290*b16650c8SDave May   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
291*b16650c8SDave May 
292*b16650c8SDave May 
293*b16650c8SDave May 
294*b16650c8SDave May   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
295*b16650c8SDave May   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
296*b16650c8SDave May   for (pk=0; pk<n_bbox_recv; pk++) {
297*b16650c8SDave May     PetscReal *array_x,*array_y;
298*b16650c8SDave May 
299*b16650c8SDave May     ierr = DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
300*b16650c8SDave May     ierr = DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
301*b16650c8SDave May 
302*b16650c8SDave May     for (p=0; p<npoints; p++) {
303*b16650c8SDave May       if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0]) ) {
304*b16650c8SDave May         if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1]) ) {
305*b16650c8SDave May           // copy point into buffer //
306*b16650c8SDave May           ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr);
307*b16650c8SDave May           // insert point buffer into DataExchanger //
308*b16650c8SDave May           ierr = DataExPackData(de,recv_bbox[pk].owner_rank,1,point_buffer);CHKERRQ(ierr);
309*b16650c8SDave May         }
310*b16650c8SDave May       }
311*b16650c8SDave May     }
312*b16650c8SDave May 
313*b16650c8SDave May     ierr = DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
314*b16650c8SDave May     ierr = DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
315*b16650c8SDave May   }
316*b16650c8SDave May 
317*b16650c8SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
318*b16650c8SDave May 
319*b16650c8SDave May   ierr = DMSwarmRestoreField(dm,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
320*b16650c8SDave May 
321*b16650c8SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
322*b16650c8SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
323*b16650c8SDave May 
324*b16650c8SDave May   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
325*b16650c8SDave May 
326*b16650c8SDave May 	ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
327*b16650c8SDave May 	ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr);
328*b16650c8SDave May 	for (p=0; p<n_points_recv; p++) {
329*b16650c8SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
330*b16650c8SDave May 
331*b16650c8SDave May     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
332*b16650c8SDave May   }
333*b16650c8SDave May 
334*b16650c8SDave May   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
335*b16650c8SDave May   PetscFree(bbox);
336*b16650c8SDave May   ierr = DataExView(de);CHKERRQ(ierr);
337*b16650c8SDave May   ierr = DataExDestroy(de);CHKERRQ(ierr);
338*b16650c8SDave May 
339*b16650c8SDave May   PetscFunctionReturn(0);
340*b16650c8SDave May }
341