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