xref: /petsc/src/dm/impls/swarm/swarm_migrate.c (revision 5627991adfa0dae03dda7a34b1ec5f8d16d5e5f7)
1dfc14de9SMatthew G. Knepley #include <petscsf.h>
208056efcSDave May #include <petscdmswarm.h>
35917a6f0SStefano Zampini #include <petscdmda.h>
4df21e3a8SDave May #include <petsc/private/dmswarmimpl.h>    /*I   "petscdmswarm.h"   I*/
5279f676cSBarry Smith #include "../src/dm/impls/swarm/data_bucket.h"
6279f676cSBarry Smith #include "../src/dm/impls/swarm/data_ex.h"
7df21e3a8SDave May 
8480eef7bSDave May /*
9480eef7bSDave May  User loads desired location (MPI rank) into field DMSwarm_rank
10480eef7bSDave May */
11df21e3a8SDave May PetscErrorCode DMSwarmMigrate_Push_Basic(DM dm,PetscBool remove_sent_points)
12df21e3a8SDave May {
13df21e3a8SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
14df21e3a8SDave May   PetscErrorCode ierr;
1577048351SPatrick Sanan   DMSwarmDataEx  de;
16df21e3a8SDave May   PetscInt       p,npoints,*rankval,n_points_recv;
17df21e3a8SDave May   PetscMPIInt    rank,nrank;
18df21e3a8SDave May   void           *point_buffer,*recv_points;
19df21e3a8SDave May   size_t         sizeof_dmswarm_point;
20df21e3a8SDave May 
21521f74f9SMatthew G. Knepley   PetscFunctionBegin;
22ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRMPI(ierr);
23df21e3a8SDave May 
2477048351SPatrick Sanan   ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
2508056efcSDave May   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
2677048351SPatrick Sanan   ierr = DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm),0, &de);CHKERRQ(ierr);
2777048351SPatrick Sanan   ierr = DMSwarmDataExTopologyInitialize(de);CHKERRQ(ierr);
28521f74f9SMatthew G. Knepley   for (p = 0; p < npoints; ++p) {
29df21e3a8SDave May     nrank = rankval[p];
30df21e3a8SDave May     if (nrank != rank) {
3177048351SPatrick Sanan       ierr = DMSwarmDataExTopologyAddNeighbour(de,nrank);CHKERRQ(ierr);
32df21e3a8SDave May     }
33df21e3a8SDave May   }
3477048351SPatrick Sanan   ierr = DMSwarmDataExTopologyFinalize(de);CHKERRQ(ierr);
3577048351SPatrick Sanan   ierr = DMSwarmDataExInitializeSendCount(de);CHKERRQ(ierr);
36df21e3a8SDave May   for (p=0; p<npoints; p++) {
37df21e3a8SDave May     nrank = rankval[p];
38df21e3a8SDave May     if (nrank != rank) {
3977048351SPatrick Sanan       ierr = DMSwarmDataExAddToSendCount(de,nrank,1);CHKERRQ(ierr);
40df21e3a8SDave May     }
41df21e3a8SDave May   }
4277048351SPatrick Sanan   ierr = DMSwarmDataExFinalizeSendCount(de);CHKERRQ(ierr);
4377048351SPatrick Sanan   ierr = DMSwarmDataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
4477048351SPatrick Sanan   ierr = DMSwarmDataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
45df21e3a8SDave May   for (p=0; p<npoints; p++) {
46df21e3a8SDave May     nrank = rankval[p];
47df21e3a8SDave May     if (nrank != rank) {
48df21e3a8SDave May       /* copy point into buffer */
4977048351SPatrick Sanan       ierr = DMSwarmDataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr);
5077048351SPatrick Sanan       /* insert point buffer into DMSwarmDataExchanger */
5177048351SPatrick Sanan       ierr = DMSwarmDataExPackData(de,nrank,1,point_buffer);CHKERRQ(ierr);
52df21e3a8SDave May     }
53df21e3a8SDave May   }
5477048351SPatrick Sanan   ierr = DMSwarmDataExPackFinalize(de);CHKERRQ(ierr);
5522a417f9SDave May   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
56df21e3a8SDave May 
57df21e3a8SDave May   if (remove_sent_points) {
5877048351SPatrick Sanan     DMSwarmDataField gfield;
5922a417f9SDave May 
6077048351SPatrick Sanan     ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,DMSwarmField_rank,&gfield);CHKERRQ(ierr);
6177048351SPatrick Sanan     ierr = DMSwarmDataFieldGetAccess(gfield);CHKERRQ(ierr);
6277048351SPatrick Sanan     ierr = DMSwarmDataFieldGetEntries(gfield,(void**)&rankval);CHKERRQ(ierr);
6322a417f9SDave May 
64df21e3a8SDave May     /* remove points which left processor */
6577048351SPatrick Sanan     ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
66df21e3a8SDave May     for (p=0; p<npoints; p++) {
67df21e3a8SDave May       nrank = rankval[p];
68df21e3a8SDave May       if (nrank != rank) {
69df21e3a8SDave May         /* kill point */
7077048351SPatrick Sanan         ierr = DMSwarmDataFieldRestoreAccess(gfield);CHKERRQ(ierr);
7122a417f9SDave May 
7277048351SPatrick Sanan         ierr = DMSwarmDataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr);
7322a417f9SDave May 
7477048351SPatrick Sanan         ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */
7577048351SPatrick Sanan         ierr = DMSwarmDataFieldGetAccess(gfield);CHKERRQ(ierr);
7677048351SPatrick Sanan         ierr = DMSwarmDataFieldGetEntries(gfield,(void**)&rankval);CHKERRQ(ierr);
77df21e3a8SDave May         p--; /* check replacement point */
78df21e3a8SDave May       }
79df21e3a8SDave May     }
8077048351SPatrick Sanan     ierr = DMSwarmDataFieldRestoreEntries(gfield,(void**)&rankval);CHKERRQ(ierr);
8177048351SPatrick Sanan     ierr = DMSwarmDataFieldRestoreAccess(gfield);CHKERRQ(ierr);
82df21e3a8SDave May   }
8377048351SPatrick Sanan   ierr = DMSwarmDataExBegin(de);CHKERRQ(ierr);
8477048351SPatrick Sanan   ierr = DMSwarmDataExEnd(de);CHKERRQ(ierr);
8577048351SPatrick Sanan   ierr = DMSwarmDataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
8677048351SPatrick Sanan   ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
8777048351SPatrick Sanan   ierr = DMSwarmDataBucketSetSizes(swarm->db,npoints + n_points_recv,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr);
88df21e3a8SDave May   for (p=0; p<n_points_recv; p++) {
89df21e3a8SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point);
90df21e3a8SDave May 
9177048351SPatrick Sanan     ierr = DMSwarmDataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
92df21e3a8SDave May   }
9377048351SPatrick Sanan   ierr = DMSwarmDataExView(de);CHKERRQ(ierr);
9477048351SPatrick Sanan   ierr = DMSwarmDataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
9577048351SPatrick Sanan   ierr = DMSwarmDataExDestroy(de);CHKERRQ(ierr);
96df21e3a8SDave May   PetscFunctionReturn(0);
97df21e3a8SDave May }
982712d1f2SDave May 
99889dbfe5SDave May PetscErrorCode DMSwarmMigrate_DMNeighborScatter(DM dm,DM dmcell,PetscBool remove_sent_points,PetscInt *npoints_prior_migration)
10040c453e9SDave May {
10140c453e9SDave May   DM_Swarm          *swarm = (DM_Swarm*)dm->data;
10240c453e9SDave May   PetscErrorCode    ierr;
10377048351SPatrick Sanan   DMSwarmDataEx     de;
10440c453e9SDave May   PetscInt          r,p,npoints,*rankval,n_points_recv;
10540c453e9SDave May   PetscMPIInt       rank,_rank;
10640c453e9SDave May   const PetscMPIInt *neighbourranks;
10740c453e9SDave May   void              *point_buffer,*recv_points;
10840c453e9SDave May   size_t            sizeof_dmswarm_point;
10940c453e9SDave May   PetscInt          nneighbors;
1107c6d1d28SDave May   PetscMPIInt       mynneigh,*myneigh;
11140c453e9SDave May 
112521f74f9SMatthew G. Knepley   PetscFunctionBegin;
113ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRMPI(ierr);
11477048351SPatrick Sanan   ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
11508056efcSDave May   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
11677048351SPatrick Sanan   ierr = DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr);
11740c453e9SDave May   ierr = DMGetNeighbors(dmcell,&nneighbors,&neighbourranks);CHKERRQ(ierr);
11877048351SPatrick Sanan   ierr = DMSwarmDataExTopologyInitialize(de);CHKERRQ(ierr);
11940c453e9SDave May   for (r=0; r<nneighbors; r++) {
12040c453e9SDave May     _rank = neighbourranks[r];
12140c453e9SDave May     if ((_rank != rank) && (_rank > 0)) {
12277048351SPatrick Sanan       ierr = DMSwarmDataExTopologyAddNeighbour(de,_rank);CHKERRQ(ierr);
12340c453e9SDave May     }
12440c453e9SDave May   }
12577048351SPatrick Sanan   ierr = DMSwarmDataExTopologyFinalize(de);CHKERRQ(ierr);
12677048351SPatrick Sanan   ierr = DMSwarmDataExTopologyGetNeighbours(de,&mynneigh,&myneigh);CHKERRQ(ierr);
12777048351SPatrick Sanan   ierr = DMSwarmDataExInitializeSendCount(de);CHKERRQ(ierr);
12840c453e9SDave May   for (p=0; p<npoints; p++) {
129f954cb40SDave May     if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
1307c6d1d28SDave May       for (r=0; r<mynneigh; r++) {
1317c6d1d28SDave May         _rank = myneigh[r];
13277048351SPatrick Sanan         ierr = DMSwarmDataExAddToSendCount(de,_rank,1);CHKERRQ(ierr);
13340c453e9SDave May       }
13440c453e9SDave May     }
13540c453e9SDave May   }
13677048351SPatrick Sanan   ierr = DMSwarmDataExFinalizeSendCount(de);CHKERRQ(ierr);
13777048351SPatrick Sanan   ierr = DMSwarmDataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
13877048351SPatrick Sanan   ierr = DMSwarmDataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
13940c453e9SDave May   for (p=0; p<npoints; p++) {
140f954cb40SDave May     if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
1417c6d1d28SDave May       for (r=0; r<mynneigh; r++) {
1427c6d1d28SDave May         _rank = myneigh[r];
14340c453e9SDave May         /* copy point into buffer */
14477048351SPatrick Sanan         ierr = DMSwarmDataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr);
14577048351SPatrick Sanan         /* insert point buffer into DMSwarmDataExchanger */
14677048351SPatrick Sanan         ierr = DMSwarmDataExPackData(de,_rank,1,point_buffer);CHKERRQ(ierr);
14740c453e9SDave May       }
14840c453e9SDave May     }
14940c453e9SDave May   }
15077048351SPatrick Sanan   ierr = DMSwarmDataExPackFinalize(de);CHKERRQ(ierr);
1517c6d1d28SDave May   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
15240c453e9SDave May   if (remove_sent_points) {
15377048351SPatrick Sanan     DMSwarmDataField PField;
1547c6d1d28SDave May 
15577048351SPatrick Sanan     ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,DMSwarmField_rank,&PField);CHKERRQ(ierr);
15677048351SPatrick Sanan     ierr = DMSwarmDataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr);
15740c453e9SDave May     /* remove points which left processor */
15877048351SPatrick Sanan     ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
15940c453e9SDave May     for (p=0; p<npoints; p++) {
160f954cb40SDave May       if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
16140c453e9SDave May         /* kill point */
16277048351SPatrick Sanan         ierr = DMSwarmDataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr);
16377048351SPatrick Sanan         ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */
16477048351SPatrick Sanan         ierr = DMSwarmDataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); /* update date point increase realloc performed */
16540c453e9SDave May         p--; /* check replacement point */
16640c453e9SDave May       }
16740c453e9SDave May     }
16840c453e9SDave May   }
16977048351SPatrick Sanan   ierr = DMSwarmDataBucketGetSizes(swarm->db,npoints_prior_migration,NULL,NULL);CHKERRQ(ierr);
17077048351SPatrick Sanan   ierr = DMSwarmDataExBegin(de);CHKERRQ(ierr);
17177048351SPatrick Sanan   ierr = DMSwarmDataExEnd(de);CHKERRQ(ierr);
17277048351SPatrick Sanan   ierr = DMSwarmDataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
17377048351SPatrick Sanan   ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
17477048351SPatrick Sanan   ierr = DMSwarmDataBucketSetSizes(swarm->db,npoints + n_points_recv,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr);
17540c453e9SDave May   for (p=0; p<n_points_recv; p++) {
17640c453e9SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point);
17740c453e9SDave May 
17877048351SPatrick Sanan     ierr = DMSwarmDataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
17940c453e9SDave May   }
18077048351SPatrick Sanan   ierr = DMSwarmDataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
18177048351SPatrick Sanan   ierr = DMSwarmDataExDestroy(de);CHKERRQ(ierr);
18240c453e9SDave May   PetscFunctionReturn(0);
18340c453e9SDave May }
184480eef7bSDave May 
18508056efcSDave May PetscErrorCode DMSwarmMigrate_CellDMScatter(DM dm,PetscBool remove_sent_points)
186480eef7bSDave May {
187480eef7bSDave May   DM_Swarm          *swarm = (DM_Swarm*)dm->data;
188480eef7bSDave May   PetscErrorCode    ierr;
1896fbf25f8SDave May   PetscInt          p,npoints,npointsg=0,npoints2,npoints2g,*rankval,npoints_prior_migration;
190bbe8250bSMatthew G. Knepley   PetscSF           sfcell = NULL;
191dfc14de9SMatthew G. Knepley   const PetscSFNode *LA_sfcell;
192480eef7bSDave May   DM                dmcell;
193480eef7bSDave May   Vec               pos;
19440c453e9SDave May   PetscBool         error_check = swarm->migrate_error_on_missing_point;
195e4fbd051SBarry Smith   PetscMPIInt       size,rank;
196480eef7bSDave May 
197521f74f9SMatthew G. Knepley   PetscFunctionBegin;
198480eef7bSDave May   ierr = DMSwarmGetCellDM(dm,&dmcell);CHKERRQ(ierr);
199480eef7bSDave May   if (!dmcell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid if cell DM provided");
200480eef7bSDave May 
201ffc4695bSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRMPI(ierr);
202ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRMPI(ierr);
2037c6d1d28SDave May 
20443a82f2bSDave May #if 1
20543a82f2bSDave May   {
20643a82f2bSDave May     PetscInt *p_cellid;
20743a82f2bSDave May     PetscInt npoints_curr,range = 0;
20843a82f2bSDave May     PetscSFNode *sf_cells;
20943a82f2bSDave May 
21077048351SPatrick Sanan     ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints_curr,NULL,NULL);CHKERRQ(ierr);
21143a82f2bSDave May     ierr = PetscMalloc1(npoints_curr, &sf_cells);CHKERRQ(ierr);
21243a82f2bSDave May 
21343a82f2bSDave May     ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
21443a82f2bSDave May     ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&p_cellid);CHKERRQ(ierr);
21543a82f2bSDave May     for (p=0; p<npoints_curr; p++) {
21643a82f2bSDave May 
21743a82f2bSDave May       sf_cells[p].rank  = 0;
21843a82f2bSDave May       sf_cells[p].index = p_cellid[p];
21943a82f2bSDave May       if (p_cellid[p] > range) {
22043a82f2bSDave May         range = p_cellid[p];
22143a82f2bSDave May       }
22243a82f2bSDave May     }
22343a82f2bSDave May     ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&p_cellid);CHKERRQ(ierr);
22443a82f2bSDave May     ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
22543a82f2bSDave May 
226c28a3c15SDave May     /*ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm),&sfcell);CHKERRQ(ierr);*/
227c28a3c15SDave May     ierr = PetscSFCreate(PETSC_COMM_SELF,&sfcell);CHKERRQ(ierr);
22843a82f2bSDave May     ierr = PetscSFSetGraph(sfcell, range, npoints_curr, NULL, PETSC_OWN_POINTER, sf_cells, PETSC_OWN_POINTER);CHKERRQ(ierr);
22943a82f2bSDave May   }
23043a82f2bSDave May #endif
23143a82f2bSDave May 
232fb1bcc12SMatthew G. Knepley   ierr = DMSwarmCreateLocalVectorFromField(dm, DMSwarmPICField_coor, &pos);CHKERRQ(ierr);
233dfc14de9SMatthew G. Knepley   ierr = DMLocatePoints(dmcell, pos, DM_POINTLOCATION_NONE, &sfcell);CHKERRQ(ierr);
234fb1bcc12SMatthew G. Knepley   ierr = DMSwarmDestroyLocalVectorFromField(dm, DMSwarmPICField_coor, &pos);CHKERRQ(ierr);
235480eef7bSDave May 
23640c453e9SDave May   if (error_check) {
23740c453e9SDave May     ierr = DMSwarmGetSize(dm,&npointsg);CHKERRQ(ierr);
23840c453e9SDave May   }
23977048351SPatrick Sanan   ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
24008056efcSDave May   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
241dfc14de9SMatthew G. Knepley   ierr = PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell);CHKERRQ(ierr);
242480eef7bSDave May   for (p=0; p<npoints; p++) {
243dfc14de9SMatthew G. Knepley     rankval[p] = LA_sfcell[p].index;
244480eef7bSDave May   }
24508056efcSDave May   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
246dfc14de9SMatthew G. Knepley   ierr = PetscSFDestroy(&sfcell);CHKERRQ(ierr);
247480eef7bSDave May 
248e4fbd051SBarry Smith   if (size > 1) {
249889dbfe5SDave May     ierr = DMSwarmMigrate_DMNeighborScatter(dm,dmcell,remove_sent_points,&npoints_prior_migration);CHKERRQ(ierr);
2506fbf25f8SDave May   } else {
25177048351SPatrick Sanan     DMSwarmDataField PField;
2520ed23c7fSDave May     PetscInt npoints_curr;
2530ed23c7fSDave May 
2540ed23c7fSDave May     /* remove points which the domain */
25577048351SPatrick Sanan     ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,DMSwarmField_rank,&PField);CHKERRQ(ierr);
25677048351SPatrick Sanan     ierr = DMSwarmDataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr);
2570ed23c7fSDave May 
25877048351SPatrick Sanan     ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints_curr,NULL,NULL);CHKERRQ(ierr);
2590ed23c7fSDave May     for (p=0; p<npoints_curr; p++) {
2600ed23c7fSDave May       if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
2610ed23c7fSDave May         /* kill point */
26277048351SPatrick Sanan         ierr = DMSwarmDataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr);
26377048351SPatrick Sanan         ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints_curr,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */
26477048351SPatrick Sanan         ierr = DMSwarmDataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); /* update date point increase realloc performed */
2650ed23c7fSDave May         p--; /* check replacement point */
2660ed23c7fSDave May       }
2670ed23c7fSDave May     }
2686fbf25f8SDave May     ierr = DMSwarmGetSize(dm,&npoints_prior_migration);CHKERRQ(ierr);
2690ed23c7fSDave May 
2706fbf25f8SDave May   }
271480eef7bSDave May 
2722d4ee042Sprj-   /* locate points newly received */
27377048351SPatrick Sanan   ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr);
274009b43efSDave May 
2757c6d1d28SDave May #if 0
2762d4ee042Sprj-   { /* safe alternative - however this performs two point locations on: (i) the initial points set and; (ii) the (initial + received) point set */
2777c6d1d28SDave May     PetscScalar *LA_coor;
2787c6d1d28SDave May     PetscInt bs;
27977048351SPatrick Sanan     DMSwarmDataField PField;
2807c6d1d28SDave May 
2817c6d1d28SDave May     ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr);
2827c6d1d28SDave May     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,bs*npoints2,(const PetscScalar*)LA_coor,&pos);CHKERRQ(ierr);
283dfc14de9SMatthew G. Knepley     ierr = DMLocatePoints(dmcell,pos,DM_POINTLOCATION_NONE,&sfcell);CHKERRQ(ierr);
2847c6d1d28SDave May 
2857c6d1d28SDave May     ierr = VecDestroy(&pos);CHKERRQ(ierr);
2867c6d1d28SDave May     ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr);
2877c6d1d28SDave May 
288dfc14de9SMatthew G. Knepley     ierr = PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell);CHKERRQ(ierr);
2897c6d1d28SDave May     ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
2907c6d1d28SDave May     for (p=0; p<npoints2; p++) {
291dfc14de9SMatthew G. Knepley       rankval[p] = LA_sfcell[p].index;
2927c6d1d28SDave May     }
293dfc14de9SMatthew G. Knepley     ierr = PetscSFDestroy(&sfcell);CHKERRQ(ierr);
29408056efcSDave May     ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
29540c453e9SDave May 
2967c6d1d28SDave May     /* remove points which left processor */
29777048351SPatrick Sanan     ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,DMSwarmField_rank,&PField);CHKERRQ(ierr);
29877048351SPatrick Sanan     ierr = DMSwarmDataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr);
2997c6d1d28SDave May 
30077048351SPatrick Sanan     ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr);
3017c6d1d28SDave May     for (p=0; p<npoints2; p++) {
302f954cb40SDave May       if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
3037c6d1d28SDave May         /* kill point */
30477048351SPatrick Sanan         ierr = DMSwarmDataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr);
30577048351SPatrick Sanan         ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */
30677048351SPatrick Sanan         ierr = DMSwarmDataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); /* update date point increase realloc performed */
3077c6d1d28SDave May         p--; /* check replacement point */
3087c6d1d28SDave May       }
3097c6d1d28SDave May     }
31040c453e9SDave May   }
311009b43efSDave May #endif
312009b43efSDave May 
313*5627991aSBarry Smith   { /* perform two point locations: (i) on the initial points set prior to communication; and (ii) on the new (received) points */
314009b43efSDave May     PetscScalar      *LA_coor;
315009b43efSDave May     PetscInt         npoints_from_neighbours,bs;
31677048351SPatrick Sanan     DMSwarmDataField PField;
317009b43efSDave May 
318009b43efSDave May     npoints_from_neighbours = npoints2 - npoints_prior_migration;
319009b43efSDave May 
320009b43efSDave May     ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr);
321009b43efSDave May     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,bs*npoints_from_neighbours,(const PetscScalar*)&LA_coor[bs*npoints_prior_migration],&pos);CHKERRQ(ierr);
322009b43efSDave May 
323009b43efSDave May     ierr = DMLocatePoints(dmcell,pos,DM_POINTLOCATION_NONE,&sfcell);CHKERRQ(ierr);
324009b43efSDave May 
325009b43efSDave May     ierr = VecDestroy(&pos);CHKERRQ(ierr);
326009b43efSDave May     ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr);
327009b43efSDave May 
328009b43efSDave May     ierr = PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell);CHKERRQ(ierr);
329009b43efSDave May     ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
330009b43efSDave May     for (p=0; p<npoints_from_neighbours; p++) {
331009b43efSDave May       rankval[npoints_prior_migration + p] = LA_sfcell[p].index;
332009b43efSDave May     }
333009b43efSDave May     ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
334009b43efSDave May     ierr = PetscSFDestroy(&sfcell);CHKERRQ(ierr);
335009b43efSDave May 
336009b43efSDave May     /* remove points which left processor */
33777048351SPatrick Sanan     ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,DMSwarmField_rank,&PField);CHKERRQ(ierr);
33877048351SPatrick Sanan     ierr = DMSwarmDataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr);
339009b43efSDave May 
34077048351SPatrick Sanan     ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr);
341009b43efSDave May     for (p=npoints_prior_migration; p<npoints2; p++) {
342009b43efSDave May       if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
343009b43efSDave May         /* kill point */
34477048351SPatrick Sanan         ierr = DMSwarmDataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr);
34577048351SPatrick Sanan         ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */
34677048351SPatrick Sanan         ierr = DMSwarmDataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); /* update date point increase realloc performed */
347009b43efSDave May         p--; /* check replacement point */
348009b43efSDave May       }
349009b43efSDave May     }
350009b43efSDave May   }
351009b43efSDave May 
3523151f1d5SDave May   {
3533151f1d5SDave May     PetscInt *p_cellid;
3543151f1d5SDave May 
35577048351SPatrick Sanan     ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr);
3563151f1d5SDave May     ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
3573151f1d5SDave May     ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&p_cellid);CHKERRQ(ierr);
3583151f1d5SDave May     for (p=0; p<npoints2; p++) {
3593151f1d5SDave May       p_cellid[p] = rankval[p];
3603151f1d5SDave May     }
3613151f1d5SDave May     ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&p_cellid);CHKERRQ(ierr);
3623151f1d5SDave May     ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
3633151f1d5SDave May   }
3643151f1d5SDave May 
36540c453e9SDave May   /* check for error on removed points */
36640c453e9SDave May   if (error_check) {
36740c453e9SDave May     ierr = DMSwarmGetSize(dm,&npoints2g);CHKERRQ(ierr);
36840c453e9SDave 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);
36940c453e9SDave May   }
370480eef7bSDave May   PetscFunctionReturn(0);
371480eef7bSDave May }
372480eef7bSDave May 
37308056efcSDave May PetscErrorCode DMSwarmMigrate_CellDMExact(DM dm,PetscBool remove_sent_points)
37408056efcSDave May {
375521f74f9SMatthew G. Knepley   PetscFunctionBegin;
37608056efcSDave May   PetscFunctionReturn(0);
37708056efcSDave May }
37808056efcSDave May 
379480eef7bSDave May /*
380480eef7bSDave May  Redundant as this assumes points can only be sent to a single rank
381480eef7bSDave May */
3822712d1f2SDave May PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm,PetscInt *globalsize)
3832712d1f2SDave May {
3842712d1f2SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
3852712d1f2SDave May   PetscErrorCode ierr;
38677048351SPatrick Sanan   DMSwarmDataEx  de;
3872712d1f2SDave May   PetscInt       p,npoints,*rankval,n_points_recv;
3882712d1f2SDave May   PetscMPIInt    rank,nrank,negrank;
3892712d1f2SDave May   void           *point_buffer,*recv_points;
3902712d1f2SDave May   size_t         sizeof_dmswarm_point;
3912712d1f2SDave May 
392521f74f9SMatthew G. Knepley   PetscFunctionBegin;
393ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRMPI(ierr);
39477048351SPatrick Sanan   ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
3952712d1f2SDave May   *globalsize = npoints;
39608056efcSDave May   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
39777048351SPatrick Sanan   ierr = DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr);
39877048351SPatrick Sanan   ierr = DMSwarmDataExTopologyInitialize(de);CHKERRQ(ierr);
3992712d1f2SDave May   for (p=0; p<npoints; p++) {
4002712d1f2SDave May     negrank = rankval[p];
4012712d1f2SDave May     if (negrank < 0) {
4022712d1f2SDave May       nrank = -negrank - 1;
40377048351SPatrick Sanan       ierr = DMSwarmDataExTopologyAddNeighbour(de,nrank);CHKERRQ(ierr);
4042712d1f2SDave May     }
4052712d1f2SDave May   }
40677048351SPatrick Sanan   ierr = DMSwarmDataExTopologyFinalize(de);CHKERRQ(ierr);
40777048351SPatrick Sanan   ierr = DMSwarmDataExInitializeSendCount(de);CHKERRQ(ierr);
4082712d1f2SDave May   for (p=0; p<npoints; p++) {
4092712d1f2SDave May     negrank = rankval[p];
4102712d1f2SDave May     if (negrank < 0) {
4112712d1f2SDave May       nrank = -negrank - 1;
41277048351SPatrick Sanan       ierr = DMSwarmDataExAddToSendCount(de,nrank,1);CHKERRQ(ierr);
4132712d1f2SDave May     }
4142712d1f2SDave May   }
41577048351SPatrick Sanan   ierr = DMSwarmDataExFinalizeSendCount(de);CHKERRQ(ierr);
41677048351SPatrick Sanan   ierr = DMSwarmDataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
41777048351SPatrick Sanan   ierr = DMSwarmDataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
4182712d1f2SDave May   for (p=0; p<npoints; p++) {
4192712d1f2SDave May     negrank = rankval[p];
4202712d1f2SDave May     if (negrank < 0) {
4212712d1f2SDave May       nrank = -negrank - 1;
4222712d1f2SDave May       rankval[p] = nrank;
4232712d1f2SDave May       /* copy point into buffer */
42477048351SPatrick Sanan       ierr = DMSwarmDataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr);
42577048351SPatrick Sanan       /* insert point buffer into DMSwarmDataExchanger */
42677048351SPatrick Sanan       ierr = DMSwarmDataExPackData(de,nrank,1,point_buffer);CHKERRQ(ierr);
4272712d1f2SDave May       rankval[p] = negrank;
4282712d1f2SDave May     }
4292712d1f2SDave May   }
43077048351SPatrick Sanan   ierr = DMSwarmDataExPackFinalize(de);CHKERRQ(ierr);
43108056efcSDave May   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
43277048351SPatrick Sanan   ierr = DMSwarmDataExBegin(de);CHKERRQ(ierr);
43377048351SPatrick Sanan   ierr = DMSwarmDataExEnd(de);CHKERRQ(ierr);
43477048351SPatrick Sanan   ierr = DMSwarmDataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
43577048351SPatrick Sanan   ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
43677048351SPatrick Sanan   ierr = DMSwarmDataBucketSetSizes(swarm->db,npoints + n_points_recv,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr);
4372712d1f2SDave May   for (p=0; p<n_points_recv; p++) {
4382712d1f2SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point);
4392712d1f2SDave May 
44077048351SPatrick Sanan     ierr = DMSwarmDataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
4412712d1f2SDave May   }
44277048351SPatrick Sanan   ierr = DMSwarmDataExView(de);CHKERRQ(ierr);
44377048351SPatrick Sanan   ierr = DMSwarmDataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
44477048351SPatrick Sanan   ierr = DMSwarmDataExDestroy(de);CHKERRQ(ierr);
4452712d1f2SDave May   PetscFunctionReturn(0);
4462712d1f2SDave May }
447b16650c8SDave May 
448b16650c8SDave May typedef struct {
449b16650c8SDave May   PetscMPIInt owner_rank;
450b16650c8SDave May   PetscReal   min[3],max[3];
451b16650c8SDave May } CollectBBox;
452b16650c8SDave May 
453fe39f135SDave May PETSC_EXTERN PetscErrorCode DMSwarmCollect_DMDABoundingBox(DM dm,PetscInt *globalsize)
454b16650c8SDave May {
455b16650c8SDave May   DM_Swarm *        swarm = (DM_Swarm*)dm->data;
456b16650c8SDave May   PetscErrorCode    ierr;
45777048351SPatrick Sanan   DMSwarmDataEx     de;
458b16650c8SDave May   PetscInt          p,pk,npoints,*rankval,n_points_recv,n_bbox_recv,dim,neighbour_cells;
459b16650c8SDave May   PetscMPIInt       rank,nrank;
460b16650c8SDave May   void              *point_buffer,*recv_points;
461b16650c8SDave May   size_t            sizeof_dmswarm_point,sizeof_bbox_ctx;
462b16650c8SDave May   PetscBool         isdmda;
463b16650c8SDave May   CollectBBox       *bbox,*recv_bbox;
464b16650c8SDave May   const PetscMPIInt *dmneighborranks;
465b16650c8SDave May   DM                dmcell;
466b16650c8SDave May 
467521f74f9SMatthew G. Knepley   PetscFunctionBegin;
468ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRMPI(ierr);
469b16650c8SDave May 
470fe39f135SDave May   ierr = DMSwarmGetCellDM(dm,&dmcell);CHKERRQ(ierr);
471fe39f135SDave May   if (!dmcell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid if cell DM provided");
472b16650c8SDave May   isdmda = PETSC_FALSE;
473b16650c8SDave May   PetscObjectTypeCompare((PetscObject)dmcell,DMDA,&isdmda);
474b16650c8SDave May   if (!isdmda) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only DMDA support for CollectBoundingBox");
475b16650c8SDave May 
4768dbd68bcSDave May   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
477b16650c8SDave May   sizeof_bbox_ctx = sizeof(CollectBBox);
478b16650c8SDave May   PetscMalloc1(1,&bbox);
479b16650c8SDave May   bbox->owner_rank = rank;
480b16650c8SDave May 
481b16650c8SDave May   /* compute the bounding box based on the overlapping / stenctil size */
482b16650c8SDave May   {
483b16650c8SDave May     Vec lcoor;
484b16650c8SDave May 
485b16650c8SDave May     ierr = DMGetCoordinatesLocal(dmcell,&lcoor);CHKERRQ(ierr);
486fe39f135SDave May     if (dim >= 1) {
487b16650c8SDave May       ierr = VecStrideMin(lcoor,0,NULL,&bbox->min[0]);CHKERRQ(ierr);
488b16650c8SDave May       ierr = VecStrideMax(lcoor,0,NULL,&bbox->max[0]);CHKERRQ(ierr);
489fe39f135SDave May     }
490fe39f135SDave May     if (dim >= 2) {
491fe39f135SDave May       ierr = VecStrideMin(lcoor,1,NULL,&bbox->min[1]);CHKERRQ(ierr);
492b16650c8SDave May       ierr = VecStrideMax(lcoor,1,NULL,&bbox->max[1]);CHKERRQ(ierr);
493b16650c8SDave May     }
494fe39f135SDave May     if (dim == 3) {
495fe39f135SDave May       ierr = VecStrideMin(lcoor,2,NULL,&bbox->min[2]);CHKERRQ(ierr);
496fe39f135SDave May       ierr = VecStrideMax(lcoor,2,NULL,&bbox->max[2]);CHKERRQ(ierr);
497fe39f135SDave May     }
498fe39f135SDave May   }
49977048351SPatrick Sanan   ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
500b16650c8SDave May   *globalsize = npoints;
50108056efcSDave May   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
50277048351SPatrick Sanan   ierr = DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr);
503b16650c8SDave May   /* use DMDA neighbours */
504b16650c8SDave May   ierr = DMDAGetNeighbors(dmcell,&dmneighborranks);CHKERRQ(ierr);
5058dbd68bcSDave May   if (dim == 1) {
5068dbd68bcSDave May     neighbour_cells = 3;
5078dbd68bcSDave May   } else if (dim == 2) {
5088dbd68bcSDave May     neighbour_cells = 9;
5098dbd68bcSDave May   } else {
5108dbd68bcSDave May     neighbour_cells = 27;
5118dbd68bcSDave May   }
51277048351SPatrick Sanan   ierr = DMSwarmDataExTopologyInitialize(de);CHKERRQ(ierr);
513b16650c8SDave May   for (p=0; p<neighbour_cells; p++) {
514b16650c8SDave May     if ((dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank)) {
51577048351SPatrick Sanan       ierr = DMSwarmDataExTopologyAddNeighbour(de,dmneighborranks[p]);CHKERRQ(ierr);
516b16650c8SDave May     }
517b16650c8SDave May   }
51877048351SPatrick Sanan   ierr = DMSwarmDataExTopologyFinalize(de);CHKERRQ(ierr);
51977048351SPatrick Sanan   ierr = DMSwarmDataExInitializeSendCount(de);CHKERRQ(ierr);
520b16650c8SDave May   for (p=0; p<neighbour_cells; p++) {
521b16650c8SDave May     if ((dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank)) {
52277048351SPatrick Sanan       ierr = DMSwarmDataExAddToSendCount(de,dmneighborranks[p],1);CHKERRQ(ierr);
523b16650c8SDave May     }
524b16650c8SDave May   }
52577048351SPatrick Sanan   ierr = DMSwarmDataExFinalizeSendCount(de);CHKERRQ(ierr);
526b16650c8SDave May   /* send bounding boxes */
52777048351SPatrick Sanan   ierr = DMSwarmDataExPackInitialize(de,sizeof_bbox_ctx);CHKERRQ(ierr);
528b16650c8SDave May   for (p=0; p<neighbour_cells; p++) {
529b16650c8SDave May     nrank = dmneighborranks[p];
530b16650c8SDave May     if ((nrank >= 0) && (nrank != rank)) {
53177048351SPatrick Sanan       /* insert bbox buffer into DMSwarmDataExchanger */
53277048351SPatrick Sanan       ierr = DMSwarmDataExPackData(de,nrank,1,bbox);CHKERRQ(ierr);
533b16650c8SDave May     }
534b16650c8SDave May   }
53577048351SPatrick Sanan   ierr = DMSwarmDataExPackFinalize(de);CHKERRQ(ierr);
536b16650c8SDave May   /* recv bounding boxes */
53777048351SPatrick Sanan   ierr = DMSwarmDataExBegin(de);CHKERRQ(ierr);
53877048351SPatrick Sanan   ierr = DMSwarmDataExEnd(de);CHKERRQ(ierr);
53977048351SPatrick Sanan   ierr = DMSwarmDataExGetRecvData(de,&n_bbox_recv,(void**)&recv_bbox);CHKERRQ(ierr);
540298827fbSBarry Smith   /*  Wrong, should not be using PETSC_COMM_WORLD */
541b16650c8SDave May   for (p=0; p<n_bbox_recv; p++) {
542298827fbSBarry Smith     ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[rank %d]: box from %d : range[%+1.4e,%+1.4e]x[%+1.4e,%+1.4e]\n",rank,recv_bbox[p].owner_rank,
543298827fbSBarry Smith            (double)recv_bbox[p].min[0],(double)recv_bbox[p].max[0],(double)recv_bbox[p].min[1],(double)recv_bbox[p].max[1]);CHKERRQ(ierr);
544b16650c8SDave May   }
545298827fbSBarry Smith   ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,stdout);CHKERRQ(ierr);
546b16650c8SDave May   /* of course this is stupid as this "generic" function should have a better way to know what the coordinates are called */
54777048351SPatrick Sanan   ierr = DMSwarmDataExInitializeSendCount(de);CHKERRQ(ierr);
548b16650c8SDave May   for (pk=0; pk<n_bbox_recv; pk++) {
549b16650c8SDave May     PetscReal *array_x,*array_y;
550b16650c8SDave May 
551b16650c8SDave May     ierr = DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
552b16650c8SDave May     ierr = DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
553b16650c8SDave May     for (p=0; p<npoints; p++) {
554b16650c8SDave May       if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0])) {
555b16650c8SDave May         if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1])) {
55677048351SPatrick Sanan           ierr = DMSwarmDataExAddToSendCount(de,recv_bbox[pk].owner_rank,1);CHKERRQ(ierr);
557b16650c8SDave May         }
558b16650c8SDave May       }
559b16650c8SDave May     }
560b16650c8SDave May     ierr = DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
561b16650c8SDave May     ierr = DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
562b16650c8SDave May   }
56377048351SPatrick Sanan   ierr = DMSwarmDataExFinalizeSendCount(de);CHKERRQ(ierr);
56477048351SPatrick Sanan   ierr = DMSwarmDataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
56577048351SPatrick Sanan   ierr = DMSwarmDataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
566b16650c8SDave May   for (pk=0; pk<n_bbox_recv; pk++) {
567b16650c8SDave May     PetscReal *array_x,*array_y;
568b16650c8SDave May 
569b16650c8SDave May     ierr = DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
570b16650c8SDave May     ierr = DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
571b16650c8SDave May     for (p=0; p<npoints; p++) {
572b16650c8SDave May       if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0])) {
573b16650c8SDave May         if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1])) {
574521f74f9SMatthew G. Knepley           /* copy point into buffer */
57577048351SPatrick Sanan           ierr = DMSwarmDataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr);
57677048351SPatrick Sanan           /* insert point buffer into DMSwarmDataExchanger */
57777048351SPatrick Sanan           ierr = DMSwarmDataExPackData(de,recv_bbox[pk].owner_rank,1,point_buffer);CHKERRQ(ierr);
578b16650c8SDave May         }
579b16650c8SDave May       }
580b16650c8SDave May     }
581b16650c8SDave May     ierr = DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
582b16650c8SDave May     ierr = DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
583b16650c8SDave May   }
58477048351SPatrick Sanan   ierr = DMSwarmDataExPackFinalize(de);CHKERRQ(ierr);
58508056efcSDave May   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
58677048351SPatrick Sanan   ierr = DMSwarmDataExBegin(de);CHKERRQ(ierr);
58777048351SPatrick Sanan   ierr = DMSwarmDataExEnd(de);CHKERRQ(ierr);
58877048351SPatrick Sanan   ierr = DMSwarmDataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
58977048351SPatrick Sanan   ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
59077048351SPatrick Sanan   ierr = DMSwarmDataBucketSetSizes(swarm->db,npoints + n_points_recv,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr);
591b16650c8SDave May   for (p=0; p<n_points_recv; p++) {
592b16650c8SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point);
593b16650c8SDave May 
59477048351SPatrick Sanan     ierr = DMSwarmDataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
595b16650c8SDave May   }
59677048351SPatrick Sanan   ierr = DMSwarmDataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
597b16650c8SDave May   PetscFree(bbox);
59877048351SPatrick Sanan   ierr = DMSwarmDataExView(de);CHKERRQ(ierr);
59977048351SPatrick Sanan   ierr = DMSwarmDataExDestroy(de);CHKERRQ(ierr);
600b16650c8SDave May   PetscFunctionReturn(0);
601b16650c8SDave May }
602a9fd7477SDave May 
603a9fd7477SDave May /* General collection when no order, or neighbour information is provided */
604a9fd7477SDave May /*
605a9fd7477SDave May  User provides context and collect() method
606a9fd7477SDave May  Broadcast user context
607a9fd7477SDave May 
608a9fd7477SDave May  for each context / rank {
609a9fd7477SDave May    collect(swarm,context,n,list)
610a9fd7477SDave May  }
611a9fd7477SDave May */
612a9fd7477SDave May PETSC_EXTERN PetscErrorCode DMSwarmCollect_General(DM dm,PetscErrorCode (*collect)(DM,void*,PetscInt*,PetscInt**),size_t ctx_size,void *ctx,PetscInt *globalsize)
613a9fd7477SDave May {
614a9fd7477SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
615a9fd7477SDave May   PetscErrorCode ierr;
61677048351SPatrick Sanan   DMSwarmDataEx  de;
617a9fd7477SDave May   PetscInt       p,r,npoints,n_points_recv;
618e4fbd051SBarry Smith   PetscMPIInt    size,rank;
619a9fd7477SDave May   void           *point_buffer,*recv_points;
620a9fd7477SDave May   void           *ctxlist;
621a9fd7477SDave May   PetscInt       *n2collect,**collectlist;
622a9fd7477SDave May   size_t         sizeof_dmswarm_point;
623a9fd7477SDave May 
624521f74f9SMatthew G. Knepley   PetscFunctionBegin;
625ffc4695bSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRMPI(ierr);
626ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRMPI(ierr);
62777048351SPatrick Sanan   ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
628a9fd7477SDave May   *globalsize = npoints;
629a9fd7477SDave May   /* Broadcast user context */
630e4fbd051SBarry Smith   PetscMalloc(ctx_size*size,&ctxlist);
631ffc4695bSBarry Smith   ierr = MPI_Allgather(ctx,ctx_size,MPI_CHAR,ctxlist,ctx_size,MPI_CHAR,PetscObjectComm((PetscObject)dm));CHKERRMPI(ierr);
632e4fbd051SBarry Smith   ierr = PetscMalloc1(size,&n2collect);CHKERRQ(ierr);
633e4fbd051SBarry Smith   ierr = PetscMalloc1(size,&collectlist);CHKERRQ(ierr);
634e4fbd051SBarry Smith   for (r=0; r<size; r++) {
635a9fd7477SDave May     PetscInt _n2collect;
636a9fd7477SDave May     PetscInt *_collectlist;
637a9fd7477SDave May     void     *_ctx_r;
638a9fd7477SDave May 
639a9fd7477SDave May     _n2collect   = 0;
640a9fd7477SDave May     _collectlist = NULL;
641a9fd7477SDave May     if (r != rank) { /* don't collect data from yourself */
642a9fd7477SDave May       _ctx_r = (void*)( (char*)ctxlist + r * ctx_size);
643a9fd7477SDave May       ierr = collect(dm,_ctx_r,&_n2collect,&_collectlist);CHKERRQ(ierr);
644a9fd7477SDave May     }
645a9fd7477SDave May     n2collect[r]   = _n2collect;
646a9fd7477SDave May     collectlist[r] = _collectlist;
647a9fd7477SDave May   }
64877048351SPatrick Sanan   ierr = DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr);
649a9fd7477SDave May   /* Define topology */
65077048351SPatrick Sanan   ierr = DMSwarmDataExTopologyInitialize(de);CHKERRQ(ierr);
651e4fbd051SBarry Smith   for (r=0; r<size; r++) {
652a9fd7477SDave May     if (n2collect[r] > 0) {
65377048351SPatrick Sanan       ierr = DMSwarmDataExTopologyAddNeighbour(de,(PetscMPIInt)r);CHKERRQ(ierr);
654a9fd7477SDave May     }
655a9fd7477SDave May   }
65677048351SPatrick Sanan   ierr = DMSwarmDataExTopologyFinalize(de);CHKERRQ(ierr);
657a9fd7477SDave May   /* Define send counts */
65877048351SPatrick Sanan   ierr = DMSwarmDataExInitializeSendCount(de);CHKERRQ(ierr);
659e4fbd051SBarry Smith   for (r=0; r<size; r++) {
660a9fd7477SDave May     if (n2collect[r] > 0) {
66177048351SPatrick Sanan       ierr = DMSwarmDataExAddToSendCount(de,r,n2collect[r]);CHKERRQ(ierr);
662a9fd7477SDave May     }
663a9fd7477SDave May   }
66477048351SPatrick Sanan   ierr = DMSwarmDataExFinalizeSendCount(de);CHKERRQ(ierr);
665a9fd7477SDave May   /* Pack data */
66677048351SPatrick Sanan   ierr = DMSwarmDataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
66777048351SPatrick Sanan   ierr = DMSwarmDataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
668e4fbd051SBarry Smith   for (r=0; r<size; r++) {
669a9fd7477SDave May     for (p=0; p<n2collect[r]; p++) {
67077048351SPatrick Sanan       ierr = DMSwarmDataBucketFillPackedArray(swarm->db,collectlist[r][p],point_buffer);CHKERRQ(ierr);
671a9fd7477SDave May       /* insert point buffer into the data exchanger */
67277048351SPatrick Sanan       ierr = DMSwarmDataExPackData(de,r,1,point_buffer);CHKERRQ(ierr);
673a9fd7477SDave May     }
674a9fd7477SDave May   }
67577048351SPatrick Sanan   ierr = DMSwarmDataExPackFinalize(de);CHKERRQ(ierr);
676a9fd7477SDave May   /* Scatter */
67777048351SPatrick Sanan   ierr = DMSwarmDataExBegin(de);CHKERRQ(ierr);
67877048351SPatrick Sanan   ierr = DMSwarmDataExEnd(de);CHKERRQ(ierr);
679a9fd7477SDave May   /* Collect data in DMSwarm container */
68077048351SPatrick Sanan   ierr = DMSwarmDataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
68177048351SPatrick Sanan   ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
68277048351SPatrick Sanan   ierr = DMSwarmDataBucketSetSizes(swarm->db,npoints + n_points_recv,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr);
683a9fd7477SDave May   for (p=0; p<n_points_recv; p++) {
684a9fd7477SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point);
685a9fd7477SDave May 
68677048351SPatrick Sanan     ierr = DMSwarmDataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
687a9fd7477SDave May   }
688a9fd7477SDave May   /* Release memory */
689e4fbd051SBarry Smith   for (r=0; r<size; r++) {
690a9fd7477SDave May     if (collectlist[r]) PetscFree(collectlist[r]);
691a9fd7477SDave May   }
692521f74f9SMatthew G. Knepley   ierr = PetscFree(collectlist);CHKERRQ(ierr);
693521f74f9SMatthew G. Knepley   ierr = PetscFree(n2collect);CHKERRQ(ierr);
694521f74f9SMatthew G. Knepley   ierr = PetscFree(ctxlist);CHKERRQ(ierr);
69577048351SPatrick Sanan   ierr = DMSwarmDataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
69677048351SPatrick Sanan   ierr = DMSwarmDataExView(de);CHKERRQ(ierr);
69777048351SPatrick Sanan   ierr = DMSwarmDataExDestroy(de);CHKERRQ(ierr);
698a9fd7477SDave May   PetscFunctionReturn(0);
699a9fd7477SDave May }
700