xref: /petsc/src/dm/impls/swarm/swarm_migrate.c (revision 08056efc59261c75231b7f4766f88349a1eb8f6d)
1df21e3a8SDave May 
2*08056efcSDave May #include <petscdmswarm.h>
3df21e3a8SDave May #include <petsc/private/dmswarmimpl.h>    /*I   "petscdmswarm.h"   I*/
4df21e3a8SDave May #include "data_bucket.h"
5df21e3a8SDave May #include "data_ex.h"
6df21e3a8SDave May 
7df21e3a8SDave May 
8480eef7bSDave May /*
9480eef7bSDave May  User loads desired location (MPI rank) into field DMSwarm_rank
10480eef7bSDave May */
11df21e3a8SDave May #undef __FUNCT__
12df21e3a8SDave May #define __FUNCT__ "DMSwarmMigrate_Push_Basic"
13df21e3a8SDave May PetscErrorCode DMSwarmMigrate_Push_Basic(DM dm,PetscBool remove_sent_points)
14df21e3a8SDave May {
15df21e3a8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
16df21e3a8SDave May   PetscErrorCode ierr;
17df21e3a8SDave May   DataEx de;
18df21e3a8SDave May   PetscInt p,npoints,*rankval,n_points_recv;
19df21e3a8SDave May   PetscMPIInt rank,nrank;
20df21e3a8SDave May   void *point_buffer,*recv_points;
21df21e3a8SDave May   size_t sizeof_dmswarm_point;
22df21e3a8SDave May 
23df21e3a8SDave May   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
24df21e3a8SDave May 
25df21e3a8SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
26*08056efcSDave May   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
27df21e3a8SDave May 
28df21e3a8SDave May   de = DataExCreate(PetscObjectComm((PetscObject)dm),0);CHKERRQ(ierr);
29df21e3a8SDave May 
30df21e3a8SDave May   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
31df21e3a8SDave May   for (p=0; p<npoints; p++) {
32df21e3a8SDave May     nrank = rankval[p];
33df21e3a8SDave May     if (nrank != rank) {
34df21e3a8SDave May       ierr = DataExTopologyAddNeighbour(de,nrank);CHKERRQ(ierr);
35df21e3a8SDave May     }
36df21e3a8SDave May   }
37df21e3a8SDave May   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
38df21e3a8SDave May 
39df21e3a8SDave May   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
40df21e3a8SDave May   for (p=0; p<npoints; p++) {
41df21e3a8SDave May     nrank = rankval[p];
42df21e3a8SDave May     if (nrank != rank) {
43df21e3a8SDave May       ierr = DataExAddToSendCount(de,nrank,1);CHKERRQ(ierr);
44df21e3a8SDave May     }
45df21e3a8SDave May   }
46df21e3a8SDave May   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
47df21e3a8SDave May 
48df21e3a8SDave May   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
49df21e3a8SDave May   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
50df21e3a8SDave May   for (p=0; p<npoints; p++) {
51df21e3a8SDave May     nrank = rankval[p];
52df21e3a8SDave May     if (nrank != rank) {
53df21e3a8SDave May       /* copy point into buffer */
54df21e3a8SDave May       ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr);
55df21e3a8SDave May       /* insert point buffer into DataExchanger */
56df21e3a8SDave May       ierr = DataExPackData(de,nrank,1,point_buffer);CHKERRQ(ierr);
57df21e3a8SDave May     }
58df21e3a8SDave May   }
59df21e3a8SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
60df21e3a8SDave May 
61df21e3a8SDave May 
62df21e3a8SDave May   if (remove_sent_points) {
63df21e3a8SDave May     /* remove points which left processor */
64df21e3a8SDave May     ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
65df21e3a8SDave May     for (p=0; p<npoints; p++) {
66df21e3a8SDave May       nrank = rankval[p];
67df21e3a8SDave May       if (nrank != rank) {
68df21e3a8SDave May         /* kill point */
69df21e3a8SDave May         ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr);
70df21e3a8SDave May         DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */
71df21e3a8SDave May         p--; /* check replacement point */
72df21e3a8SDave May       }
73df21e3a8SDave May     }
74df21e3a8SDave May   }
75*08056efcSDave May   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
76df21e3a8SDave May 
77df21e3a8SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
78df21e3a8SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
79df21e3a8SDave May 
80df21e3a8SDave May   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
81df21e3a8SDave May 
82df21e3a8SDave May 	ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
83df21e3a8SDave May 	ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr);
84df21e3a8SDave May 	for (p=0; p<n_points_recv; p++) {
85df21e3a8SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
86df21e3a8SDave May 
87df21e3a8SDave May     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
88df21e3a8SDave May   }
89df21e3a8SDave May   ierr = DataExView(de);CHKERRQ(ierr);
90df21e3a8SDave May   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
91df21e3a8SDave May   ierr = DataExDestroy(de);CHKERRQ(ierr);
92df21e3a8SDave May 
93df21e3a8SDave May   PetscFunctionReturn(0);
94df21e3a8SDave May }
952712d1f2SDave May 
9640c453e9SDave May #undef __FUNCT__
97889dbfe5SDave May #define __FUNCT__ "DMSwarmMigrate_DMNeighborScatter"
98889dbfe5SDave May PetscErrorCode DMSwarmMigrate_DMNeighborScatter(DM dm,DM dmcell,PetscBool remove_sent_points,PetscInt *npoints_prior_migration)
9940c453e9SDave May {
10040c453e9SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
10140c453e9SDave May   PetscErrorCode ierr;
10240c453e9SDave May   DataEx de;
10340c453e9SDave May   PetscInt r,p,npoints,*rankval,n_points_recv;
10440c453e9SDave May   PetscMPIInt rank,_rank;
10540c453e9SDave May   const PetscMPIInt *neighbourranks;
10640c453e9SDave May   void *point_buffer,*recv_points;
10740c453e9SDave May   size_t sizeof_dmswarm_point;
10840c453e9SDave May   PetscInt nneighbors;
10940c453e9SDave May 
11040c453e9SDave May   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
11140c453e9SDave May 
11240c453e9SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
113*08056efcSDave May   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
11440c453e9SDave May 
11540c453e9SDave May   de = DataExCreate(PetscObjectComm((PetscObject)dm),0);CHKERRQ(ierr);
11640c453e9SDave May 
11740c453e9SDave May   ierr = DMGetNeighbors(dmcell,&nneighbors,&neighbourranks);CHKERRQ(ierr);
11840c453e9SDave May   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
11940c453e9SDave May   for (r=0; r<nneighbors; r++) {
12040c453e9SDave May     _rank = neighbourranks[r];
12140c453e9SDave May     if ((_rank != rank) && (_rank > 0)) {
12240c453e9SDave May       ierr = DataExTopologyAddNeighbour(de,_rank);CHKERRQ(ierr);
12340c453e9SDave May     }
12440c453e9SDave May   }
12540c453e9SDave May   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
12640c453e9SDave May 
12740c453e9SDave May   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
12840c453e9SDave May   for (p=0; p<npoints; p++) {
12940c453e9SDave May     if (rankval[p] == -1) {
13040c453e9SDave May       for (r=0; r<nneighbors; r++) {
13140c453e9SDave May         _rank = neighbourranks[r];
13240c453e9SDave May         if ((_rank != rank) && (_rank > 0)) {
13340c453e9SDave May           ierr = DataExAddToSendCount(de,_rank,1);CHKERRQ(ierr);
13440c453e9SDave May         }
13540c453e9SDave May       }
13640c453e9SDave May     }
13740c453e9SDave May   }
13840c453e9SDave May   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
13940c453e9SDave May 
14040c453e9SDave May   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
14140c453e9SDave May   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
14240c453e9SDave May   for (p=0; p<npoints; p++) {
14340c453e9SDave May     if (rankval[p] == -1) {
14440c453e9SDave May       for (r=0; r<nneighbors; r++) {
14540c453e9SDave May         _rank = neighbourranks[r];
14640c453e9SDave May         if ((_rank != rank) && (_rank > 0)) {
14740c453e9SDave May           /* copy point into buffer */
14840c453e9SDave May           ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr);
14940c453e9SDave May           /* insert point buffer into DataExchanger */
15040c453e9SDave May           ierr = DataExPackData(de,_rank,1,point_buffer);CHKERRQ(ierr);
15140c453e9SDave May         }
15240c453e9SDave May       }
15340c453e9SDave May     }
15440c453e9SDave May   }
15540c453e9SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
15640c453e9SDave May 
15740c453e9SDave May 
15840c453e9SDave May   if (remove_sent_points) {
15940c453e9SDave May     /* remove points which left processor */
16040c453e9SDave May     ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
16140c453e9SDave May     for (p=0; p<npoints; p++) {
16240c453e9SDave May       if (rankval[p] == -1) {
16340c453e9SDave May         /* kill point */
16440c453e9SDave May         ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr);
16540c453e9SDave May         DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */
16640c453e9SDave May         p--; /* check replacement point */
16740c453e9SDave May       }
16840c453e9SDave May     }
16940c453e9SDave May   }
170*08056efcSDave May   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
17140c453e9SDave May   ierr = DataBucketGetSizes(swarm->db,npoints_prior_migration,NULL,NULL);CHKERRQ(ierr);
17240c453e9SDave May 
17340c453e9SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
17440c453e9SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
17540c453e9SDave May 
17640c453e9SDave May   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
17740c453e9SDave May 
17840c453e9SDave May 	ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
17940c453e9SDave May 	ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr);
18040c453e9SDave May 	for (p=0; p<n_points_recv; p++) {
18140c453e9SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
18240c453e9SDave May 
18340c453e9SDave May     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
18440c453e9SDave May   }
18540c453e9SDave May   ierr = DataExView(de);CHKERRQ(ierr);
18640c453e9SDave May   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
18740c453e9SDave May   ierr = DataExDestroy(de);CHKERRQ(ierr);
18840c453e9SDave May 
18940c453e9SDave May   PetscFunctionReturn(0);
19040c453e9SDave May }
191480eef7bSDave May 
192480eef7bSDave May #undef __FUNCT__
193*08056efcSDave May #define __FUNCT__ "DMSwarmMigrate_CellDMScatter"
194*08056efcSDave May PetscErrorCode DMSwarmMigrate_CellDMScatter(DM dm,PetscBool remove_sent_points)
195480eef7bSDave May {
196480eef7bSDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
197480eef7bSDave May   PetscErrorCode ierr;
19840c453e9SDave May   PetscInt p,npoints,npointsg,npoints2,npoints2g,len,*rankval,npoints_prior_migration;
199480eef7bSDave May   const PetscInt *LA_iscell;
200480eef7bSDave May   DM dmcell;
201480eef7bSDave May   IS iscell;
202480eef7bSDave May   Vec pos;
20340c453e9SDave May   PetscBool error_check = swarm->migrate_error_on_missing_point;
204480eef7bSDave May 
205480eef7bSDave May   ierr = DMSwarmGetCellDM(dm,&dmcell);CHKERRQ(ierr);
206480eef7bSDave May   if (!dmcell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid if cell DM provided");
207480eef7bSDave May 
208*08056efcSDave May   ierr = DMSwarmCreateGlobalVectorFromField(dm,DMSwarmPICField_coor,&pos);CHKERRQ(ierr);
209480eef7bSDave May   ierr = DMLocatePoints(dmcell,pos,&iscell);CHKERRQ(ierr);
210*08056efcSDave May   ierr = DMSwarmDestroyGlobalVectorFromField(dm,DMSwarmPICField_coor,&pos);CHKERRQ(ierr);
211480eef7bSDave May 
21240c453e9SDave May   if (error_check) {
21340c453e9SDave May     ierr = DMSwarmGetSize(dm,&npointsg);CHKERRQ(ierr);
21440c453e9SDave May   }
215480eef7bSDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
216*08056efcSDave May   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
217480eef7bSDave May   ierr = ISGetIndices(iscell,&LA_iscell);CHKERRQ(ierr);
218480eef7bSDave May   for (p=0; p<npoints; p++) {
219480eef7bSDave May     if (LA_iscell[p] == -1) {
220480eef7bSDave May       rankval[p] = -1;
221480eef7bSDave May     }
222480eef7bSDave May   }
223480eef7bSDave May   ierr = ISRestoreIndices(iscell,&LA_iscell);CHKERRQ(ierr);
224*08056efcSDave May   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
225480eef7bSDave May   ierr = ISDestroy(&iscell);CHKERRQ(ierr);
226480eef7bSDave May 
227889dbfe5SDave May   ierr = DMSwarmMigrate_DMNeighborScatter(dm,dmcell,remove_sent_points,&npoints_prior_migration);CHKERRQ(ierr);
228480eef7bSDave May 
22940c453e9SDave May   /* locate points newly recevied */
23040c453e9SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr);
23140c453e9SDave May   len = npoints2 - npoints_prior_migration;
23240c453e9SDave May   if (len > 0) {
23340c453e9SDave May     PetscScalar *LA_coor;
23440c453e9SDave May     PetscInt bs = 1;
23540c453e9SDave May 
236*08056efcSDave May     ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&LA_coor);CHKERRQ(ierr);
23740c453e9SDave May     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,len,(const PetscScalar*)&LA_coor[npoints_prior_migration],&pos);CHKERRQ(ierr);
23840c453e9SDave May 
23940c453e9SDave May     ierr = DMLocatePoints(dmcell,pos,&iscell);CHKERRQ(ierr);
24040c453e9SDave May 
24140c453e9SDave May     ierr = VecDestroy(&pos);CHKERRQ(ierr);
242*08056efcSDave May     ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&LA_coor);CHKERRQ(ierr);
24340c453e9SDave May 
24440c453e9SDave May     ierr = ISGetIndices(iscell,&LA_iscell);CHKERRQ(ierr);
245*08056efcSDave May     ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
24640c453e9SDave May     for (p=0; p<len; p++) {
24740c453e9SDave May       if (LA_iscell[p] == -1) {
24840c453e9SDave May         rankval[npoints_prior_migration+p] = -1;
24940c453e9SDave May       }
25040c453e9SDave May     }
25140c453e9SDave May     ierr = ISRestoreIndices(iscell,&LA_iscell);CHKERRQ(ierr);
25240c453e9SDave May     ierr = ISDestroy(&iscell);CHKERRQ(ierr);
25340c453e9SDave May 
25440c453e9SDave May     /* remove points which left processor */
25540c453e9SDave May     ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
25640c453e9SDave May     for (p=npoints_prior_migration; p<npoints2; p++) {
25740c453e9SDave May       if (rankval[p] == -1) {
25840c453e9SDave May         /* kill point */
25940c453e9SDave May         ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr);
26040c453e9SDave May         DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */
26140c453e9SDave May         p--; /* check replacement point */
26240c453e9SDave May       }
26340c453e9SDave May     }
264*08056efcSDave May     ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
26540c453e9SDave May 
26640c453e9SDave May   }
26740c453e9SDave May 
26840c453e9SDave May   /* check for error on removed points */
26940c453e9SDave May   if (error_check) {
27040c453e9SDave May     ierr = DMSwarmGetSize(dm,&npoints2g);CHKERRQ(ierr);
27140c453e9SDave May     if (npointsg != npoints2g) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Points from the DMSwarm must remain constant during migration (initial %D - final %D)",npointsg,npoints2g);
27240c453e9SDave May   }
273480eef7bSDave May   PetscFunctionReturn(0);
274480eef7bSDave May }
275480eef7bSDave May 
276*08056efcSDave May #undef __FUNCT__
277*08056efcSDave May #define __FUNCT__ "DMSwarmMigrate_CellDMExact"
278*08056efcSDave May PetscErrorCode DMSwarmMigrate_CellDMExact(DM dm,PetscBool remove_sent_points)
279*08056efcSDave May {
280*08056efcSDave May //  DM_Swarm *swarm = (DM_Swarm*)dm->data;
281*08056efcSDave May //  PetscErrorCode ierr;
282*08056efcSDave May 
283*08056efcSDave May   PetscFunctionReturn(0);
284*08056efcSDave May }
285*08056efcSDave May 
286480eef7bSDave May /*
287480eef7bSDave May  Redundant as this assumes points can only be sent to a single rank
288480eef7bSDave May */
2892712d1f2SDave May #undef __FUNCT__
2902712d1f2SDave May #define __FUNCT__ "DMSwarmMigrate_GlobalToLocal_Basic"
2912712d1f2SDave May PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm,PetscInt *globalsize)
2922712d1f2SDave May {
2932712d1f2SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
2942712d1f2SDave May   PetscErrorCode ierr;
2952712d1f2SDave May   DataEx de;
2962712d1f2SDave May   PetscInt p,npoints,*rankval,n_points_recv;
2972712d1f2SDave May   PetscMPIInt rank,nrank,negrank;
2982712d1f2SDave May   void *point_buffer,*recv_points;
2992712d1f2SDave May   size_t sizeof_dmswarm_point;
3002712d1f2SDave May 
3012712d1f2SDave May   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
3022712d1f2SDave May 
3032712d1f2SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
3042712d1f2SDave May   *globalsize = npoints;
305*08056efcSDave May   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
3062712d1f2SDave May 
3072712d1f2SDave May   de = DataExCreate(PetscObjectComm((PetscObject)dm),0);CHKERRQ(ierr);
3082712d1f2SDave May 
3092712d1f2SDave May   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
3102712d1f2SDave May   for (p=0; p<npoints; p++) {
3112712d1f2SDave May     negrank = rankval[p];
3122712d1f2SDave May     if (negrank < 0) {
3132712d1f2SDave May       nrank = -negrank - 1;
3142712d1f2SDave May       ierr = DataExTopologyAddNeighbour(de,nrank);CHKERRQ(ierr);
3152712d1f2SDave May     }
3162712d1f2SDave May   }
3172712d1f2SDave May   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
3182712d1f2SDave May 
3192712d1f2SDave May   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
3202712d1f2SDave May   for (p=0; p<npoints; p++) {
3212712d1f2SDave May     negrank = rankval[p];
3222712d1f2SDave May     if (negrank < 0) {
3232712d1f2SDave May       nrank = -negrank - 1;
3242712d1f2SDave May       ierr = DataExAddToSendCount(de,nrank,1);CHKERRQ(ierr);
3252712d1f2SDave May     }
3262712d1f2SDave May   }
3272712d1f2SDave May   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
3282712d1f2SDave May 
3292712d1f2SDave May   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
3302712d1f2SDave May   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
3312712d1f2SDave May   for (p=0; p<npoints; p++) {
3322712d1f2SDave May     negrank = rankval[p];
3332712d1f2SDave May     if (negrank < 0) {
3342712d1f2SDave May       nrank = -negrank - 1;
3352712d1f2SDave May       rankval[p] = nrank;
3362712d1f2SDave May       /* copy point into buffer */
3372712d1f2SDave May       ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr);
3382712d1f2SDave May       /* insert point buffer into DataExchanger */
3392712d1f2SDave May       ierr = DataExPackData(de,nrank,1,point_buffer);CHKERRQ(ierr);
3402712d1f2SDave May       rankval[p] = negrank;
3412712d1f2SDave May     }
3422712d1f2SDave May   }
3432712d1f2SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
3442712d1f2SDave May 
345*08056efcSDave May   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
3462712d1f2SDave May 
3472712d1f2SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
3482712d1f2SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
3492712d1f2SDave May 
3502712d1f2SDave May   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
3512712d1f2SDave May 
3522712d1f2SDave May 	ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
3532712d1f2SDave May 	ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr);
3542712d1f2SDave May 	for (p=0; p<n_points_recv; p++) {
3552712d1f2SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
3562712d1f2SDave May 
3572712d1f2SDave May     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
3582712d1f2SDave May   }
3592712d1f2SDave May   ierr = DataExView(de);CHKERRQ(ierr);
3602712d1f2SDave May   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
3612712d1f2SDave May   ierr = DataExDestroy(de);CHKERRQ(ierr);
3622712d1f2SDave May 
3632712d1f2SDave May   PetscFunctionReturn(0);
3642712d1f2SDave May }
365b16650c8SDave May 
366b16650c8SDave May typedef struct {
367b16650c8SDave May   PetscMPIInt owner_rank;
368b16650c8SDave May   PetscReal min[3],max[3];
369b16650c8SDave May } CollectBBox;
370b16650c8SDave May 
371b16650c8SDave May #undef __FUNCT__
372fe39f135SDave May #define __FUNCT__ "DMSwarmCollect_DMDABoundingBox"
373fe39f135SDave May PETSC_EXTERN PetscErrorCode DMSwarmCollect_DMDABoundingBox(DM dm,PetscInt *globalsize)
374b16650c8SDave May {
375b16650c8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
376b16650c8SDave May   PetscErrorCode ierr;
377b16650c8SDave May   DataEx de;
378b16650c8SDave May   PetscInt p,pk,npoints,*rankval,n_points_recv,n_bbox_recv,dim,neighbour_cells;
379b16650c8SDave May   PetscMPIInt rank,nrank;
380b16650c8SDave May   void *point_buffer,*recv_points;
381b16650c8SDave May   size_t sizeof_dmswarm_point,sizeof_bbox_ctx;
382b16650c8SDave May   PetscBool isdmda;
383b16650c8SDave May   CollectBBox *bbox,*recv_bbox;
384b16650c8SDave May   const PetscMPIInt *dmneighborranks;
385b16650c8SDave May   DM dmcell;
386b16650c8SDave May 
387b16650c8SDave May   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
388b16650c8SDave May 
389fe39f135SDave May   ierr = DMSwarmGetCellDM(dm,&dmcell);CHKERRQ(ierr);
390fe39f135SDave May   if (!dmcell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid if cell DM provided");
391b16650c8SDave May 
392b16650c8SDave May   isdmda = PETSC_FALSE;
393b16650c8SDave May   PetscObjectTypeCompare((PetscObject)dmcell,DMDA,&isdmda);
394b16650c8SDave May   if (!isdmda) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only DMDA support for CollectBoundingBox");
395b16650c8SDave May 
3968dbd68bcSDave May   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
397fe39f135SDave May 
398b16650c8SDave May   sizeof_bbox_ctx = sizeof(CollectBBox);
399b16650c8SDave May   PetscMalloc1(1,&bbox);
400b16650c8SDave May   bbox->owner_rank = rank;
401b16650c8SDave May 
402b16650c8SDave May   /* compute the bounding box based on the overlapping / stenctil size */
403b16650c8SDave May   //ierr = DMDAGetLocalBoundingBox(dmcell,bbox->min,bbox->max);CHKERRQ(ierr);
404b16650c8SDave May   {
405b16650c8SDave May     Vec lcoor;
406b16650c8SDave May 
407b16650c8SDave May     ierr = DMGetCoordinatesLocal(dmcell,&lcoor);CHKERRQ(ierr);
408b16650c8SDave May 
409fe39f135SDave May     if (dim >= 1) {
410b16650c8SDave May       ierr = VecStrideMin(lcoor,0,NULL,&bbox->min[0]);CHKERRQ(ierr);
411b16650c8SDave May       ierr = VecStrideMax(lcoor,0,NULL,&bbox->max[0]);CHKERRQ(ierr);
412fe39f135SDave May     }
413fe39f135SDave May 
414fe39f135SDave May     if (dim >= 2) {
415fe39f135SDave May       ierr = VecStrideMin(lcoor,1,NULL,&bbox->min[1]);CHKERRQ(ierr);
416b16650c8SDave May       ierr = VecStrideMax(lcoor,1,NULL,&bbox->max[1]);CHKERRQ(ierr);
417b16650c8SDave May     }
418b16650c8SDave May 
419fe39f135SDave May     if (dim == 3) {
420fe39f135SDave May       ierr = VecStrideMin(lcoor,2,NULL,&bbox->min[2]);CHKERRQ(ierr);
421fe39f135SDave May       ierr = VecStrideMax(lcoor,2,NULL,&bbox->max[2]);CHKERRQ(ierr);
422fe39f135SDave May     }
423fe39f135SDave May }
424fe39f135SDave May 
425b16650c8SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
426b16650c8SDave May   *globalsize = npoints;
427*08056efcSDave May   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
428b16650c8SDave May 
429b16650c8SDave May   de = DataExCreate(PetscObjectComm((PetscObject)dm),0);CHKERRQ(ierr);
430b16650c8SDave May 
431b16650c8SDave May   /* use DMDA neighbours */
432b16650c8SDave May 	ierr = DMDAGetNeighbors(dmcell,&dmneighborranks);CHKERRQ(ierr);
4338dbd68bcSDave May   if (dim == 1) {
4348dbd68bcSDave May     neighbour_cells = 3;
4358dbd68bcSDave May   } else if (dim == 2) {
4368dbd68bcSDave May     neighbour_cells = 9;
4378dbd68bcSDave May   } else {
4388dbd68bcSDave May     neighbour_cells = 27;
4398dbd68bcSDave May   }
440b16650c8SDave May 
441b16650c8SDave May   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
442b16650c8SDave May   for (p=0; p<neighbour_cells; p++) {
443b16650c8SDave May 		if ( (dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank) ) {
444b16650c8SDave May       ierr = DataExTopologyAddNeighbour(de,dmneighborranks[p]);CHKERRQ(ierr);
445b16650c8SDave May     }
446b16650c8SDave May   }
447b16650c8SDave May   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
448b16650c8SDave May 
449b16650c8SDave May   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
450b16650c8SDave May   for (p=0; p<neighbour_cells; p++) {
451b16650c8SDave May 		if ( (dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank) ) {
452b16650c8SDave May       ierr = DataExAddToSendCount(de,dmneighborranks[p],1);CHKERRQ(ierr);
453b16650c8SDave May     }
454b16650c8SDave May   }
455b16650c8SDave May   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
456b16650c8SDave May 
457b16650c8SDave May 
458b16650c8SDave May   /* send bounding boxes */
459b16650c8SDave May   ierr = DataExPackInitialize(de,sizeof_bbox_ctx);CHKERRQ(ierr);
460b16650c8SDave May   for (p=0; p<neighbour_cells; p++) {
461b16650c8SDave May     nrank = dmneighborranks[p];
462b16650c8SDave May 		if ( (nrank >= 0) && (nrank != rank) ) {
463b16650c8SDave May       /* insert bbox buffer into DataExchanger */
464b16650c8SDave May       ierr = DataExPackData(de,nrank,1,bbox);CHKERRQ(ierr);
465b16650c8SDave May     }
466b16650c8SDave May   }
467b16650c8SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
468b16650c8SDave May 
469b16650c8SDave May   /* recv bounding boxes */
470b16650c8SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
471b16650c8SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
472b16650c8SDave May 
473b16650c8SDave May   ierr = DataExGetRecvData(de,&n_bbox_recv,(void**)&recv_bbox);CHKERRQ(ierr);
474b16650c8SDave May 
475b16650c8SDave May 
476b16650c8SDave May   for (p=0; p<n_bbox_recv; p++) {
477b16650c8SDave May     printf("[rank %d]: box from %d : range[%+1.4e,%+1.4e]x[%+1.4e,%+1.4e]\n",rank,recv_bbox[p].owner_rank,
478b16650c8SDave May            recv_bbox[p].min[0],recv_bbox[p].max[0],recv_bbox[p].min[1],recv_bbox[p].max[1]);
479b16650c8SDave May   }
480b16650c8SDave May 
481b16650c8SDave May   /* of course this is stupid as this "generic" function should have a better way to know what the coordinates are called */
482b16650c8SDave May   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
483b16650c8SDave May   for (pk=0; pk<n_bbox_recv; pk++) {
484b16650c8SDave May     PetscReal *array_x,*array_y;
485b16650c8SDave May 
486b16650c8SDave May     ierr = DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
487b16650c8SDave May     ierr = DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
488b16650c8SDave May 
489b16650c8SDave May     for (p=0; p<npoints; p++) {
490b16650c8SDave May       if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0]) ) {
491b16650c8SDave May         if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1]) ) {
492b16650c8SDave May 
493b16650c8SDave May           ierr = DataExAddToSendCount(de,recv_bbox[pk].owner_rank,1);CHKERRQ(ierr);
494b16650c8SDave May 
495b16650c8SDave May         }
496b16650c8SDave May       }
497b16650c8SDave May     }
498b16650c8SDave May     ierr = DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
499b16650c8SDave May     ierr = DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
500b16650c8SDave May   }
501b16650c8SDave May   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
502b16650c8SDave May 
503b16650c8SDave May 
504b16650c8SDave May 
505b16650c8SDave May   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
506b16650c8SDave May   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
507b16650c8SDave May   for (pk=0; pk<n_bbox_recv; pk++) {
508b16650c8SDave May     PetscReal *array_x,*array_y;
509b16650c8SDave May 
510b16650c8SDave May     ierr = DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
511b16650c8SDave May     ierr = DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
512b16650c8SDave May 
513b16650c8SDave May     for (p=0; p<npoints; p++) {
514b16650c8SDave May       if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0]) ) {
515b16650c8SDave May         if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1]) ) {
516b16650c8SDave May           // copy point into buffer //
517b16650c8SDave May           ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr);
518b16650c8SDave May           // insert point buffer into DataExchanger //
519b16650c8SDave May           ierr = DataExPackData(de,recv_bbox[pk].owner_rank,1,point_buffer);CHKERRQ(ierr);
520b16650c8SDave May         }
521b16650c8SDave May       }
522b16650c8SDave May     }
523b16650c8SDave May 
524b16650c8SDave May     ierr = DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
525b16650c8SDave May     ierr = DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
526b16650c8SDave May   }
527b16650c8SDave May 
528b16650c8SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
529b16650c8SDave May 
530*08056efcSDave May   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
531b16650c8SDave May 
532b16650c8SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
533b16650c8SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
534b16650c8SDave May 
535b16650c8SDave May   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
536b16650c8SDave May 
537b16650c8SDave May 	ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
538b16650c8SDave May 	ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr);
539b16650c8SDave May 	for (p=0; p<n_points_recv; p++) {
540b16650c8SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
541b16650c8SDave May 
542b16650c8SDave May     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
543b16650c8SDave May   }
544b16650c8SDave May 
545b16650c8SDave May   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
546b16650c8SDave May   PetscFree(bbox);
547b16650c8SDave May   ierr = DataExView(de);CHKERRQ(ierr);
548b16650c8SDave May   ierr = DataExDestroy(de);CHKERRQ(ierr);
549b16650c8SDave May 
550b16650c8SDave May   PetscFunctionReturn(0);
551b16650c8SDave May }
552a9fd7477SDave May 
553a9fd7477SDave May 
554a9fd7477SDave May /* General collection when no order, or neighbour information is provided */
555a9fd7477SDave May /*
556a9fd7477SDave May  User provides context and collect() method
557a9fd7477SDave May  Broadcast user context
558a9fd7477SDave May 
559a9fd7477SDave May  for each context / rank {
560a9fd7477SDave May    collect(swarm,context,n,list)
561a9fd7477SDave May  }
562a9fd7477SDave May 
563a9fd7477SDave May */
564a9fd7477SDave May #undef __FUNCT__
565a9fd7477SDave May #define __FUNCT__ "DMSwarmCollect_General"
566a9fd7477SDave May PETSC_EXTERN PetscErrorCode DMSwarmCollect_General(DM dm,PetscErrorCode (*collect)(DM,void*,PetscInt*,PetscInt**),size_t ctx_size,void *ctx,PetscInt *globalsize)
567a9fd7477SDave May {
568a9fd7477SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
569a9fd7477SDave May   PetscErrorCode ierr;
570a9fd7477SDave May   DataEx         de;
571a9fd7477SDave May   PetscInt       p,r,npoints,n_points_recv;
572a9fd7477SDave May   PetscMPIInt    commsize,rank;
573a9fd7477SDave May   void           *point_buffer,*recv_points;
574a9fd7477SDave May   void           *ctxlist;
575a9fd7477SDave May   PetscInt       *n2collect,**collectlist;
576a9fd7477SDave May   size_t         sizeof_dmswarm_point;
577a9fd7477SDave May 
578a9fd7477SDave May   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&commsize);CHKERRQ(ierr);
579a9fd7477SDave May   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
580a9fd7477SDave May 
581a9fd7477SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
582a9fd7477SDave May   *globalsize = npoints;
583a9fd7477SDave May 
584a9fd7477SDave May   /* Broadcast user context */
585a9fd7477SDave May   PetscMalloc(ctx_size*commsize,&ctxlist);
586a9fd7477SDave May   ierr = MPI_Allgather(ctx,ctx_size,MPI_CHAR,ctxlist,ctx_size,MPI_CHAR,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
587a9fd7477SDave May 
588a9fd7477SDave May   PetscMalloc1(commsize,&n2collect);
589a9fd7477SDave May   PetscMalloc1(commsize,&collectlist);
590a9fd7477SDave May 
591a9fd7477SDave May   for (r=0; r<commsize; r++) {
592a9fd7477SDave May     PetscInt _n2collect;
593a9fd7477SDave May     PetscInt *_collectlist;
594a9fd7477SDave May     void     *_ctx_r;
595a9fd7477SDave May 
596a9fd7477SDave May     _n2collect   = 0;
597a9fd7477SDave May     _collectlist = NULL;
598a9fd7477SDave May 
599a9fd7477SDave May     if (r != rank) { /* don't collect data from yourself */
600a9fd7477SDave May       _ctx_r = (void*)( (char*)ctxlist + r * ctx_size );
601a9fd7477SDave May       ierr = collect(dm,_ctx_r,&_n2collect,&_collectlist);CHKERRQ(ierr);
602a9fd7477SDave May     }
603a9fd7477SDave May 
604a9fd7477SDave May     n2collect[r]   = _n2collect;
605a9fd7477SDave May     collectlist[r] = _collectlist;
606a9fd7477SDave May   }
607a9fd7477SDave May 
608a9fd7477SDave May   de = DataExCreate(PetscObjectComm((PetscObject)dm),0);CHKERRQ(ierr);
609a9fd7477SDave May 
610a9fd7477SDave May   /* Define topology */
611a9fd7477SDave May   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
612a9fd7477SDave May   for (r=0; r<commsize; r++) {
613a9fd7477SDave May     if (n2collect[r] > 0) {
614a9fd7477SDave May       ierr = DataExTopologyAddNeighbour(de,(PetscMPIInt)r);CHKERRQ(ierr);
615a9fd7477SDave May     }
616a9fd7477SDave May   }
617a9fd7477SDave May   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
618a9fd7477SDave May 
619a9fd7477SDave May   /* Define send counts */
620a9fd7477SDave May   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
621a9fd7477SDave May   for (r=0; r<commsize; r++) {
622a9fd7477SDave May     if (n2collect[r] > 0) {
623a9fd7477SDave May       ierr = DataExAddToSendCount(de,r,n2collect[r]);CHKERRQ(ierr);
624a9fd7477SDave May     }
625a9fd7477SDave May   }
626a9fd7477SDave May   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
627a9fd7477SDave May 
628a9fd7477SDave May   /* Pack data */
629a9fd7477SDave May   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
630a9fd7477SDave May 
631a9fd7477SDave May   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
632a9fd7477SDave May   for (r=0; r<commsize; r++) {
633a9fd7477SDave May 
634a9fd7477SDave May     for (p=0; p<n2collect[r]; p++) {
635a9fd7477SDave May       ierr = DataBucketFillPackedArray(swarm->db,collectlist[r][p],point_buffer);CHKERRQ(ierr);
636a9fd7477SDave May       /* insert point buffer into the data exchanger */
637a9fd7477SDave May       ierr = DataExPackData(de,r,1,point_buffer);CHKERRQ(ierr);
638a9fd7477SDave May     }
639a9fd7477SDave May   }
640a9fd7477SDave May 
641a9fd7477SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
642a9fd7477SDave May 
643a9fd7477SDave May   /* Scatter */
644a9fd7477SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
645a9fd7477SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
646a9fd7477SDave May 
647a9fd7477SDave May   /* Collect data in DMSwarm container */
648a9fd7477SDave May   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
649a9fd7477SDave May 
650a9fd7477SDave May 	ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
651a9fd7477SDave May 	ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr);
652a9fd7477SDave May 	for (p=0; p<n_points_recv; p++) {
653a9fd7477SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
654a9fd7477SDave May 
655a9fd7477SDave May     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
656a9fd7477SDave May   }
657a9fd7477SDave May 
658a9fd7477SDave May   /* Release memory */
659a9fd7477SDave May   for (r=0; r<commsize; r++) {
660a9fd7477SDave May     if (collectlist[r]) PetscFree(collectlist[r]);
661a9fd7477SDave May   }
662a9fd7477SDave May   PetscFree(collectlist);
663a9fd7477SDave May   PetscFree(n2collect);
664a9fd7477SDave May   PetscFree(ctxlist);
665a9fd7477SDave May   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
666a9fd7477SDave May   ierr = DataExView(de);CHKERRQ(ierr);
667a9fd7477SDave May   ierr = DataExDestroy(de);CHKERRQ(ierr);
668a9fd7477SDave May 
669a9fd7477SDave May   PetscFunctionReturn(0);
670a9fd7477SDave May }
671a9fd7477SDave May 
672