xref: /petsc/src/dm/impls/swarm/swarm_migrate.c (revision a9fd74774b7c05f3763235eff6df2ab103b4033b)
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 }
168b16650c8SDave May 
169b16650c8SDave May typedef struct {
170b16650c8SDave May   PetscMPIInt owner_rank;
171b16650c8SDave May   PetscReal min[3],max[3];
172b16650c8SDave May } CollectBBox;
173b16650c8SDave May 
174b16650c8SDave May #undef __FUNCT__
175b16650c8SDave May #define __FUNCT__ "DMSwarmMigrate_GlobalToLocal_BoundingBox"
176b16650c8SDave May PETSC_EXTERN PetscErrorCode DMSwarmMigrate_GlobalToLocal_BoundingBox(DM dm,PetscInt *globalsize)
177b16650c8SDave May {
178b16650c8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
179b16650c8SDave May   PetscErrorCode ierr;
180b16650c8SDave May   DataEx de;
181b16650c8SDave May   PetscInt p,pk,npoints,*rankval,n_points_recv,n_bbox_recv,dim,neighbour_cells;
182b16650c8SDave May   PetscMPIInt rank,nrank;
183b16650c8SDave May   void *point_buffer,*recv_points;
184b16650c8SDave May   size_t sizeof_dmswarm_point,sizeof_bbox_ctx;
185b16650c8SDave May   PetscBool isdmda;
186b16650c8SDave May   CollectBBox *bbox,*recv_bbox;
187b16650c8SDave May   const PetscMPIInt *dmneighborranks;
188b16650c8SDave May   DM dmcell;
189b16650c8SDave May 
190b16650c8SDave May   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
191b16650c8SDave May 
192b16650c8SDave May   //ierr = DMSwarmGetDMCell(swarm,&dmcell);CHKERRQ(ierr);
193b16650c8SDave May   dmcell = swarm->dmcell;
194b16650c8SDave May 
195b16650c8SDave May   dim = 2;
196b16650c8SDave May   if (dim == 2) {
197b16650c8SDave May     neighbour_cells = 9;
198b16650c8SDave May   }
199b16650c8SDave May 
200b16650c8SDave May   isdmda = PETSC_FALSE;
201b16650c8SDave May   PetscObjectTypeCompare((PetscObject)dmcell,DMDA,&isdmda);
202b16650c8SDave May   if (!isdmda) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only DMDA support for CollectBoundingBox");
203b16650c8SDave May 
204b16650c8SDave May   sizeof_bbox_ctx = sizeof(CollectBBox);
205b16650c8SDave May   PetscMalloc1(1,&bbox);
206b16650c8SDave May   bbox->owner_rank = rank;
207b16650c8SDave May 
208b16650c8SDave May   /* compute the bounding box based on the overlapping / stenctil size */
209b16650c8SDave May   //ierr = DMDAGetLocalBoundingBox(dmcell,bbox->min,bbox->max);CHKERRQ(ierr);
210b16650c8SDave May   {
211b16650c8SDave May     Vec lcoor;
212b16650c8SDave May 
213b16650c8SDave May     ierr = DMGetCoordinatesLocal(dmcell,&lcoor);CHKERRQ(ierr);
214b16650c8SDave May 
215b16650c8SDave May     ierr = VecStrideMin(lcoor,0,NULL,&bbox->min[0]);CHKERRQ(ierr);
216b16650c8SDave May     ierr = VecStrideMin(lcoor,1,NULL,&bbox->min[1]);CHKERRQ(ierr);
217b16650c8SDave May     ierr = VecStrideMax(lcoor,0,NULL,&bbox->max[0]);CHKERRQ(ierr);
218b16650c8SDave May     ierr = VecStrideMax(lcoor,1,NULL,&bbox->max[1]);CHKERRQ(ierr);
219b16650c8SDave May   }
220b16650c8SDave May 
221b16650c8SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
222b16650c8SDave May   *globalsize = npoints;
223b16650c8SDave May   ierr = DMSwarmGetField(dm,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
224b16650c8SDave May 
225b16650c8SDave May   de = DataExCreate(PetscObjectComm((PetscObject)dm),0);CHKERRQ(ierr);
226b16650c8SDave May 
227b16650c8SDave May   /* use DMDA neighbours */
228b16650c8SDave May 	ierr = DMDAGetNeighbors(dmcell,&dmneighborranks);CHKERRQ(ierr);
229b16650c8SDave May 
230b16650c8SDave May   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
231b16650c8SDave May   for (p=0; p<neighbour_cells; p++) {
232b16650c8SDave May 		if ( (dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank) ) {
233b16650c8SDave May       ierr = DataExTopologyAddNeighbour(de,dmneighborranks[p]);CHKERRQ(ierr);
234b16650c8SDave May     }
235b16650c8SDave May   }
236b16650c8SDave May   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
237b16650c8SDave May 
238b16650c8SDave May   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
239b16650c8SDave May   for (p=0; p<neighbour_cells; p++) {
240b16650c8SDave May 		if ( (dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank) ) {
241b16650c8SDave May       ierr = DataExAddToSendCount(de,dmneighborranks[p],1);CHKERRQ(ierr);
242b16650c8SDave May     }
243b16650c8SDave May   }
244b16650c8SDave May   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
245b16650c8SDave May 
246b16650c8SDave May 
247b16650c8SDave May   /* send bounding boxes */
248b16650c8SDave May   ierr = DataExPackInitialize(de,sizeof_bbox_ctx);CHKERRQ(ierr);
249b16650c8SDave May   for (p=0; p<neighbour_cells; p++) {
250b16650c8SDave May     nrank = dmneighborranks[p];
251b16650c8SDave May 		if ( (nrank >= 0) && (nrank != rank) ) {
252b16650c8SDave May       /* insert bbox buffer into DataExchanger */
253b16650c8SDave May       ierr = DataExPackData(de,nrank,1,bbox);CHKERRQ(ierr);
254b16650c8SDave May     }
255b16650c8SDave May   }
256b16650c8SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
257b16650c8SDave May 
258b16650c8SDave May   /* recv bounding boxes */
259b16650c8SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
260b16650c8SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
261b16650c8SDave May 
262b16650c8SDave May   ierr = DataExGetRecvData(de,&n_bbox_recv,(void**)&recv_bbox);CHKERRQ(ierr);
263b16650c8SDave May 
264b16650c8SDave May 
265b16650c8SDave May   for (p=0; p<n_bbox_recv; p++) {
266b16650c8SDave May     printf("[rank %d]: box from %d : range[%+1.4e,%+1.4e]x[%+1.4e,%+1.4e]\n",rank,recv_bbox[p].owner_rank,
267b16650c8SDave May            recv_bbox[p].min[0],recv_bbox[p].max[0],recv_bbox[p].min[1],recv_bbox[p].max[1]);
268b16650c8SDave May   }
269b16650c8SDave May 
270b16650c8SDave May   /* of course this is stupid as this "generic" function should have a better way to know what the coordinates are called */
271b16650c8SDave May   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
272b16650c8SDave May   for (pk=0; pk<n_bbox_recv; pk++) {
273b16650c8SDave May     PetscReal *array_x,*array_y;
274b16650c8SDave May 
275b16650c8SDave May     ierr = DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
276b16650c8SDave May     ierr = DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
277b16650c8SDave May 
278b16650c8SDave May     for (p=0; p<npoints; p++) {
279b16650c8SDave May       if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0]) ) {
280b16650c8SDave May         if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1]) ) {
281b16650c8SDave May 
282b16650c8SDave May           ierr = DataExAddToSendCount(de,recv_bbox[pk].owner_rank,1);CHKERRQ(ierr);
283b16650c8SDave May 
284b16650c8SDave May         }
285b16650c8SDave May       }
286b16650c8SDave May     }
287b16650c8SDave May     ierr = DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
288b16650c8SDave May     ierr = DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
289b16650c8SDave May   }
290b16650c8SDave May   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
291b16650c8SDave May 
292b16650c8SDave May 
293b16650c8SDave May 
294b16650c8SDave May   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
295b16650c8SDave May   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
296b16650c8SDave May   for (pk=0; pk<n_bbox_recv; pk++) {
297b16650c8SDave May     PetscReal *array_x,*array_y;
298b16650c8SDave May 
299b16650c8SDave May     ierr = DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
300b16650c8SDave May     ierr = DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
301b16650c8SDave May 
302b16650c8SDave May     for (p=0; p<npoints; p++) {
303b16650c8SDave May       if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0]) ) {
304b16650c8SDave May         if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1]) ) {
305b16650c8SDave May           // copy point into buffer //
306b16650c8SDave May           ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr);
307b16650c8SDave May           // insert point buffer into DataExchanger //
308b16650c8SDave May           ierr = DataExPackData(de,recv_bbox[pk].owner_rank,1,point_buffer);CHKERRQ(ierr);
309b16650c8SDave May         }
310b16650c8SDave May       }
311b16650c8SDave May     }
312b16650c8SDave May 
313b16650c8SDave May     ierr = DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
314b16650c8SDave May     ierr = DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
315b16650c8SDave May   }
316b16650c8SDave May 
317b16650c8SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
318b16650c8SDave May 
319b16650c8SDave May   ierr = DMSwarmRestoreField(dm,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
320b16650c8SDave May 
321b16650c8SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
322b16650c8SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
323b16650c8SDave May 
324b16650c8SDave May   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
325b16650c8SDave May 
326b16650c8SDave May 	ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
327b16650c8SDave May 	ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr);
328b16650c8SDave May 	for (p=0; p<n_points_recv; p++) {
329b16650c8SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
330b16650c8SDave May 
331b16650c8SDave May     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
332b16650c8SDave May   }
333b16650c8SDave May 
334b16650c8SDave May   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
335b16650c8SDave May   PetscFree(bbox);
336b16650c8SDave May   ierr = DataExView(de);CHKERRQ(ierr);
337b16650c8SDave May   ierr = DataExDestroy(de);CHKERRQ(ierr);
338b16650c8SDave May 
339b16650c8SDave May   PetscFunctionReturn(0);
340b16650c8SDave May }
341*a9fd7477SDave May 
342*a9fd7477SDave May 
343*a9fd7477SDave May /* General collection when no order, or neighbour information is provided */
344*a9fd7477SDave May /*
345*a9fd7477SDave May  User provides context and collect() method
346*a9fd7477SDave May  Broadcast user context
347*a9fd7477SDave May 
348*a9fd7477SDave May  for each context / rank {
349*a9fd7477SDave May    collect(swarm,context,n,list)
350*a9fd7477SDave May  }
351*a9fd7477SDave May 
352*a9fd7477SDave May */
353*a9fd7477SDave May #undef __FUNCT__
354*a9fd7477SDave May #define __FUNCT__ "DMSwarmCollect_General"
355*a9fd7477SDave May PETSC_EXTERN PetscErrorCode DMSwarmCollect_General(DM dm,PetscErrorCode (*collect)(DM,void*,PetscInt*,PetscInt**),size_t ctx_size,void *ctx,PetscInt *globalsize)
356*a9fd7477SDave May {
357*a9fd7477SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
358*a9fd7477SDave May   PetscErrorCode ierr;
359*a9fd7477SDave May   DataEx         de;
360*a9fd7477SDave May   PetscInt       p,r,npoints,n_points_recv;
361*a9fd7477SDave May   PetscMPIInt    commsize,rank;
362*a9fd7477SDave May   void           *point_buffer,*recv_points;
363*a9fd7477SDave May   void           *ctxlist;
364*a9fd7477SDave May   PetscInt       *n2collect,**collectlist;
365*a9fd7477SDave May   size_t         sizeof_dmswarm_point;
366*a9fd7477SDave May 
367*a9fd7477SDave May   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&commsize);CHKERRQ(ierr);
368*a9fd7477SDave May   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
369*a9fd7477SDave May 
370*a9fd7477SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
371*a9fd7477SDave May   *globalsize = npoints;
372*a9fd7477SDave May 
373*a9fd7477SDave May   /* Broadcast user context */
374*a9fd7477SDave May   PetscMalloc(ctx_size*commsize,&ctxlist);
375*a9fd7477SDave May   ierr = MPI_Allgather(ctx,ctx_size,MPI_CHAR,ctxlist,ctx_size,MPI_CHAR,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
376*a9fd7477SDave May 
377*a9fd7477SDave May   PetscMalloc1(commsize,&n2collect);
378*a9fd7477SDave May   PetscMalloc1(commsize,&collectlist);
379*a9fd7477SDave May 
380*a9fd7477SDave May   for (r=0; r<commsize; r++) {
381*a9fd7477SDave May     PetscInt _n2collect;
382*a9fd7477SDave May     PetscInt *_collectlist;
383*a9fd7477SDave May     void     *_ctx_r;
384*a9fd7477SDave May 
385*a9fd7477SDave May     _n2collect   = 0;
386*a9fd7477SDave May     _collectlist = NULL;
387*a9fd7477SDave May 
388*a9fd7477SDave May     if (r != rank) { /* don't collect data from yourself */
389*a9fd7477SDave May       _ctx_r = (void*)( (char*)ctxlist + r * ctx_size );
390*a9fd7477SDave May       ierr = collect(dm,_ctx_r,&_n2collect,&_collectlist);CHKERRQ(ierr);
391*a9fd7477SDave May     }
392*a9fd7477SDave May 
393*a9fd7477SDave May     n2collect[r]   = _n2collect;
394*a9fd7477SDave May     collectlist[r] = _collectlist;
395*a9fd7477SDave May   }
396*a9fd7477SDave May 
397*a9fd7477SDave May   de = DataExCreate(PetscObjectComm((PetscObject)dm),0);CHKERRQ(ierr);
398*a9fd7477SDave May 
399*a9fd7477SDave May   /* Define topology */
400*a9fd7477SDave May   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
401*a9fd7477SDave May   for (r=0; r<commsize; r++) {
402*a9fd7477SDave May     if (n2collect[r] > 0) {
403*a9fd7477SDave May       ierr = DataExTopologyAddNeighbour(de,(PetscMPIInt)r);CHKERRQ(ierr);
404*a9fd7477SDave May     }
405*a9fd7477SDave May   }
406*a9fd7477SDave May   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
407*a9fd7477SDave May 
408*a9fd7477SDave May   /* Define send counts */
409*a9fd7477SDave May   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
410*a9fd7477SDave May   for (r=0; r<commsize; r++) {
411*a9fd7477SDave May     if (n2collect[r] > 0) {
412*a9fd7477SDave May       ierr = DataExAddToSendCount(de,r,n2collect[r]);CHKERRQ(ierr);
413*a9fd7477SDave May     }
414*a9fd7477SDave May   }
415*a9fd7477SDave May   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
416*a9fd7477SDave May 
417*a9fd7477SDave May   /* Pack data */
418*a9fd7477SDave May   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
419*a9fd7477SDave May 
420*a9fd7477SDave May   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
421*a9fd7477SDave May   for (r=0; r<commsize; r++) {
422*a9fd7477SDave May 
423*a9fd7477SDave May     for (p=0; p<n2collect[r]; p++) {
424*a9fd7477SDave May       ierr = DataBucketFillPackedArray(swarm->db,collectlist[r][p],point_buffer);CHKERRQ(ierr);
425*a9fd7477SDave May       /* insert point buffer into the data exchanger */
426*a9fd7477SDave May       ierr = DataExPackData(de,r,1,point_buffer);CHKERRQ(ierr);
427*a9fd7477SDave May     }
428*a9fd7477SDave May   }
429*a9fd7477SDave May 
430*a9fd7477SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
431*a9fd7477SDave May 
432*a9fd7477SDave May   /* Scatter */
433*a9fd7477SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
434*a9fd7477SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
435*a9fd7477SDave May 
436*a9fd7477SDave May   /* Collect data in DMSwarm container */
437*a9fd7477SDave May   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
438*a9fd7477SDave May 
439*a9fd7477SDave May 	ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
440*a9fd7477SDave May 	ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr);
441*a9fd7477SDave May 	for (p=0; p<n_points_recv; p++) {
442*a9fd7477SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
443*a9fd7477SDave May 
444*a9fd7477SDave May     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
445*a9fd7477SDave May   }
446*a9fd7477SDave May 
447*a9fd7477SDave May   /* Release memory */
448*a9fd7477SDave May   for (r=0; r<commsize; r++) {
449*a9fd7477SDave May     if (collectlist[r]) PetscFree(collectlist[r]);
450*a9fd7477SDave May   }
451*a9fd7477SDave May   PetscFree(collectlist);
452*a9fd7477SDave May   PetscFree(n2collect);
453*a9fd7477SDave May   PetscFree(ctxlist);
454*a9fd7477SDave May   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
455*a9fd7477SDave May   ierr = DataExView(de);CHKERRQ(ierr);
456*a9fd7477SDave May   ierr = DataExDestroy(de);CHKERRQ(ierr);
457*a9fd7477SDave May 
458*a9fd7477SDave May   PetscFunctionReturn(0);
459*a9fd7477SDave May }
460*a9fd7477SDave May 
461