xref: /petsc/src/dm/impls/swarm/swarm_migrate.c (revision 65141ba87fb79262f677d83274e7164876c1be87)
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 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 
22521f74f9SMatthew G. Knepley   PetscFunctionBegin;
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);
2608056efcSDave May   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
27521f74f9SMatthew G. Knepley   ierr = DataExCreate(PetscObjectComm((PetscObject)dm),0, &de);CHKERRQ(ierr);
28df21e3a8SDave May   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
29521f74f9SMatthew G. Knepley   for (p = 0; p < npoints; ++p) {
30df21e3a8SDave May     nrank = rankval[p];
31df21e3a8SDave May     if (nrank != rank) {
32df21e3a8SDave May       ierr = DataExTopologyAddNeighbour(de,nrank);CHKERRQ(ierr);
33df21e3a8SDave May     }
34df21e3a8SDave May   }
35df21e3a8SDave May   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
36df21e3a8SDave May   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
37df21e3a8SDave May   for (p=0; p<npoints; p++) {
38df21e3a8SDave May     nrank = rankval[p];
39df21e3a8SDave May     if (nrank != rank) {
40df21e3a8SDave May       ierr = DataExAddToSendCount(de,nrank,1);CHKERRQ(ierr);
41df21e3a8SDave May     }
42df21e3a8SDave May   }
43df21e3a8SDave May   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
44df21e3a8SDave May   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
45df21e3a8SDave May   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
46df21e3a8SDave May   for (p=0; p<npoints; p++) {
47df21e3a8SDave May     nrank = rankval[p];
48df21e3a8SDave May     if (nrank != rank) {
49df21e3a8SDave May       /* copy point into buffer */
50df21e3a8SDave May       ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr);
51df21e3a8SDave May       /* insert point buffer into DataExchanger */
52df21e3a8SDave May       ierr = DataExPackData(de,nrank,1,point_buffer);CHKERRQ(ierr);
53df21e3a8SDave May     }
54df21e3a8SDave May   }
55df21e3a8SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
5622a417f9SDave May   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
57df21e3a8SDave May 
58df21e3a8SDave May   if (remove_sent_points) {
5922a417f9SDave May     DataField gfield;
6022a417f9SDave May 
6122a417f9SDave May     ierr = DataBucketGetDataFieldByName(swarm->db,DMSwarmField_rank,&gfield);CHKERRQ(ierr);
6222a417f9SDave May     ierr = DataFieldGetAccess(gfield);CHKERRQ(ierr);
6322a417f9SDave May     ierr = DataFieldGetEntries(gfield,(void**)&rankval);CHKERRQ(ierr);
6422a417f9SDave May 
65df21e3a8SDave May     /* remove points which left processor */
66df21e3a8SDave May     ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
67df21e3a8SDave May     for (p=0; p<npoints; p++) {
68df21e3a8SDave May       nrank = rankval[p];
69df21e3a8SDave May       if (nrank != rank) {
70df21e3a8SDave May         /* kill point */
7122a417f9SDave May         ierr = DataFieldRestoreAccess(gfield);CHKERRQ(ierr);
7222a417f9SDave May 
73df21e3a8SDave May         ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr);
7422a417f9SDave May 
7522a417f9SDave May         ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */
7622a417f9SDave May         ierr = DataFieldGetAccess(gfield);CHKERRQ(ierr);
7722a417f9SDave May         ierr = DataFieldGetEntries(gfield,(void**)&rankval);CHKERRQ(ierr);
78df21e3a8SDave May         p--; /* check replacement point */
79df21e3a8SDave May       }
80df21e3a8SDave May     }
8122a417f9SDave May     ierr = DataFieldRestoreEntries(gfield,(void**)&rankval);CHKERRQ(ierr);
8222a417f9SDave May     ierr = DataFieldRestoreAccess(gfield);CHKERRQ(ierr);
83df21e3a8SDave May   }
84df21e3a8SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
85df21e3a8SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
86df21e3a8SDave May   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
87df21e3a8SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
88*65141ba8SDave May   ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr);
89df21e3a8SDave May   for (p=0; p<n_points_recv; p++) {
90df21e3a8SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
91df21e3a8SDave May 
92df21e3a8SDave May     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
93df21e3a8SDave May   }
94df21e3a8SDave May   ierr = DataExView(de);CHKERRQ(ierr);
95df21e3a8SDave May   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
96df21e3a8SDave May   ierr = DataExDestroy(de);CHKERRQ(ierr);
97df21e3a8SDave May   PetscFunctionReturn(0);
98df21e3a8SDave May }
992712d1f2SDave May 
100889dbfe5SDave May PetscErrorCode DMSwarmMigrate_DMNeighborScatter(DM dm,DM dmcell,PetscBool remove_sent_points,PetscInt *npoints_prior_migration)
10140c453e9SDave May {
10240c453e9SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
10340c453e9SDave May   PetscErrorCode ierr;
10440c453e9SDave May   DataEx de;
10540c453e9SDave May   PetscInt r,p,npoints,*rankval,n_points_recv;
10640c453e9SDave May   PetscMPIInt rank,_rank;
10740c453e9SDave May   const PetscMPIInt *neighbourranks;
10840c453e9SDave May   void *point_buffer,*recv_points;
10940c453e9SDave May   size_t sizeof_dmswarm_point;
11040c453e9SDave May   PetscInt nneighbors;
1117c6d1d28SDave May   PetscMPIInt mynneigh,*myneigh;
11240c453e9SDave May 
113521f74f9SMatthew G. Knepley   PetscFunctionBegin;
11440c453e9SDave May   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
11540c453e9SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
11608056efcSDave May   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
117521f74f9SMatthew G. Knepley   ierr = DataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr);
11840c453e9SDave May   ierr = DMGetNeighbors(dmcell,&nneighbors,&neighbourranks);CHKERRQ(ierr);
11940c453e9SDave May   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
12040c453e9SDave May   for (r=0; r<nneighbors; r++) {
12140c453e9SDave May     _rank = neighbourranks[r];
12240c453e9SDave May     if ((_rank != rank) && (_rank > 0)) {
12340c453e9SDave May       ierr = DataExTopologyAddNeighbour(de,_rank);CHKERRQ(ierr);
12440c453e9SDave May     }
12540c453e9SDave May   }
12640c453e9SDave May   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
1277c6d1d28SDave May   ierr = DataExTopologyGetNeighbours(de,&mynneigh,&myneigh);CHKERRQ(ierr);
12840c453e9SDave May   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
12940c453e9SDave May   for (p=0; p<npoints; p++) {
130f954cb40SDave May     if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
1317c6d1d28SDave May       for (r=0; r<mynneigh; r++) {
1327c6d1d28SDave May         _rank = myneigh[r];
13340c453e9SDave May         ierr = DataExAddToSendCount(de,_rank,1);CHKERRQ(ierr);
13440c453e9SDave May       }
13540c453e9SDave May     }
13640c453e9SDave May   }
13740c453e9SDave May   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
13840c453e9SDave May   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
13940c453e9SDave May   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
14040c453e9SDave May   for (p=0; p<npoints; p++) {
141f954cb40SDave May     if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
1427c6d1d28SDave May       for (r=0; r<mynneigh; r++) {
1437c6d1d28SDave May         _rank = myneigh[r];
14440c453e9SDave May         /* copy point into buffer */
14540c453e9SDave May         ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr);
14640c453e9SDave May         /* insert point buffer into DataExchanger */
14740c453e9SDave May         ierr = DataExPackData(de,_rank,1,point_buffer);CHKERRQ(ierr);
14840c453e9SDave May       }
14940c453e9SDave May     }
15040c453e9SDave May   }
15140c453e9SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
1527c6d1d28SDave May   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
15340c453e9SDave May   if (remove_sent_points) {
1547c6d1d28SDave May     DataField PField;
1557c6d1d28SDave May 
1567c6d1d28SDave May     ierr = DataBucketGetDataFieldByName(swarm->db,DMSwarmField_rank,&PField);CHKERRQ(ierr);
1577c6d1d28SDave May     ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr);
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++) {
161f954cb40SDave May       if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
16240c453e9SDave May         /* kill point */
16340c453e9SDave May         ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr);
1647c6d1d28SDave May         ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */
1657c6d1d28SDave May         ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); /* update date point increase realloc performed */
16640c453e9SDave May         p--; /* check replacement point */
16740c453e9SDave May       }
16840c453e9SDave May     }
16940c453e9SDave May   }
17040c453e9SDave May   ierr = DataBucketGetSizes(swarm->db,npoints_prior_migration,NULL,NULL);CHKERRQ(ierr);
17140c453e9SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
17240c453e9SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
17340c453e9SDave May   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
17440c453e9SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
175*65141ba8SDave May   ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr);
17640c453e9SDave May   for (p=0; p<n_points_recv; p++) {
17740c453e9SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
17840c453e9SDave May 
17940c453e9SDave May     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
18040c453e9SDave May   }
18140c453e9SDave May   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
18240c453e9SDave May   ierr = DataExDestroy(de);CHKERRQ(ierr);
18340c453e9SDave May   PetscFunctionReturn(0);
18440c453e9SDave May }
185480eef7bSDave May 
18608056efcSDave May PetscErrorCode DMSwarmMigrate_CellDMScatter(DM dm,PetscBool remove_sent_points)
187480eef7bSDave May {
188480eef7bSDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
189480eef7bSDave May   PetscErrorCode ierr;
1906fbf25f8SDave May   PetscInt p,npoints,npointsg=0,npoints2,npoints2g,*rankval,npoints_prior_migration;
191bbe8250bSMatthew G. Knepley   PetscSF sfcell = NULL;
192dfc14de9SMatthew G. Knepley   const PetscSFNode *LA_sfcell;
193480eef7bSDave May   DM dmcell;
194480eef7bSDave May   Vec pos;
19540c453e9SDave May   PetscBool error_check = swarm->migrate_error_on_missing_point;
1967c6d1d28SDave May   PetscMPIInt commsize,rank;
197480eef7bSDave May 
198521f74f9SMatthew G. Knepley   PetscFunctionBegin;
199480eef7bSDave May   ierr = DMSwarmGetCellDM(dm,&dmcell);CHKERRQ(ierr);
200480eef7bSDave May   if (!dmcell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid if cell DM provided");
201480eef7bSDave May 
2027c6d1d28SDave May   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&commsize);CHKERRQ(ierr);
2037c6d1d28SDave May   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
2047c6d1d28SDave May 
205fb1bcc12SMatthew G. Knepley   ierr = DMSwarmCreateLocalVectorFromField(dm, DMSwarmPICField_coor, &pos);CHKERRQ(ierr);
206dfc14de9SMatthew G. Knepley   ierr = DMLocatePoints(dmcell, pos, DM_POINTLOCATION_NONE, &sfcell);CHKERRQ(ierr);
207fb1bcc12SMatthew G. Knepley   ierr = DMSwarmDestroyLocalVectorFromField(dm, DMSwarmPICField_coor, &pos);CHKERRQ(ierr);
208480eef7bSDave May 
20940c453e9SDave May   if (error_check) {
21040c453e9SDave May     ierr = DMSwarmGetSize(dm,&npointsg);CHKERRQ(ierr);
21140c453e9SDave May   }
212480eef7bSDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
21308056efcSDave May   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
214dfc14de9SMatthew G. Knepley   ierr = PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell);CHKERRQ(ierr);
215480eef7bSDave May   for (p=0; p<npoints; p++) {
216dfc14de9SMatthew G. Knepley     rankval[p] = LA_sfcell[p].index;
217480eef7bSDave May   }
21808056efcSDave May   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
219dfc14de9SMatthew G. Knepley   ierr = PetscSFDestroy(&sfcell);CHKERRQ(ierr);
220480eef7bSDave May 
2216fbf25f8SDave May   if (commsize > 1) {
222889dbfe5SDave May     ierr = DMSwarmMigrate_DMNeighborScatter(dm,dmcell,remove_sent_points,&npoints_prior_migration);CHKERRQ(ierr);
2236fbf25f8SDave May   } else {
2246fbf25f8SDave May     ierr = DMSwarmGetSize(dm,&npoints_prior_migration);CHKERRQ(ierr);
2256fbf25f8SDave May   }
226480eef7bSDave May 
22740c453e9SDave May   /* locate points newly recevied */
22840c453e9SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr);
2297c6d1d28SDave May #if 0
23040c453e9SDave May   len = npoints2 - npoints_prior_migration;
23140c453e9SDave May   if (len > 0) {
23240c453e9SDave May     PetscScalar *LA_coor;
2337c6d1d28SDave May     PetscInt bs;
23440c453e9SDave May 
2357c6d1d28SDave May     ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr);
2367c6d1d28SDave May     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,bs*len,(const PetscScalar*)&LA_coor[bs*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);
2417c6d1d28SDave May     ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr);
24240c453e9SDave May 
24340c453e9SDave May     ierr = ISGetIndices(iscell,&LA_iscell);CHKERRQ(ierr);
24408056efcSDave May     ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
24540c453e9SDave May     for (p=0; p<len; p++) {
246f954cb40SDave May       rankval[npoints_prior_migration+p] = LA_iscell[p];
24740c453e9SDave May     }
24840c453e9SDave May     ierr = ISRestoreIndices(iscell,&LA_iscell);CHKERRQ(ierr);
24940c453e9SDave May     ierr = ISDestroy(&iscell);CHKERRQ(ierr);
2507c6d1d28SDave May     ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
25140c453e9SDave May 
25240c453e9SDave May     /* remove points which left processor */
2537c6d1d28SDave May     ierr = DataBucketGetDataFieldByName(swarm->db,DMSwarmField_rank,&PField);CHKERRQ(ierr);
2547c6d1d28SDave May     ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr);
2557c6d1d28SDave May 
25640c453e9SDave May     ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
25740c453e9SDave May     for (p=npoints_prior_migration; p<npoints2; p++) {
258f954cb40SDave May       if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
25940c453e9SDave May         /* kill point */
26040c453e9SDave May         ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr);
2617c6d1d28SDave May         ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */
2627c6d1d28SDave May         ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); /* update date point increase realloc performed */
26340c453e9SDave May         p--; /* check replacement point */
26440c453e9SDave May       }
26540c453e9SDave May     }
2667c6d1d28SDave May 
2677c6d1d28SDave May   }
2687c6d1d28SDave May #endif
2697c6d1d28SDave May 
2707c6d1d28SDave May   {
2717c6d1d28SDave May     PetscScalar *LA_coor;
2727c6d1d28SDave May     PetscInt bs;
2737c6d1d28SDave May     DataField PField;
2747c6d1d28SDave May 
2757c6d1d28SDave May     ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr);
2767c6d1d28SDave May     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,bs*npoints2,(const PetscScalar*)LA_coor,&pos);CHKERRQ(ierr);
277dfc14de9SMatthew G. Knepley     ierr = DMLocatePoints(dmcell,pos,DM_POINTLOCATION_NONE,&sfcell);CHKERRQ(ierr);
2787c6d1d28SDave May 
2797c6d1d28SDave May     ierr = VecDestroy(&pos);CHKERRQ(ierr);
2807c6d1d28SDave May     ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr);
2817c6d1d28SDave May 
282dfc14de9SMatthew G. Knepley     ierr = PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell);CHKERRQ(ierr);
2837c6d1d28SDave May     ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
2847c6d1d28SDave May     for (p=0; p<npoints2; p++) {
285dfc14de9SMatthew G. Knepley       rankval[p] = LA_sfcell[p].index;
2867c6d1d28SDave May     }
287dfc14de9SMatthew G. Knepley     ierr = PetscSFDestroy(&sfcell);CHKERRQ(ierr);
28808056efcSDave May     ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
28940c453e9SDave May 
2907c6d1d28SDave May     /* remove points which left processor */
2917c6d1d28SDave May     ierr = DataBucketGetDataFieldByName(swarm->db,DMSwarmField_rank,&PField);CHKERRQ(ierr);
2927c6d1d28SDave May     ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr);
2937c6d1d28SDave May 
2947c6d1d28SDave May     ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr);
2957c6d1d28SDave May     for (p=0; p<npoints2; p++) {
296f954cb40SDave May       if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
2977c6d1d28SDave May         /* kill point */
2987c6d1d28SDave May         ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr);
2997c6d1d28SDave May         ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */
3007c6d1d28SDave May         ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); /* update date point increase realloc performed */
3017c6d1d28SDave May         p--; /* check replacement point */
3027c6d1d28SDave May       }
3037c6d1d28SDave May     }
30440c453e9SDave May   }
30540c453e9SDave May   /* check for error on removed points */
30640c453e9SDave May   if (error_check) {
30740c453e9SDave May     ierr = DMSwarmGetSize(dm,&npoints2g);CHKERRQ(ierr);
30840c453e9SDave 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);
30940c453e9SDave May   }
310480eef7bSDave May   PetscFunctionReturn(0);
311480eef7bSDave May }
312480eef7bSDave May 
31308056efcSDave May PetscErrorCode DMSwarmMigrate_CellDMExact(DM dm,PetscBool remove_sent_points)
31408056efcSDave May {
315521f74f9SMatthew G. Knepley   PetscFunctionBegin;
31608056efcSDave May   PetscFunctionReturn(0);
31708056efcSDave May }
31808056efcSDave May 
319480eef7bSDave May /*
320480eef7bSDave May  Redundant as this assumes points can only be sent to a single rank
321480eef7bSDave May */
3222712d1f2SDave May PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm,PetscInt *globalsize)
3232712d1f2SDave May {
3242712d1f2SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
3252712d1f2SDave May   PetscErrorCode ierr;
3262712d1f2SDave May   DataEx de;
3272712d1f2SDave May   PetscInt p,npoints,*rankval,n_points_recv;
3282712d1f2SDave May   PetscMPIInt rank,nrank,negrank;
3292712d1f2SDave May   void *point_buffer,*recv_points;
3302712d1f2SDave May   size_t sizeof_dmswarm_point;
3312712d1f2SDave May 
332521f74f9SMatthew G. Knepley   PetscFunctionBegin;
3332712d1f2SDave May   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
3342712d1f2SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
3352712d1f2SDave May   *globalsize = npoints;
33608056efcSDave May   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
337521f74f9SMatthew G. Knepley   ierr = DataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr);
3382712d1f2SDave May   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
3392712d1f2SDave May   for (p=0; p<npoints; p++) {
3402712d1f2SDave May     negrank = rankval[p];
3412712d1f2SDave May     if (negrank < 0) {
3422712d1f2SDave May       nrank = -negrank - 1;
3432712d1f2SDave May       ierr = DataExTopologyAddNeighbour(de,nrank);CHKERRQ(ierr);
3442712d1f2SDave May     }
3452712d1f2SDave May   }
3462712d1f2SDave May   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
3472712d1f2SDave May   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
3482712d1f2SDave May   for (p=0; p<npoints; p++) {
3492712d1f2SDave May     negrank = rankval[p];
3502712d1f2SDave May     if (negrank < 0) {
3512712d1f2SDave May       nrank = -negrank - 1;
3522712d1f2SDave May       ierr = DataExAddToSendCount(de,nrank,1);CHKERRQ(ierr);
3532712d1f2SDave May     }
3542712d1f2SDave May   }
3552712d1f2SDave May   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
3562712d1f2SDave May   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
3572712d1f2SDave May   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
3582712d1f2SDave May   for (p=0; p<npoints; p++) {
3592712d1f2SDave May     negrank = rankval[p];
3602712d1f2SDave May     if (negrank < 0) {
3612712d1f2SDave May       nrank = -negrank - 1;
3622712d1f2SDave May       rankval[p] = nrank;
3632712d1f2SDave May       /* copy point into buffer */
3642712d1f2SDave May       ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr);
3652712d1f2SDave May       /* insert point buffer into DataExchanger */
3662712d1f2SDave May       ierr = DataExPackData(de,nrank,1,point_buffer);CHKERRQ(ierr);
3672712d1f2SDave May       rankval[p] = negrank;
3682712d1f2SDave May     }
3692712d1f2SDave May   }
3702712d1f2SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
37108056efcSDave May   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
3722712d1f2SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
3732712d1f2SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
3742712d1f2SDave May   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
3752712d1f2SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
376*65141ba8SDave May   ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr);
3772712d1f2SDave May   for (p=0; p<n_points_recv; p++) {
3782712d1f2SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
3792712d1f2SDave May 
3802712d1f2SDave May     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
3812712d1f2SDave May   }
3822712d1f2SDave May   ierr = DataExView(de);CHKERRQ(ierr);
3832712d1f2SDave May   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
3842712d1f2SDave May   ierr = DataExDestroy(de);CHKERRQ(ierr);
3852712d1f2SDave May   PetscFunctionReturn(0);
3862712d1f2SDave May }
387b16650c8SDave May 
388b16650c8SDave May typedef struct {
389b16650c8SDave May   PetscMPIInt owner_rank;
390b16650c8SDave May   PetscReal min[3],max[3];
391b16650c8SDave May } CollectBBox;
392b16650c8SDave May 
393fe39f135SDave May PETSC_EXTERN PetscErrorCode DMSwarmCollect_DMDABoundingBox(DM dm,PetscInt *globalsize)
394b16650c8SDave May {
395b16650c8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
396b16650c8SDave May   PetscErrorCode ierr;
397b16650c8SDave May   DataEx de;
398b16650c8SDave May   PetscInt p,pk,npoints,*rankval,n_points_recv,n_bbox_recv,dim,neighbour_cells;
399b16650c8SDave May   PetscMPIInt rank,nrank;
400b16650c8SDave May   void *point_buffer,*recv_points;
401b16650c8SDave May   size_t sizeof_dmswarm_point,sizeof_bbox_ctx;
402b16650c8SDave May   PetscBool isdmda;
403b16650c8SDave May   CollectBBox *bbox,*recv_bbox;
404b16650c8SDave May   const PetscMPIInt *dmneighborranks;
405b16650c8SDave May   DM dmcell;
406b16650c8SDave May 
407521f74f9SMatthew G. Knepley   PetscFunctionBegin;
408b16650c8SDave May   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
409b16650c8SDave May 
410fe39f135SDave May   ierr = DMSwarmGetCellDM(dm,&dmcell);CHKERRQ(ierr);
411fe39f135SDave May   if (!dmcell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid if cell DM provided");
412b16650c8SDave May   isdmda = PETSC_FALSE;
413b16650c8SDave May   PetscObjectTypeCompare((PetscObject)dmcell,DMDA,&isdmda);
414b16650c8SDave May   if (!isdmda) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only DMDA support for CollectBoundingBox");
415b16650c8SDave May 
4168dbd68bcSDave May   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
417b16650c8SDave May   sizeof_bbox_ctx = sizeof(CollectBBox);
418b16650c8SDave May   PetscMalloc1(1,&bbox);
419b16650c8SDave May   bbox->owner_rank = rank;
420b16650c8SDave May 
421b16650c8SDave May   /* compute the bounding box based on the overlapping / stenctil size */
422b16650c8SDave May   {
423b16650c8SDave May     Vec lcoor;
424b16650c8SDave May 
425b16650c8SDave May     ierr = DMGetCoordinatesLocal(dmcell,&lcoor);CHKERRQ(ierr);
426fe39f135SDave May     if (dim >= 1) {
427b16650c8SDave May       ierr = VecStrideMin(lcoor,0,NULL,&bbox->min[0]);CHKERRQ(ierr);
428b16650c8SDave May       ierr = VecStrideMax(lcoor,0,NULL,&bbox->max[0]);CHKERRQ(ierr);
429fe39f135SDave May     }
430fe39f135SDave May     if (dim >= 2) {
431fe39f135SDave May       ierr = VecStrideMin(lcoor,1,NULL,&bbox->min[1]);CHKERRQ(ierr);
432b16650c8SDave May       ierr = VecStrideMax(lcoor,1,NULL,&bbox->max[1]);CHKERRQ(ierr);
433b16650c8SDave May     }
434fe39f135SDave May     if (dim == 3) {
435fe39f135SDave May       ierr = VecStrideMin(lcoor,2,NULL,&bbox->min[2]);CHKERRQ(ierr);
436fe39f135SDave May       ierr = VecStrideMax(lcoor,2,NULL,&bbox->max[2]);CHKERRQ(ierr);
437fe39f135SDave May     }
438fe39f135SDave May   }
439b16650c8SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
440b16650c8SDave May   *globalsize = npoints;
44108056efcSDave May   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
442521f74f9SMatthew G. Knepley   ierr = DataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr);
443b16650c8SDave May   /* use DMDA neighbours */
444b16650c8SDave May   ierr = DMDAGetNeighbors(dmcell,&dmneighborranks);CHKERRQ(ierr);
4458dbd68bcSDave May   if (dim == 1) {
4468dbd68bcSDave May     neighbour_cells = 3;
4478dbd68bcSDave May   } else if (dim == 2) {
4488dbd68bcSDave May     neighbour_cells = 9;
4498dbd68bcSDave May   } else {
4508dbd68bcSDave May     neighbour_cells = 27;
4518dbd68bcSDave May   }
452b16650c8SDave May   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
453b16650c8SDave May   for (p=0; p<neighbour_cells; p++) {
454b16650c8SDave May     if ( (dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank) ) {
455b16650c8SDave May       ierr = DataExTopologyAddNeighbour(de,dmneighborranks[p]);CHKERRQ(ierr);
456b16650c8SDave May     }
457b16650c8SDave May   }
458b16650c8SDave May   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
459b16650c8SDave May   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
460b16650c8SDave May   for (p=0; p<neighbour_cells; p++) {
461b16650c8SDave May     if ( (dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank) ) {
462b16650c8SDave May       ierr = DataExAddToSendCount(de,dmneighborranks[p],1);CHKERRQ(ierr);
463b16650c8SDave May     }
464b16650c8SDave May   }
465b16650c8SDave May   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
466b16650c8SDave May   /* send bounding boxes */
467b16650c8SDave May   ierr = DataExPackInitialize(de,sizeof_bbox_ctx);CHKERRQ(ierr);
468b16650c8SDave May   for (p=0; p<neighbour_cells; p++) {
469b16650c8SDave May     nrank = dmneighborranks[p];
470b16650c8SDave May     if ( (nrank >= 0) && (nrank != rank) ) {
471b16650c8SDave May       /* insert bbox buffer into DataExchanger */
472b16650c8SDave May       ierr = DataExPackData(de,nrank,1,bbox);CHKERRQ(ierr);
473b16650c8SDave May     }
474b16650c8SDave May   }
475b16650c8SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
476b16650c8SDave May   /* recv bounding boxes */
477b16650c8SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
478b16650c8SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
479b16650c8SDave May   ierr = DataExGetRecvData(de,&n_bbox_recv,(void**)&recv_bbox);CHKERRQ(ierr);
480b16650c8SDave May   for (p=0; p<n_bbox_recv; p++) {
481b4b02483SDave 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,
482b4b02483SDave 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]);
483b16650c8SDave May   }
484b16650c8SDave May   /* of course this is stupid as this "generic" function should have a better way to know what the coordinates are called */
485b16650c8SDave May   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
486b16650c8SDave May   for (pk=0; pk<n_bbox_recv; pk++) {
487b16650c8SDave May     PetscReal *array_x,*array_y;
488b16650c8SDave May 
489b16650c8SDave May     ierr = DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
490b16650c8SDave May     ierr = DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
491b16650c8SDave May     for (p=0; p<npoints; p++) {
492b16650c8SDave May       if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0]) ) {
493b16650c8SDave May         if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1]) ) {
494b16650c8SDave May           ierr = DataExAddToSendCount(de,recv_bbox[pk].owner_rank,1);CHKERRQ(ierr);
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   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
503b16650c8SDave May   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
504b16650c8SDave May   for (pk=0; pk<n_bbox_recv; pk++) {
505b16650c8SDave May     PetscReal *array_x,*array_y;
506b16650c8SDave May 
507b16650c8SDave May     ierr = DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
508b16650c8SDave May     ierr = DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
509b16650c8SDave May     for (p=0; p<npoints; p++) {
510b16650c8SDave May       if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0]) ) {
511b16650c8SDave May         if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1]) ) {
512521f74f9SMatthew G. Knepley           /* copy point into buffer */
513b16650c8SDave May           ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr);
514521f74f9SMatthew G. Knepley           /* insert point buffer into DataExchanger */
515b16650c8SDave May           ierr = DataExPackData(de,recv_bbox[pk].owner_rank,1,point_buffer);CHKERRQ(ierr);
516b16650c8SDave May         }
517b16650c8SDave May       }
518b16650c8SDave May     }
519b16650c8SDave May     ierr = DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
520b16650c8SDave May     ierr = DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
521b16650c8SDave May   }
522b16650c8SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
52308056efcSDave May   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
524b16650c8SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
525b16650c8SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
526b16650c8SDave May   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
527b16650c8SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
528*65141ba8SDave May   ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr);
529b16650c8SDave May   for (p=0; p<n_points_recv; p++) {
530b16650c8SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
531b16650c8SDave May 
532b16650c8SDave May     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
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   PetscFunctionReturn(0);
539b16650c8SDave May }
540a9fd7477SDave May 
541a9fd7477SDave May 
542a9fd7477SDave May /* General collection when no order, or neighbour information is provided */
543a9fd7477SDave May /*
544a9fd7477SDave May  User provides context and collect() method
545a9fd7477SDave May  Broadcast user context
546a9fd7477SDave May 
547a9fd7477SDave May  for each context / rank {
548a9fd7477SDave May    collect(swarm,context,n,list)
549a9fd7477SDave May  }
550a9fd7477SDave May */
551a9fd7477SDave May PETSC_EXTERN PetscErrorCode DMSwarmCollect_General(DM dm,PetscErrorCode (*collect)(DM,void*,PetscInt*,PetscInt**),size_t ctx_size,void *ctx,PetscInt *globalsize)
552a9fd7477SDave May {
553a9fd7477SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
554a9fd7477SDave May   PetscErrorCode ierr;
555a9fd7477SDave May   DataEx         de;
556a9fd7477SDave May   PetscInt       p,r,npoints,n_points_recv;
557a9fd7477SDave May   PetscMPIInt    commsize,rank;
558a9fd7477SDave May   void           *point_buffer,*recv_points;
559a9fd7477SDave May   void           *ctxlist;
560a9fd7477SDave May   PetscInt       *n2collect,**collectlist;
561a9fd7477SDave May   size_t         sizeof_dmswarm_point;
562a9fd7477SDave May 
563521f74f9SMatthew G. Knepley   PetscFunctionBegin;
564a9fd7477SDave May   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&commsize);CHKERRQ(ierr);
565a9fd7477SDave May   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
566a9fd7477SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
567a9fd7477SDave May   *globalsize = npoints;
568a9fd7477SDave May   /* Broadcast user context */
569a9fd7477SDave May   PetscMalloc(ctx_size*commsize,&ctxlist);
570a9fd7477SDave May   ierr = MPI_Allgather(ctx,ctx_size,MPI_CHAR,ctxlist,ctx_size,MPI_CHAR,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
571521f74f9SMatthew G. Knepley   ierr = PetscMalloc1(commsize,&n2collect);CHKERRQ(ierr);
572521f74f9SMatthew G. Knepley   ierr = PetscMalloc1(commsize,&collectlist);CHKERRQ(ierr);
573a9fd7477SDave May   for (r=0; r<commsize; r++) {
574a9fd7477SDave May     PetscInt _n2collect;
575a9fd7477SDave May     PetscInt *_collectlist;
576a9fd7477SDave May     void     *_ctx_r;
577a9fd7477SDave May 
578a9fd7477SDave May     _n2collect   = 0;
579a9fd7477SDave May     _collectlist = NULL;
580a9fd7477SDave May     if (r != rank) { /* don't collect data from yourself */
581a9fd7477SDave May       _ctx_r = (void*)( (char*)ctxlist + r * ctx_size );
582a9fd7477SDave May       ierr = collect(dm,_ctx_r,&_n2collect,&_collectlist);CHKERRQ(ierr);
583a9fd7477SDave May     }
584a9fd7477SDave May     n2collect[r]   = _n2collect;
585a9fd7477SDave May     collectlist[r] = _collectlist;
586a9fd7477SDave May   }
587521f74f9SMatthew G. Knepley   ierr = DataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr);
588a9fd7477SDave May   /* Define topology */
589a9fd7477SDave May   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
590a9fd7477SDave May   for (r=0; r<commsize; r++) {
591a9fd7477SDave May     if (n2collect[r] > 0) {
592a9fd7477SDave May       ierr = DataExTopologyAddNeighbour(de,(PetscMPIInt)r);CHKERRQ(ierr);
593a9fd7477SDave May     }
594a9fd7477SDave May   }
595a9fd7477SDave May   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
596a9fd7477SDave May   /* Define send counts */
597a9fd7477SDave May   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
598a9fd7477SDave May   for (r=0; r<commsize; r++) {
599a9fd7477SDave May     if (n2collect[r] > 0) {
600a9fd7477SDave May       ierr = DataExAddToSendCount(de,r,n2collect[r]);CHKERRQ(ierr);
601a9fd7477SDave May     }
602a9fd7477SDave May   }
603a9fd7477SDave May   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
604a9fd7477SDave May   /* Pack data */
605a9fd7477SDave May   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
606a9fd7477SDave May   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
607a9fd7477SDave May   for (r=0; r<commsize; r++) {
608a9fd7477SDave May     for (p=0; p<n2collect[r]; p++) {
609a9fd7477SDave May       ierr = DataBucketFillPackedArray(swarm->db,collectlist[r][p],point_buffer);CHKERRQ(ierr);
610a9fd7477SDave May       /* insert point buffer into the data exchanger */
611a9fd7477SDave May       ierr = DataExPackData(de,r,1,point_buffer);CHKERRQ(ierr);
612a9fd7477SDave May     }
613a9fd7477SDave May   }
614a9fd7477SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
615a9fd7477SDave May   /* Scatter */
616a9fd7477SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
617a9fd7477SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
618a9fd7477SDave May   /* Collect data in DMSwarm container */
619a9fd7477SDave May   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
620a9fd7477SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
621*65141ba8SDave May   ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr);
622a9fd7477SDave May   for (p=0; p<n_points_recv; p++) {
623a9fd7477SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
624a9fd7477SDave May 
625a9fd7477SDave May     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
626a9fd7477SDave May   }
627a9fd7477SDave May   /* Release memory */
628a9fd7477SDave May   for (r=0; r<commsize; r++) {
629a9fd7477SDave May     if (collectlist[r]) PetscFree(collectlist[r]);
630a9fd7477SDave May   }
631521f74f9SMatthew G. Knepley   ierr = PetscFree(collectlist);CHKERRQ(ierr);
632521f74f9SMatthew G. Knepley   ierr = PetscFree(n2collect);CHKERRQ(ierr);
633521f74f9SMatthew G. Knepley   ierr = PetscFree(ctxlist);CHKERRQ(ierr);
634a9fd7477SDave May   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
635a9fd7477SDave May   ierr = DataExView(de);CHKERRQ(ierr);
636a9fd7477SDave May   ierr = DataExDestroy(de);CHKERRQ(ierr);
637a9fd7477SDave May   PetscFunctionReturn(0);
638a9fd7477SDave May }
639a9fd7477SDave May 
640