xref: /petsc/src/dm/impls/swarm/swarm_migrate.c (revision b4b02483d4b6e972a56a98fc1cecebc3b91c36d3)
1df21e3a8SDave May 
2dfc14de9SMatthew G. Knepley #include <petscsf.h>
308056efcSDave May #include <petscdmswarm.h>
4df21e3a8SDave May #include <petsc/private/dmswarmimpl.h>    /*I   "petscdmswarm.h"   I*/
5df21e3a8SDave May #include "data_bucket.h"
6df21e3a8SDave May #include "data_ex.h"
7df21e3a8SDave May 
8df21e3a8SDave May 
9480eef7bSDave May /*
10480eef7bSDave May  User loads desired location (MPI rank) into field DMSwarm_rank
11480eef7bSDave May */
12df21e3a8SDave May #undef __FUNCT__
13df21e3a8SDave May #define __FUNCT__ "DMSwarmMigrate_Push_Basic"
14df21e3a8SDave May PetscErrorCode DMSwarmMigrate_Push_Basic(DM dm,PetscBool remove_sent_points)
15df21e3a8SDave May {
16df21e3a8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
17df21e3a8SDave May   PetscErrorCode ierr;
18df21e3a8SDave May   DataEx de;
19df21e3a8SDave May   PetscInt p,npoints,*rankval,n_points_recv;
20df21e3a8SDave May   PetscMPIInt rank,nrank;
21df21e3a8SDave May   void *point_buffer,*recv_points;
22df21e3a8SDave May   size_t sizeof_dmswarm_point;
23df21e3a8SDave May 
24521f74f9SMatthew G. Knepley   PetscFunctionBegin;
25df21e3a8SDave May   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
26df21e3a8SDave May 
27df21e3a8SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
2808056efcSDave May   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
29521f74f9SMatthew G. Knepley   ierr = DataExCreate(PetscObjectComm((PetscObject)dm),0, &de);CHKERRQ(ierr);
30df21e3a8SDave May   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
31521f74f9SMatthew G. Knepley   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   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   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
47df21e3a8SDave May   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
48df21e3a8SDave May   for (p=0; p<npoints; p++) {
49df21e3a8SDave May     nrank = rankval[p];
50df21e3a8SDave May     if (nrank != rank) {
51df21e3a8SDave May       /* copy point into buffer */
52df21e3a8SDave May       ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr);
53df21e3a8SDave May       /* insert point buffer into DataExchanger */
54df21e3a8SDave May       ierr = DataExPackData(de,nrank,1,point_buffer);CHKERRQ(ierr);
55df21e3a8SDave May     }
56df21e3a8SDave May   }
57df21e3a8SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
58df21e3a8SDave May 
59df21e3a8SDave May   if (remove_sent_points) {
60df21e3a8SDave May     /* remove points which left processor */
61df21e3a8SDave May     ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
62df21e3a8SDave May     for (p=0; p<npoints; p++) {
63df21e3a8SDave May       nrank = rankval[p];
64df21e3a8SDave May       if (nrank != rank) {
65df21e3a8SDave May         /* kill point */
66df21e3a8SDave May         ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr);
67df21e3a8SDave May         DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */
68df21e3a8SDave May         p--; /* check replacement point */
69df21e3a8SDave May       }
70df21e3a8SDave May     }
71df21e3a8SDave May   }
7208056efcSDave May   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
73df21e3a8SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
74df21e3a8SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
75df21e3a8SDave May   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
76df21e3a8SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
77df21e3a8SDave May   ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr);
78df21e3a8SDave May   for (p=0; p<n_points_recv; p++) {
79df21e3a8SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
80df21e3a8SDave May 
81df21e3a8SDave May     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
82df21e3a8SDave May   }
83df21e3a8SDave May   ierr = DataExView(de);CHKERRQ(ierr);
84df21e3a8SDave May   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
85df21e3a8SDave May   ierr = DataExDestroy(de);CHKERRQ(ierr);
86df21e3a8SDave May   PetscFunctionReturn(0);
87df21e3a8SDave May }
882712d1f2SDave May 
8940c453e9SDave May #undef __FUNCT__
90889dbfe5SDave May #define __FUNCT__ "DMSwarmMigrate_DMNeighborScatter"
91889dbfe5SDave May PetscErrorCode DMSwarmMigrate_DMNeighborScatter(DM dm,DM dmcell,PetscBool remove_sent_points,PetscInt *npoints_prior_migration)
9240c453e9SDave May {
9340c453e9SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
9440c453e9SDave May   PetscErrorCode ierr;
9540c453e9SDave May   DataEx de;
9640c453e9SDave May   PetscInt r,p,npoints,*rankval,n_points_recv;
9740c453e9SDave May   PetscMPIInt rank,_rank;
9840c453e9SDave May   const PetscMPIInt *neighbourranks;
9940c453e9SDave May   void *point_buffer,*recv_points;
10040c453e9SDave May   size_t sizeof_dmswarm_point;
10140c453e9SDave May   PetscInt nneighbors;
1027c6d1d28SDave May   PetscMPIInt mynneigh,*myneigh;
10340c453e9SDave May 
104521f74f9SMatthew G. Knepley   PetscFunctionBegin;
10540c453e9SDave May   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
10640c453e9SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
10708056efcSDave May   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
108521f74f9SMatthew G. Knepley   ierr = DataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr);
10940c453e9SDave May   ierr = DMGetNeighbors(dmcell,&nneighbors,&neighbourranks);CHKERRQ(ierr);
11040c453e9SDave May   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
11140c453e9SDave May   for (r=0; r<nneighbors; r++) {
11240c453e9SDave May     _rank = neighbourranks[r];
11340c453e9SDave May     if ((_rank != rank) && (_rank > 0)) {
11440c453e9SDave May       ierr = DataExTopologyAddNeighbour(de,_rank);CHKERRQ(ierr);
11540c453e9SDave May     }
11640c453e9SDave May   }
11740c453e9SDave May   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
1187c6d1d28SDave May   ierr = DataExTopologyGetNeighbours(de,&mynneigh,&myneigh);CHKERRQ(ierr);
11940c453e9SDave May   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
12040c453e9SDave May   for (p=0; p<npoints; p++) {
121f954cb40SDave May     if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
1227c6d1d28SDave May       for (r=0; r<mynneigh; r++) {
1237c6d1d28SDave May         _rank = myneigh[r];
12440c453e9SDave May         ierr = DataExAddToSendCount(de,_rank,1);CHKERRQ(ierr);
12540c453e9SDave May       }
12640c453e9SDave May     }
12740c453e9SDave May   }
12840c453e9SDave May   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
12940c453e9SDave May   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
13040c453e9SDave May   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
13140c453e9SDave May   for (p=0; p<npoints; p++) {
132f954cb40SDave May     if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
1337c6d1d28SDave May       for (r=0; r<mynneigh; r++) {
1347c6d1d28SDave May         _rank = myneigh[r];
13540c453e9SDave May         /* copy point into buffer */
13640c453e9SDave May         ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr);
13740c453e9SDave May         /* insert point buffer into DataExchanger */
13840c453e9SDave May         ierr = DataExPackData(de,_rank,1,point_buffer);CHKERRQ(ierr);
13940c453e9SDave May       }
14040c453e9SDave May     }
14140c453e9SDave May   }
14240c453e9SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
1437c6d1d28SDave May   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
14440c453e9SDave May   if (remove_sent_points) {
1457c6d1d28SDave May     DataField PField;
1467c6d1d28SDave May 
1477c6d1d28SDave May     ierr = DataBucketGetDataFieldByName(swarm->db,DMSwarmField_rank,&PField);CHKERRQ(ierr);
1487c6d1d28SDave May     ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr);
14940c453e9SDave May     /* remove points which left processor */
15040c453e9SDave May     ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
15140c453e9SDave May     for (p=0; p<npoints; p++) {
152f954cb40SDave May       if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
15340c453e9SDave May         /* kill point */
15440c453e9SDave May         ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr);
1557c6d1d28SDave May         ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */
1567c6d1d28SDave May         ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); /* update date point increase realloc performed */
15740c453e9SDave May         p--; /* check replacement point */
15840c453e9SDave May       }
15940c453e9SDave May     }
16040c453e9SDave May   }
16140c453e9SDave May   ierr = DataBucketGetSizes(swarm->db,npoints_prior_migration,NULL,NULL);CHKERRQ(ierr);
16240c453e9SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
16340c453e9SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
16440c453e9SDave May   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
16540c453e9SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
16640c453e9SDave May   ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr);
16740c453e9SDave May   for (p=0; p<n_points_recv; p++) {
16840c453e9SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
16940c453e9SDave May 
17040c453e9SDave May     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
17140c453e9SDave May   }
17240c453e9SDave May   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
17340c453e9SDave May   ierr = DataExDestroy(de);CHKERRQ(ierr);
17440c453e9SDave May   PetscFunctionReturn(0);
17540c453e9SDave May }
176480eef7bSDave May 
177480eef7bSDave May #undef __FUNCT__
17808056efcSDave May #define __FUNCT__ "DMSwarmMigrate_CellDMScatter"
17908056efcSDave May PetscErrorCode DMSwarmMigrate_CellDMScatter(DM dm,PetscBool remove_sent_points)
180480eef7bSDave May {
181480eef7bSDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
182480eef7bSDave May   PetscErrorCode ierr;
1836fbf25f8SDave May   PetscInt p,npoints,npointsg=0,npoints2,npoints2g,*rankval,npoints_prior_migration;
184bbe8250bSMatthew G. Knepley   PetscSF sfcell = NULL;
185dfc14de9SMatthew G. Knepley   const PetscSFNode *LA_sfcell;
186480eef7bSDave May   DM dmcell;
187480eef7bSDave May   Vec pos;
18840c453e9SDave May   PetscBool error_check = swarm->migrate_error_on_missing_point;
1897c6d1d28SDave May   PetscMPIInt commsize,rank;
190480eef7bSDave May 
191521f74f9SMatthew G. Knepley   PetscFunctionBegin;
192480eef7bSDave May   ierr = DMSwarmGetCellDM(dm,&dmcell);CHKERRQ(ierr);
193480eef7bSDave May   if (!dmcell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid if cell DM provided");
194480eef7bSDave May 
1957c6d1d28SDave May   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&commsize);CHKERRQ(ierr);
1967c6d1d28SDave May   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
1977c6d1d28SDave May 
198fb1bcc12SMatthew G. Knepley   ierr = DMSwarmCreateLocalVectorFromField(dm, DMSwarmPICField_coor, &pos);CHKERRQ(ierr);
199dfc14de9SMatthew G. Knepley   ierr = DMLocatePoints(dmcell, pos, DM_POINTLOCATION_NONE, &sfcell);CHKERRQ(ierr);
200fb1bcc12SMatthew G. Knepley   ierr = DMSwarmDestroyLocalVectorFromField(dm, DMSwarmPICField_coor, &pos);CHKERRQ(ierr);
201480eef7bSDave May 
20240c453e9SDave May   if (error_check) {
20340c453e9SDave May     ierr = DMSwarmGetSize(dm,&npointsg);CHKERRQ(ierr);
20440c453e9SDave May   }
205480eef7bSDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
20608056efcSDave May   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
207dfc14de9SMatthew G. Knepley   ierr = PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell);CHKERRQ(ierr);
208480eef7bSDave May   for (p=0; p<npoints; p++) {
209dfc14de9SMatthew G. Knepley     rankval[p] = LA_sfcell[p].index;
210480eef7bSDave May   }
21108056efcSDave May   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
212dfc14de9SMatthew G. Knepley   ierr = PetscSFDestroy(&sfcell);CHKERRQ(ierr);
213480eef7bSDave May 
2146fbf25f8SDave May   if (commsize > 1) {
215889dbfe5SDave May     ierr = DMSwarmMigrate_DMNeighborScatter(dm,dmcell,remove_sent_points,&npoints_prior_migration);CHKERRQ(ierr);
2166fbf25f8SDave May   } else {
2176fbf25f8SDave May     ierr = DMSwarmGetSize(dm,&npoints_prior_migration);CHKERRQ(ierr);
2186fbf25f8SDave May   }
219480eef7bSDave May 
22040c453e9SDave May   /* locate points newly recevied */
22140c453e9SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr);
2227c6d1d28SDave May #if 0
22340c453e9SDave May   len = npoints2 - npoints_prior_migration;
22440c453e9SDave May   if (len > 0) {
22540c453e9SDave May     PetscScalar *LA_coor;
2267c6d1d28SDave May     PetscInt bs;
22740c453e9SDave May 
2287c6d1d28SDave May     ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr);
2297c6d1d28SDave May     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,bs*len,(const PetscScalar*)&LA_coor[bs*npoints_prior_migration],&pos);CHKERRQ(ierr);
23040c453e9SDave May 
23140c453e9SDave May     ierr = DMLocatePoints(dmcell,pos,&iscell);CHKERRQ(ierr);
23240c453e9SDave May 
23340c453e9SDave May     ierr = VecDestroy(&pos);CHKERRQ(ierr);
2347c6d1d28SDave May     ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr);
23540c453e9SDave May 
23640c453e9SDave May     ierr = ISGetIndices(iscell,&LA_iscell);CHKERRQ(ierr);
23708056efcSDave May     ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
23840c453e9SDave May     for (p=0; p<len; p++) {
239f954cb40SDave May       rankval[npoints_prior_migration+p] = LA_iscell[p];
24040c453e9SDave May     }
24140c453e9SDave May     ierr = ISRestoreIndices(iscell,&LA_iscell);CHKERRQ(ierr);
24240c453e9SDave May     ierr = ISDestroy(&iscell);CHKERRQ(ierr);
2437c6d1d28SDave May     ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
24440c453e9SDave May 
24540c453e9SDave May     /* remove points which left processor */
2467c6d1d28SDave May     ierr = DataBucketGetDataFieldByName(swarm->db,DMSwarmField_rank,&PField);CHKERRQ(ierr);
2477c6d1d28SDave May     ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr);
2487c6d1d28SDave May 
24940c453e9SDave May     ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
25040c453e9SDave May     for (p=npoints_prior_migration; p<npoints2; p++) {
251f954cb40SDave May       if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
25240c453e9SDave May         /* kill point */
25340c453e9SDave May         ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr);
2547c6d1d28SDave May         ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */
2557c6d1d28SDave May         ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); /* update date point increase realloc performed */
25640c453e9SDave May         p--; /* check replacement point */
25740c453e9SDave May       }
25840c453e9SDave May     }
2597c6d1d28SDave May 
2607c6d1d28SDave May   }
2617c6d1d28SDave May #endif
2627c6d1d28SDave May 
2637c6d1d28SDave May   {
2647c6d1d28SDave May     PetscScalar *LA_coor;
2657c6d1d28SDave May     PetscInt bs;
2667c6d1d28SDave May     DataField PField;
2677c6d1d28SDave May 
2687c6d1d28SDave May     ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr);
2697c6d1d28SDave May     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,bs*npoints2,(const PetscScalar*)LA_coor,&pos);CHKERRQ(ierr);
270dfc14de9SMatthew G. Knepley     ierr = DMLocatePoints(dmcell,pos,DM_POINTLOCATION_NONE,&sfcell);CHKERRQ(ierr);
2717c6d1d28SDave May 
2727c6d1d28SDave May     ierr = VecDestroy(&pos);CHKERRQ(ierr);
2737c6d1d28SDave May     ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr);
2747c6d1d28SDave May 
275dfc14de9SMatthew G. Knepley     ierr = PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell);CHKERRQ(ierr);
2767c6d1d28SDave May     ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
2777c6d1d28SDave May     for (p=0; p<npoints2; p++) {
278dfc14de9SMatthew G. Knepley       rankval[p] = LA_sfcell[p].index;
2797c6d1d28SDave May     }
280dfc14de9SMatthew G. Knepley     ierr = PetscSFDestroy(&sfcell);CHKERRQ(ierr);
28108056efcSDave May     ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
28240c453e9SDave May 
2837c6d1d28SDave May     /* remove points which left processor */
2847c6d1d28SDave May     ierr = DataBucketGetDataFieldByName(swarm->db,DMSwarmField_rank,&PField);CHKERRQ(ierr);
2857c6d1d28SDave May     ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr);
2867c6d1d28SDave May 
2877c6d1d28SDave May     ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr);
2887c6d1d28SDave May     for (p=0; p<npoints2; p++) {
289f954cb40SDave May       if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
2907c6d1d28SDave May         /* kill point */
2917c6d1d28SDave May         ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr);
2927c6d1d28SDave May         ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */
2937c6d1d28SDave May         ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); /* update date point increase realloc performed */
2947c6d1d28SDave May         p--; /* check replacement point */
2957c6d1d28SDave May       }
2967c6d1d28SDave May     }
29740c453e9SDave May   }
29840c453e9SDave May   /* check for error on removed points */
29940c453e9SDave May   if (error_check) {
30040c453e9SDave May     ierr = DMSwarmGetSize(dm,&npoints2g);CHKERRQ(ierr);
30140c453e9SDave 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);
30240c453e9SDave May   }
303480eef7bSDave May   PetscFunctionReturn(0);
304480eef7bSDave May }
305480eef7bSDave May 
30608056efcSDave May #undef __FUNCT__
30708056efcSDave May #define __FUNCT__ "DMSwarmMigrate_CellDMExact"
30808056efcSDave May PetscErrorCode DMSwarmMigrate_CellDMExact(DM dm,PetscBool remove_sent_points)
30908056efcSDave May {
310521f74f9SMatthew G. Knepley   PetscFunctionBegin;
31108056efcSDave May   PetscFunctionReturn(0);
31208056efcSDave May }
31308056efcSDave May 
314480eef7bSDave May /*
315480eef7bSDave May  Redundant as this assumes points can only be sent to a single rank
316480eef7bSDave May */
3172712d1f2SDave May #undef __FUNCT__
3182712d1f2SDave May #define __FUNCT__ "DMSwarmMigrate_GlobalToLocal_Basic"
3192712d1f2SDave May PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm,PetscInt *globalsize)
3202712d1f2SDave May {
3212712d1f2SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
3222712d1f2SDave May   PetscErrorCode ierr;
3232712d1f2SDave May   DataEx de;
3242712d1f2SDave May   PetscInt p,npoints,*rankval,n_points_recv;
3252712d1f2SDave May   PetscMPIInt rank,nrank,negrank;
3262712d1f2SDave May   void *point_buffer,*recv_points;
3272712d1f2SDave May   size_t sizeof_dmswarm_point;
3282712d1f2SDave May 
329521f74f9SMatthew G. Knepley   PetscFunctionBegin;
3302712d1f2SDave May   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
3312712d1f2SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
3322712d1f2SDave May   *globalsize = npoints;
33308056efcSDave May   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
334521f74f9SMatthew G. Knepley   ierr = DataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr);
3352712d1f2SDave May   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
3362712d1f2SDave May   for (p=0; p<npoints; p++) {
3372712d1f2SDave May     negrank = rankval[p];
3382712d1f2SDave May     if (negrank < 0) {
3392712d1f2SDave May       nrank = -negrank - 1;
3402712d1f2SDave May       ierr = DataExTopologyAddNeighbour(de,nrank);CHKERRQ(ierr);
3412712d1f2SDave May     }
3422712d1f2SDave May   }
3432712d1f2SDave May   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
3442712d1f2SDave May   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
3452712d1f2SDave May   for (p=0; p<npoints; p++) {
3462712d1f2SDave May     negrank = rankval[p];
3472712d1f2SDave May     if (negrank < 0) {
3482712d1f2SDave May       nrank = -negrank - 1;
3492712d1f2SDave May       ierr = DataExAddToSendCount(de,nrank,1);CHKERRQ(ierr);
3502712d1f2SDave May     }
3512712d1f2SDave May   }
3522712d1f2SDave May   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
3532712d1f2SDave May   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
3542712d1f2SDave May   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
3552712d1f2SDave May   for (p=0; p<npoints; p++) {
3562712d1f2SDave May     negrank = rankval[p];
3572712d1f2SDave May     if (negrank < 0) {
3582712d1f2SDave May       nrank = -negrank - 1;
3592712d1f2SDave May       rankval[p] = nrank;
3602712d1f2SDave May       /* copy point into buffer */
3612712d1f2SDave May       ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr);
3622712d1f2SDave May       /* insert point buffer into DataExchanger */
3632712d1f2SDave May       ierr = DataExPackData(de,nrank,1,point_buffer);CHKERRQ(ierr);
3642712d1f2SDave May       rankval[p] = negrank;
3652712d1f2SDave May     }
3662712d1f2SDave May   }
3672712d1f2SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
36808056efcSDave May   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
3692712d1f2SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
3702712d1f2SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
3712712d1f2SDave May   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
3722712d1f2SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
3732712d1f2SDave May   ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr);
3742712d1f2SDave May   for (p=0; p<n_points_recv; p++) {
3752712d1f2SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
3762712d1f2SDave May 
3772712d1f2SDave May     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
3782712d1f2SDave May   }
3792712d1f2SDave May   ierr = DataExView(de);CHKERRQ(ierr);
3802712d1f2SDave May   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
3812712d1f2SDave May   ierr = DataExDestroy(de);CHKERRQ(ierr);
3822712d1f2SDave May   PetscFunctionReturn(0);
3832712d1f2SDave May }
384b16650c8SDave May 
385b16650c8SDave May typedef struct {
386b16650c8SDave May   PetscMPIInt owner_rank;
387b16650c8SDave May   PetscReal min[3],max[3];
388b16650c8SDave May } CollectBBox;
389b16650c8SDave May 
390b16650c8SDave May #undef __FUNCT__
391fe39f135SDave May #define __FUNCT__ "DMSwarmCollect_DMDABoundingBox"
392fe39f135SDave May PETSC_EXTERN PetscErrorCode DMSwarmCollect_DMDABoundingBox(DM dm,PetscInt *globalsize)
393b16650c8SDave May {
394b16650c8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
395b16650c8SDave May   PetscErrorCode ierr;
396b16650c8SDave May   DataEx de;
397b16650c8SDave May   PetscInt p,pk,npoints,*rankval,n_points_recv,n_bbox_recv,dim,neighbour_cells;
398b16650c8SDave May   PetscMPIInt rank,nrank;
399b16650c8SDave May   void *point_buffer,*recv_points;
400b16650c8SDave May   size_t sizeof_dmswarm_point,sizeof_bbox_ctx;
401b16650c8SDave May   PetscBool isdmda;
402b16650c8SDave May   CollectBBox *bbox,*recv_bbox;
403b16650c8SDave May   const PetscMPIInt *dmneighborranks;
404b16650c8SDave May   DM dmcell;
405b16650c8SDave May 
406521f74f9SMatthew G. Knepley   PetscFunctionBegin;
407b16650c8SDave May   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
408b16650c8SDave May 
409fe39f135SDave May   ierr = DMSwarmGetCellDM(dm,&dmcell);CHKERRQ(ierr);
410fe39f135SDave May   if (!dmcell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid if cell DM provided");
411b16650c8SDave May   isdmda = PETSC_FALSE;
412b16650c8SDave May   PetscObjectTypeCompare((PetscObject)dmcell,DMDA,&isdmda);
413b16650c8SDave May   if (!isdmda) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only DMDA support for CollectBoundingBox");
414b16650c8SDave May 
4158dbd68bcSDave May   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
416b16650c8SDave May   sizeof_bbox_ctx = sizeof(CollectBBox);
417b16650c8SDave May   PetscMalloc1(1,&bbox);
418b16650c8SDave May   bbox->owner_rank = rank;
419b16650c8SDave May 
420b16650c8SDave May   /* compute the bounding box based on the overlapping / stenctil size */
421b16650c8SDave May   {
422b16650c8SDave May     Vec lcoor;
423b16650c8SDave May 
424b16650c8SDave May     ierr = DMGetCoordinatesLocal(dmcell,&lcoor);CHKERRQ(ierr);
425fe39f135SDave May     if (dim >= 1) {
426b16650c8SDave May       ierr = VecStrideMin(lcoor,0,NULL,&bbox->min[0]);CHKERRQ(ierr);
427b16650c8SDave May       ierr = VecStrideMax(lcoor,0,NULL,&bbox->max[0]);CHKERRQ(ierr);
428fe39f135SDave May     }
429fe39f135SDave May     if (dim >= 2) {
430fe39f135SDave May       ierr = VecStrideMin(lcoor,1,NULL,&bbox->min[1]);CHKERRQ(ierr);
431b16650c8SDave May       ierr = VecStrideMax(lcoor,1,NULL,&bbox->max[1]);CHKERRQ(ierr);
432b16650c8SDave May     }
433fe39f135SDave May     if (dim == 3) {
434fe39f135SDave May       ierr = VecStrideMin(lcoor,2,NULL,&bbox->min[2]);CHKERRQ(ierr);
435fe39f135SDave May       ierr = VecStrideMax(lcoor,2,NULL,&bbox->max[2]);CHKERRQ(ierr);
436fe39f135SDave May     }
437fe39f135SDave May   }
438b16650c8SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
439b16650c8SDave May   *globalsize = npoints;
44008056efcSDave May   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
441521f74f9SMatthew G. Knepley   ierr = DataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr);
442b16650c8SDave May   /* use DMDA neighbours */
443b16650c8SDave May   ierr = DMDAGetNeighbors(dmcell,&dmneighborranks);CHKERRQ(ierr);
4448dbd68bcSDave May   if (dim == 1) {
4458dbd68bcSDave May     neighbour_cells = 3;
4468dbd68bcSDave May   } else if (dim == 2) {
4478dbd68bcSDave May     neighbour_cells = 9;
4488dbd68bcSDave May   } else {
4498dbd68bcSDave May     neighbour_cells = 27;
4508dbd68bcSDave May   }
451b16650c8SDave May   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
452b16650c8SDave May   for (p=0; p<neighbour_cells; p++) {
453b16650c8SDave May     if ( (dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank) ) {
454b16650c8SDave May       ierr = DataExTopologyAddNeighbour(de,dmneighborranks[p]);CHKERRQ(ierr);
455b16650c8SDave May     }
456b16650c8SDave May   }
457b16650c8SDave May   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
458b16650c8SDave May   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
459b16650c8SDave May   for (p=0; p<neighbour_cells; p++) {
460b16650c8SDave May     if ( (dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank) ) {
461b16650c8SDave May       ierr = DataExAddToSendCount(de,dmneighborranks[p],1);CHKERRQ(ierr);
462b16650c8SDave May     }
463b16650c8SDave May   }
464b16650c8SDave May   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
465b16650c8SDave May   /* send bounding boxes */
466b16650c8SDave May   ierr = DataExPackInitialize(de,sizeof_bbox_ctx);CHKERRQ(ierr);
467b16650c8SDave May   for (p=0; p<neighbour_cells; p++) {
468b16650c8SDave May     nrank = dmneighborranks[p];
469b16650c8SDave May     if ( (nrank >= 0) && (nrank != rank) ) {
470b16650c8SDave May       /* insert bbox buffer into DataExchanger */
471b16650c8SDave May       ierr = DataExPackData(de,nrank,1,bbox);CHKERRQ(ierr);
472b16650c8SDave May     }
473b16650c8SDave May   }
474b16650c8SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
475b16650c8SDave May   /* recv bounding boxes */
476b16650c8SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
477b16650c8SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
478b16650c8SDave May   ierr = DataExGetRecvData(de,&n_bbox_recv,(void**)&recv_bbox);CHKERRQ(ierr);
479b16650c8SDave May   for (p=0; p<n_bbox_recv; p++) {
480*b4b02483SDave May     PetscPrintf(PETSC_COMM_SELF,"[rank %d]: box from %d : range[%+1.4e,%+1.4e]x[%+1.4e,%+1.4e]\n",rank,recv_bbox[p].owner_rank,
481*b4b02483SDave May            (double)recv_bbox[p].min[0],(double)recv_bbox[p].max[0],(double)recv_bbox[p].min[1],(double)recv_bbox[p].max[1]);
482b16650c8SDave May   }
483b16650c8SDave May   /* of course this is stupid as this "generic" function should have a better way to know what the coordinates are called */
484b16650c8SDave May   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
485b16650c8SDave May   for (pk=0; pk<n_bbox_recv; pk++) {
486b16650c8SDave May     PetscReal *array_x,*array_y;
487b16650c8SDave May 
488b16650c8SDave May     ierr = DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
489b16650c8SDave May     ierr = DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
490b16650c8SDave May     for (p=0; p<npoints; p++) {
491b16650c8SDave May       if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0]) ) {
492b16650c8SDave May         if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1]) ) {
493b16650c8SDave May           ierr = DataExAddToSendCount(de,recv_bbox[pk].owner_rank,1);CHKERRQ(ierr);
494b16650c8SDave May         }
495b16650c8SDave May       }
496b16650c8SDave May     }
497b16650c8SDave May     ierr = DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
498b16650c8SDave May     ierr = DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
499b16650c8SDave May   }
500b16650c8SDave May   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
501b16650c8SDave May   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
502b16650c8SDave May   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
503b16650c8SDave May   for (pk=0; pk<n_bbox_recv; pk++) {
504b16650c8SDave May     PetscReal *array_x,*array_y;
505b16650c8SDave May 
506b16650c8SDave May     ierr = DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
507b16650c8SDave May     ierr = DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
508b16650c8SDave May     for (p=0; p<npoints; p++) {
509b16650c8SDave May       if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0]) ) {
510b16650c8SDave May         if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1]) ) {
511521f74f9SMatthew G. Knepley           /* copy point into buffer */
512b16650c8SDave May           ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr);
513521f74f9SMatthew G. Knepley           /* insert point buffer into DataExchanger */
514b16650c8SDave May           ierr = DataExPackData(de,recv_bbox[pk].owner_rank,1,point_buffer);CHKERRQ(ierr);
515b16650c8SDave May         }
516b16650c8SDave May       }
517b16650c8SDave May     }
518b16650c8SDave May     ierr = DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
519b16650c8SDave May     ierr = DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
520b16650c8SDave May   }
521b16650c8SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
52208056efcSDave May   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
523b16650c8SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
524b16650c8SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
525b16650c8SDave May   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
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   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
534b16650c8SDave May   PetscFree(bbox);
535b16650c8SDave May   ierr = DataExView(de);CHKERRQ(ierr);
536b16650c8SDave May   ierr = DataExDestroy(de);CHKERRQ(ierr);
537b16650c8SDave May   PetscFunctionReturn(0);
538b16650c8SDave May }
539a9fd7477SDave May 
540a9fd7477SDave May 
541a9fd7477SDave May /* General collection when no order, or neighbour information is provided */
542a9fd7477SDave May /*
543a9fd7477SDave May  User provides context and collect() method
544a9fd7477SDave May  Broadcast user context
545a9fd7477SDave May 
546a9fd7477SDave May  for each context / rank {
547a9fd7477SDave May    collect(swarm,context,n,list)
548a9fd7477SDave May  }
549a9fd7477SDave May */
550a9fd7477SDave May #undef __FUNCT__
551a9fd7477SDave May #define __FUNCT__ "DMSwarmCollect_General"
552a9fd7477SDave May PETSC_EXTERN PetscErrorCode DMSwarmCollect_General(DM dm,PetscErrorCode (*collect)(DM,void*,PetscInt*,PetscInt**),size_t ctx_size,void *ctx,PetscInt *globalsize)
553a9fd7477SDave May {
554a9fd7477SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
555a9fd7477SDave May   PetscErrorCode ierr;
556a9fd7477SDave May   DataEx         de;
557a9fd7477SDave May   PetscInt       p,r,npoints,n_points_recv;
558a9fd7477SDave May   PetscMPIInt    commsize,rank;
559a9fd7477SDave May   void           *point_buffer,*recv_points;
560a9fd7477SDave May   void           *ctxlist;
561a9fd7477SDave May   PetscInt       *n2collect,**collectlist;
562a9fd7477SDave May   size_t         sizeof_dmswarm_point;
563a9fd7477SDave May 
564521f74f9SMatthew G. Knepley   PetscFunctionBegin;
565a9fd7477SDave May   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&commsize);CHKERRQ(ierr);
566a9fd7477SDave May   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
567a9fd7477SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
568a9fd7477SDave May   *globalsize = npoints;
569a9fd7477SDave May   /* Broadcast user context */
570a9fd7477SDave May   PetscMalloc(ctx_size*commsize,&ctxlist);
571a9fd7477SDave May   ierr = MPI_Allgather(ctx,ctx_size,MPI_CHAR,ctxlist,ctx_size,MPI_CHAR,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
572521f74f9SMatthew G. Knepley   ierr = PetscMalloc1(commsize,&n2collect);CHKERRQ(ierr);
573521f74f9SMatthew G. Knepley   ierr = PetscMalloc1(commsize,&collectlist);CHKERRQ(ierr);
574a9fd7477SDave May   for (r=0; r<commsize; r++) {
575a9fd7477SDave May     PetscInt _n2collect;
576a9fd7477SDave May     PetscInt *_collectlist;
577a9fd7477SDave May     void     *_ctx_r;
578a9fd7477SDave May 
579a9fd7477SDave May     _n2collect   = 0;
580a9fd7477SDave May     _collectlist = NULL;
581a9fd7477SDave May     if (r != rank) { /* don't collect data from yourself */
582a9fd7477SDave May       _ctx_r = (void*)( (char*)ctxlist + r * ctx_size );
583a9fd7477SDave May       ierr = collect(dm,_ctx_r,&_n2collect,&_collectlist);CHKERRQ(ierr);
584a9fd7477SDave May     }
585a9fd7477SDave May     n2collect[r]   = _n2collect;
586a9fd7477SDave May     collectlist[r] = _collectlist;
587a9fd7477SDave May   }
588521f74f9SMatthew G. Knepley   ierr = DataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr);
589a9fd7477SDave May   /* Define topology */
590a9fd7477SDave May   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
591a9fd7477SDave May   for (r=0; r<commsize; r++) {
592a9fd7477SDave May     if (n2collect[r] > 0) {
593a9fd7477SDave May       ierr = DataExTopologyAddNeighbour(de,(PetscMPIInt)r);CHKERRQ(ierr);
594a9fd7477SDave May     }
595a9fd7477SDave May   }
596a9fd7477SDave May   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
597a9fd7477SDave May   /* Define send counts */
598a9fd7477SDave May   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
599a9fd7477SDave May   for (r=0; r<commsize; r++) {
600a9fd7477SDave May     if (n2collect[r] > 0) {
601a9fd7477SDave May       ierr = DataExAddToSendCount(de,r,n2collect[r]);CHKERRQ(ierr);
602a9fd7477SDave May     }
603a9fd7477SDave May   }
604a9fd7477SDave May   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
605a9fd7477SDave May   /* Pack data */
606a9fd7477SDave May   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
607a9fd7477SDave May   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
608a9fd7477SDave May   for (r=0; r<commsize; r++) {
609a9fd7477SDave May     for (p=0; p<n2collect[r]; p++) {
610a9fd7477SDave May       ierr = DataBucketFillPackedArray(swarm->db,collectlist[r][p],point_buffer);CHKERRQ(ierr);
611a9fd7477SDave May       /* insert point buffer into the data exchanger */
612a9fd7477SDave May       ierr = DataExPackData(de,r,1,point_buffer);CHKERRQ(ierr);
613a9fd7477SDave May     }
614a9fd7477SDave May   }
615a9fd7477SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
616a9fd7477SDave May   /* Scatter */
617a9fd7477SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
618a9fd7477SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
619a9fd7477SDave May   /* Collect data in DMSwarm container */
620a9fd7477SDave May   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
621a9fd7477SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
622a9fd7477SDave May   ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr);
623a9fd7477SDave May   for (p=0; p<n_points_recv; p++) {
624a9fd7477SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
625a9fd7477SDave May 
626a9fd7477SDave May     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
627a9fd7477SDave May   }
628a9fd7477SDave May   /* Release memory */
629a9fd7477SDave May   for (r=0; r<commsize; r++) {
630a9fd7477SDave May     if (collectlist[r]) PetscFree(collectlist[r]);
631a9fd7477SDave May   }
632521f74f9SMatthew G. Knepley   ierr = PetscFree(collectlist);CHKERRQ(ierr);
633521f74f9SMatthew G. Knepley   ierr = PetscFree(n2collect);CHKERRQ(ierr);
634521f74f9SMatthew G. Knepley   ierr = PetscFree(ctxlist);CHKERRQ(ierr);
635a9fd7477SDave May   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
636a9fd7477SDave May   ierr = DataExView(de);CHKERRQ(ierr);
637a9fd7477SDave May   ierr = DataExDestroy(de);CHKERRQ(ierr);
638a9fd7477SDave May   PetscFunctionReturn(0);
639a9fd7477SDave May }
640a9fd7477SDave May 
641