xref: /petsc/src/dm/impls/swarm/swarm_migrate.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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;
1477048351SPatrick Sanan   DMSwarmDataEx  de;
15df21e3a8SDave May   PetscInt       p,npoints,*rankval,n_points_recv;
16df21e3a8SDave May   PetscMPIInt    rank,nrank;
17df21e3a8SDave May   void           *point_buffer,*recv_points;
18df21e3a8SDave May   size_t         sizeof_dmswarm_point;
19df21e3a8SDave May 
20521f74f9SMatthew G. Knepley   PetscFunctionBegin;
21*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank));
22df21e3a8SDave May 
23*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL));
24*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval));
25*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm),0, &de));
26*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExTopologyInitialize(de));
27521f74f9SMatthew G. Knepley   for (p = 0; p < npoints; ++p) {
28df21e3a8SDave May     nrank = rankval[p];
29df21e3a8SDave May     if (nrank != rank) {
30*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSwarmDataExTopologyAddNeighbour(de,nrank));
31df21e3a8SDave May     }
32df21e3a8SDave May   }
33*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExTopologyFinalize(de));
34*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExInitializeSendCount(de));
35df21e3a8SDave May   for (p=0; p<npoints; p++) {
36df21e3a8SDave May     nrank = rankval[p];
37df21e3a8SDave May     if (nrank != rank) {
38*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSwarmDataExAddToSendCount(de,nrank,1));
39df21e3a8SDave May     }
40df21e3a8SDave May   }
41*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExFinalizeSendCount(de));
42*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer));
43*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExPackInitialize(de,sizeof_dmswarm_point));
44df21e3a8SDave May   for (p=0; p<npoints; p++) {
45df21e3a8SDave May     nrank = rankval[p];
46df21e3a8SDave May     if (nrank != rank) {
47df21e3a8SDave May       /* copy point into buffer */
48*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSwarmDataBucketFillPackedArray(swarm->db,p,point_buffer));
4977048351SPatrick Sanan       /* insert point buffer into DMSwarmDataExchanger */
50*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSwarmDataExPackData(de,nrank,1,point_buffer));
51df21e3a8SDave May     }
52df21e3a8SDave May   }
53*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExPackFinalize(de));
54*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval));
55df21e3a8SDave May 
56df21e3a8SDave May   if (remove_sent_points) {
5777048351SPatrick Sanan     DMSwarmDataField gfield;
5822a417f9SDave May 
59*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,DMSwarmField_rank,&gfield));
60*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataFieldGetAccess(gfield));
61*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataFieldGetEntries(gfield,(void**)&rankval));
6222a417f9SDave May 
63df21e3a8SDave May     /* remove points which left processor */
64*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL));
65df21e3a8SDave May     for (p=0; p<npoints; p++) {
66df21e3a8SDave May       nrank = rankval[p];
67df21e3a8SDave May       if (nrank != rank) {
68df21e3a8SDave May         /* kill point */
69*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMSwarmDataFieldRestoreAccess(gfield));
7022a417f9SDave May 
71*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMSwarmDataBucketRemovePointAtIndex(swarm->db,p));
7222a417f9SDave May 
73*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL)); /* you need to update npoints as the list size decreases! */
74*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMSwarmDataFieldGetAccess(gfield));
75*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMSwarmDataFieldGetEntries(gfield,(void**)&rankval));
76df21e3a8SDave May         p--; /* check replacement point */
77df21e3a8SDave May       }
78df21e3a8SDave May     }
79*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataFieldRestoreEntries(gfield,(void**)&rankval));
80*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataFieldRestoreAccess(gfield));
81df21e3a8SDave May   }
82*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExBegin(de));
83*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExEnd(de));
84*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExGetRecvData(de,&n_points_recv,(void**)&recv_points));
85*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL));
86*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketSetSizes(swarm->db,npoints + n_points_recv,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
87df21e3a8SDave May   for (p=0; p<n_points_recv; p++) {
88df21e3a8SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point);
89df21e3a8SDave May 
90*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataBucketInsertPackedArray(swarm->db,npoints+p,data_p));
91df21e3a8SDave May   }
92*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExView(de));
93*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketDestroyPackedArray(swarm->db,&point_buffer));
94*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExDestroy(de));
95df21e3a8SDave May   PetscFunctionReturn(0);
96df21e3a8SDave May }
972712d1f2SDave May 
98889dbfe5SDave May PetscErrorCode DMSwarmMigrate_DMNeighborScatter(DM dm,DM dmcell,PetscBool remove_sent_points,PetscInt *npoints_prior_migration)
9940c453e9SDave May {
10040c453e9SDave May   DM_Swarm          *swarm = (DM_Swarm*)dm->data;
10177048351SPatrick Sanan   DMSwarmDataEx     de;
10240c453e9SDave May   PetscInt          r,p,npoints,*rankval,n_points_recv;
10340c453e9SDave May   PetscMPIInt       rank,_rank;
10440c453e9SDave May   const PetscMPIInt *neighbourranks;
10540c453e9SDave May   void              *point_buffer,*recv_points;
10640c453e9SDave May   size_t            sizeof_dmswarm_point;
10740c453e9SDave May   PetscInt          nneighbors;
1087c6d1d28SDave May   PetscMPIInt       mynneigh,*myneigh;
10940c453e9SDave May 
110521f74f9SMatthew G. Knepley   PetscFunctionBegin;
111*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank));
112*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL));
113*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval));
114*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm),0,&de));
115*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetNeighbors(dmcell,&nneighbors,&neighbourranks));
116*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExTopologyInitialize(de));
11740c453e9SDave May   for (r=0; r<nneighbors; r++) {
11840c453e9SDave May     _rank = neighbourranks[r];
11940c453e9SDave May     if ((_rank != rank) && (_rank > 0)) {
120*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSwarmDataExTopologyAddNeighbour(de,_rank));
12140c453e9SDave May     }
12240c453e9SDave May   }
123*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExTopologyFinalize(de));
124*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExTopologyGetNeighbours(de,&mynneigh,&myneigh));
125*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExInitializeSendCount(de));
12640c453e9SDave May   for (p=0; p<npoints; p++) {
127f954cb40SDave May     if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
1287c6d1d28SDave May       for (r=0; r<mynneigh; r++) {
1297c6d1d28SDave May         _rank = myneigh[r];
130*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMSwarmDataExAddToSendCount(de,_rank,1));
13140c453e9SDave May       }
13240c453e9SDave May     }
13340c453e9SDave May   }
134*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExFinalizeSendCount(de));
135*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer));
136*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExPackInitialize(de,sizeof_dmswarm_point));
13740c453e9SDave May   for (p=0; p<npoints; p++) {
138f954cb40SDave May     if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
1397c6d1d28SDave May       for (r=0; r<mynneigh; r++) {
1407c6d1d28SDave May         _rank = myneigh[r];
14140c453e9SDave May         /* copy point into buffer */
142*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMSwarmDataBucketFillPackedArray(swarm->db,p,point_buffer));
14377048351SPatrick Sanan         /* insert point buffer into DMSwarmDataExchanger */
144*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMSwarmDataExPackData(de,_rank,1,point_buffer));
14540c453e9SDave May       }
14640c453e9SDave May     }
14740c453e9SDave May   }
148*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExPackFinalize(de));
149*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval));
15040c453e9SDave May   if (remove_sent_points) {
15177048351SPatrick Sanan     DMSwarmDataField PField;
1527c6d1d28SDave May 
153*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,DMSwarmField_rank,&PField));
154*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataFieldGetEntries(PField,(void**)&rankval));
15540c453e9SDave May     /* remove points which left processor */
156*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL));
15740c453e9SDave May     for (p=0; p<npoints; p++) {
158f954cb40SDave May       if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
15940c453e9SDave May         /* kill point */
160*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMSwarmDataBucketRemovePointAtIndex(swarm->db,p));
161*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL)); /* you need to update npoints as the list size decreases! */
162*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMSwarmDataFieldGetEntries(PField,(void**)&rankval)); /* update date point increase realloc performed */
16340c453e9SDave May         p--; /* check replacement point */
16440c453e9SDave May       }
16540c453e9SDave May     }
16640c453e9SDave May   }
167*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db,npoints_prior_migration,NULL,NULL));
168*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExBegin(de));
169*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExEnd(de));
170*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExGetRecvData(de,&n_points_recv,(void**)&recv_points));
171*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL));
172*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketSetSizes(swarm->db,npoints + n_points_recv,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
17340c453e9SDave May   for (p=0; p<n_points_recv; p++) {
17440c453e9SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point);
17540c453e9SDave May 
176*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataBucketInsertPackedArray(swarm->db,npoints+p,data_p));
17740c453e9SDave May   }
178*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketDestroyPackedArray(swarm->db,&point_buffer));
179*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExDestroy(de));
18040c453e9SDave May   PetscFunctionReturn(0);
18140c453e9SDave May }
182480eef7bSDave May 
18308056efcSDave May PetscErrorCode DMSwarmMigrate_CellDMScatter(DM dm,PetscBool remove_sent_points)
184480eef7bSDave May {
185480eef7bSDave May   DM_Swarm          *swarm = (DM_Swarm*)dm->data;
1866fbf25f8SDave May   PetscInt          p,npoints,npointsg=0,npoints2,npoints2g,*rankval,npoints_prior_migration;
187bbe8250bSMatthew G. Knepley   PetscSF           sfcell = NULL;
188dfc14de9SMatthew G. Knepley   const PetscSFNode *LA_sfcell;
189480eef7bSDave May   DM                dmcell;
190480eef7bSDave May   Vec               pos;
19140c453e9SDave May   PetscBool         error_check = swarm->migrate_error_on_missing_point;
192e4fbd051SBarry Smith   PetscMPIInt       size,rank;
193480eef7bSDave May 
194521f74f9SMatthew G. Knepley   PetscFunctionBegin;
195*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetCellDM(dm,&dmcell));
1962c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!dmcell,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid if cell DM provided");
197480eef7bSDave May 
198*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size));
199*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank));
2007c6d1d28SDave May 
20143a82f2bSDave May #if 1
20243a82f2bSDave May   {
20343a82f2bSDave May     PetscInt *p_cellid;
20443a82f2bSDave May     PetscInt npoints_curr,range = 0;
20543a82f2bSDave May     PetscSFNode *sf_cells;
20643a82f2bSDave May 
207*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db,&npoints_curr,NULL,NULL));
208*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(npoints_curr, &sf_cells));
20943a82f2bSDave May 
210*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval));
211*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&p_cellid));
21243a82f2bSDave May     for (p=0; p<npoints_curr; p++) {
21343a82f2bSDave May 
21443a82f2bSDave May       sf_cells[p].rank  = 0;
21543a82f2bSDave May       sf_cells[p].index = p_cellid[p];
21643a82f2bSDave May       if (p_cellid[p] > range) {
21743a82f2bSDave May         range = p_cellid[p];
21843a82f2bSDave May       }
21943a82f2bSDave May     }
220*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&p_cellid));
221*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval));
22243a82f2bSDave May 
223*5f80ce2aSJacob Faibussowitsch     /*CHKERRQ(PetscSFCreate(PetscObjectComm((PetscObject)dm),&sfcell));*/
224*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFCreate(PETSC_COMM_SELF,&sfcell));
225*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFSetGraph(sfcell, range, npoints_curr, NULL, PETSC_OWN_POINTER, sf_cells, PETSC_OWN_POINTER));
22643a82f2bSDave May   }
22743a82f2bSDave May #endif
22843a82f2bSDave May 
229*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmCreateLocalVectorFromField(dm, DMSwarmPICField_coor, &pos));
230*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLocatePoints(dmcell, pos, DM_POINTLOCATION_NONE, &sfcell));
231*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDestroyLocalVectorFromField(dm, DMSwarmPICField_coor, &pos));
232480eef7bSDave May 
23340c453e9SDave May   if (error_check) {
234*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmGetSize(dm,&npointsg));
23540c453e9SDave May   }
236*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL));
237*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval));
238*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell));
239480eef7bSDave May   for (p=0; p<npoints; p++) {
240dfc14de9SMatthew G. Knepley     rankval[p] = LA_sfcell[p].index;
241480eef7bSDave May   }
242*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval));
243*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFDestroy(&sfcell));
244480eef7bSDave May 
245e4fbd051SBarry Smith   if (size > 1) {
246*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmMigrate_DMNeighborScatter(dm,dmcell,remove_sent_points,&npoints_prior_migration));
2476fbf25f8SDave May   } else {
24877048351SPatrick Sanan     DMSwarmDataField PField;
2490ed23c7fSDave May     PetscInt npoints_curr;
2500ed23c7fSDave May 
2510ed23c7fSDave May     /* remove points which the domain */
252*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,DMSwarmField_rank,&PField));
253*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataFieldGetEntries(PField,(void**)&rankval));
2540ed23c7fSDave May 
255*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db,&npoints_curr,NULL,NULL));
2560ed23c7fSDave May     for (p=0; p<npoints_curr; p++) {
2570ed23c7fSDave May       if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
2580ed23c7fSDave May         /* kill point */
259*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMSwarmDataBucketRemovePointAtIndex(swarm->db,p));
260*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db,&npoints_curr,NULL,NULL)); /* you need to update npoints as the list size decreases! */
261*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMSwarmDataFieldGetEntries(PField,(void**)&rankval)); /* update date point increase realloc performed */
2620ed23c7fSDave May         p--; /* check replacement point */
2630ed23c7fSDave May       }
2640ed23c7fSDave May     }
265*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmGetSize(dm,&npoints_prior_migration));
2660ed23c7fSDave May 
2676fbf25f8SDave May   }
268480eef7bSDave May 
2692d4ee042Sprj-   /* locate points newly received */
270*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db,&npoints2,NULL,NULL));
271009b43efSDave May 
2727c6d1d28SDave May #if 0
2732d4ee042Sprj-   { /* safe alternative - however this performs two point locations on: (i) the initial points set and; (ii) the (initial + received) point set */
2747c6d1d28SDave May     PetscScalar *LA_coor;
2757c6d1d28SDave May     PetscInt bs;
27677048351SPatrick Sanan     DMSwarmDataField PField;
2777c6d1d28SDave May 
278*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmGetField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor));
279*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreateSeqWithArray(PETSC_COMM_SELF,bs,bs*npoints2,(const PetscScalar*)LA_coor,&pos));
280*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLocatePoints(dmcell,pos,DM_POINTLOCATION_NONE,&sfcell));
2817c6d1d28SDave May 
282*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&pos));
283*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmRestoreField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor));
2847c6d1d28SDave May 
285*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell));
286*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval));
2877c6d1d28SDave May     for (p=0; p<npoints2; p++) {
288dfc14de9SMatthew G. Knepley       rankval[p] = LA_sfcell[p].index;
2897c6d1d28SDave May     }
290*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFDestroy(&sfcell));
291*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval));
29240c453e9SDave May 
2937c6d1d28SDave May     /* remove points which left processor */
294*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,DMSwarmField_rank,&PField));
295*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataFieldGetEntries(PField,(void**)&rankval));
2967c6d1d28SDave May 
297*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db,&npoints2,NULL,NULL));
2987c6d1d28SDave May     for (p=0; p<npoints2; p++) {
299f954cb40SDave May       if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
3007c6d1d28SDave May         /* kill point */
301*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMSwarmDataBucketRemovePointAtIndex(swarm->db,p));
302*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db,&npoints2,NULL,NULL)); /* you need to update npoints as the list size decreases! */
303*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMSwarmDataFieldGetEntries(PField,(void**)&rankval)); /* update date point increase realloc performed */
3047c6d1d28SDave May         p--; /* check replacement point */
3057c6d1d28SDave May       }
3067c6d1d28SDave May     }
30740c453e9SDave May   }
308009b43efSDave May #endif
309009b43efSDave May 
3105627991aSBarry Smith   { /* perform two point locations: (i) on the initial points set prior to communication; and (ii) on the new (received) points */
311009b43efSDave May     PetscScalar      *LA_coor;
312009b43efSDave May     PetscInt         npoints_from_neighbours,bs;
31377048351SPatrick Sanan     DMSwarmDataField PField;
314009b43efSDave May 
315009b43efSDave May     npoints_from_neighbours = npoints2 - npoints_prior_migration;
316009b43efSDave May 
317*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmGetField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor));
318*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreateSeqWithArray(PETSC_COMM_SELF,bs,bs*npoints_from_neighbours,(const PetscScalar*)&LA_coor[bs*npoints_prior_migration],&pos));
319009b43efSDave May 
320*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLocatePoints(dmcell,pos,DM_POINTLOCATION_NONE,&sfcell));
321009b43efSDave May 
322*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&pos));
323*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmRestoreField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor));
324009b43efSDave May 
325*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell));
326*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval));
327009b43efSDave May     for (p=0; p<npoints_from_neighbours; p++) {
328009b43efSDave May       rankval[npoints_prior_migration + p] = LA_sfcell[p].index;
329009b43efSDave May     }
330*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval));
331*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFDestroy(&sfcell));
332009b43efSDave May 
333009b43efSDave May     /* remove points which left processor */
334*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,DMSwarmField_rank,&PField));
335*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataFieldGetEntries(PField,(void**)&rankval));
336009b43efSDave May 
337*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db,&npoints2,NULL,NULL));
338009b43efSDave May     for (p=npoints_prior_migration; p<npoints2; p++) {
339009b43efSDave May       if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
340009b43efSDave May         /* kill point */
341*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMSwarmDataBucketRemovePointAtIndex(swarm->db,p));
342*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db,&npoints2,NULL,NULL)); /* you need to update npoints as the list size decreases! */
343*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMSwarmDataFieldGetEntries(PField,(void**)&rankval)); /* update date point increase realloc performed */
344009b43efSDave May         p--; /* check replacement point */
345009b43efSDave May       }
346009b43efSDave May     }
347009b43efSDave May   }
348009b43efSDave May 
3493151f1d5SDave May   {
3503151f1d5SDave May     PetscInt *p_cellid;
3513151f1d5SDave May 
352*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db,&npoints2,NULL,NULL));
353*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval));
354*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&p_cellid));
3553151f1d5SDave May     for (p=0; p<npoints2; p++) {
3563151f1d5SDave May       p_cellid[p] = rankval[p];
3573151f1d5SDave May     }
358*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&p_cellid));
359*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval));
3603151f1d5SDave May   }
3613151f1d5SDave May 
36240c453e9SDave May   /* check for error on removed points */
36340c453e9SDave May   if (error_check) {
364*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmGetSize(dm,&npoints2g));
3652c71b3e2SJacob Faibussowitsch     PetscCheckFalse(npointsg != npoints2g,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Points from the DMSwarm must remain constant during migration (initial %D - final %D)",npointsg,npoints2g);
36640c453e9SDave May   }
367480eef7bSDave May   PetscFunctionReturn(0);
368480eef7bSDave May }
369480eef7bSDave May 
37008056efcSDave May PetscErrorCode DMSwarmMigrate_CellDMExact(DM dm,PetscBool remove_sent_points)
37108056efcSDave May {
372521f74f9SMatthew G. Knepley   PetscFunctionBegin;
37308056efcSDave May   PetscFunctionReturn(0);
37408056efcSDave May }
37508056efcSDave May 
376480eef7bSDave May /*
377480eef7bSDave May  Redundant as this assumes points can only be sent to a single rank
378480eef7bSDave May */
3792712d1f2SDave May PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm,PetscInt *globalsize)
3802712d1f2SDave May {
3812712d1f2SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
38277048351SPatrick Sanan   DMSwarmDataEx  de;
3832712d1f2SDave May   PetscInt       p,npoints,*rankval,n_points_recv;
3842712d1f2SDave May   PetscMPIInt    rank,nrank,negrank;
3852712d1f2SDave May   void           *point_buffer,*recv_points;
3862712d1f2SDave May   size_t         sizeof_dmswarm_point;
3872712d1f2SDave May 
388521f74f9SMatthew G. Knepley   PetscFunctionBegin;
389*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank));
390*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL));
3912712d1f2SDave May   *globalsize = npoints;
392*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval));
393*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm),0,&de));
394*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExTopologyInitialize(de));
3952712d1f2SDave May   for (p=0; p<npoints; p++) {
3962712d1f2SDave May     negrank = rankval[p];
3972712d1f2SDave May     if (negrank < 0) {
3982712d1f2SDave May       nrank = -negrank - 1;
399*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSwarmDataExTopologyAddNeighbour(de,nrank));
4002712d1f2SDave May     }
4012712d1f2SDave May   }
402*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExTopologyFinalize(de));
403*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExInitializeSendCount(de));
4042712d1f2SDave May   for (p=0; p<npoints; p++) {
4052712d1f2SDave May     negrank = rankval[p];
4062712d1f2SDave May     if (negrank < 0) {
4072712d1f2SDave May       nrank = -negrank - 1;
408*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSwarmDataExAddToSendCount(de,nrank,1));
4092712d1f2SDave May     }
4102712d1f2SDave May   }
411*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExFinalizeSendCount(de));
412*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer));
413*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExPackInitialize(de,sizeof_dmswarm_point));
4142712d1f2SDave May   for (p=0; p<npoints; p++) {
4152712d1f2SDave May     negrank = rankval[p];
4162712d1f2SDave May     if (negrank < 0) {
4172712d1f2SDave May       nrank = -negrank - 1;
4182712d1f2SDave May       rankval[p] = nrank;
4192712d1f2SDave May       /* copy point into buffer */
420*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSwarmDataBucketFillPackedArray(swarm->db,p,point_buffer));
42177048351SPatrick Sanan       /* insert point buffer into DMSwarmDataExchanger */
422*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSwarmDataExPackData(de,nrank,1,point_buffer));
4232712d1f2SDave May       rankval[p] = negrank;
4242712d1f2SDave May     }
4252712d1f2SDave May   }
426*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExPackFinalize(de));
427*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval));
428*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExBegin(de));
429*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExEnd(de));
430*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExGetRecvData(de,&n_points_recv,(void**)&recv_points));
431*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL));
432*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketSetSizes(swarm->db,npoints + n_points_recv,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
4332712d1f2SDave May   for (p=0; p<n_points_recv; p++) {
4342712d1f2SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point);
4352712d1f2SDave May 
436*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataBucketInsertPackedArray(swarm->db,npoints+p,data_p));
4372712d1f2SDave May   }
438*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExView(de));
439*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketDestroyPackedArray(swarm->db,&point_buffer));
440*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExDestroy(de));
4412712d1f2SDave May   PetscFunctionReturn(0);
4422712d1f2SDave May }
443b16650c8SDave May 
444b16650c8SDave May typedef struct {
445b16650c8SDave May   PetscMPIInt owner_rank;
446b16650c8SDave May   PetscReal   min[3],max[3];
447b16650c8SDave May } CollectBBox;
448b16650c8SDave May 
449fe39f135SDave May PETSC_EXTERN PetscErrorCode DMSwarmCollect_DMDABoundingBox(DM dm,PetscInt *globalsize)
450b16650c8SDave May {
451b16650c8SDave May   DM_Swarm *        swarm = (DM_Swarm*)dm->data;
452b16650c8SDave May   PetscErrorCode    ierr;
45377048351SPatrick Sanan   DMSwarmDataEx     de;
454b16650c8SDave May   PetscInt          p,pk,npoints,*rankval,n_points_recv,n_bbox_recv,dim,neighbour_cells;
455b16650c8SDave May   PetscMPIInt       rank,nrank;
456b16650c8SDave May   void              *point_buffer,*recv_points;
457b16650c8SDave May   size_t            sizeof_dmswarm_point,sizeof_bbox_ctx;
458b16650c8SDave May   PetscBool         isdmda;
459b16650c8SDave May   CollectBBox       *bbox,*recv_bbox;
460b16650c8SDave May   const PetscMPIInt *dmneighborranks;
461b16650c8SDave May   DM                dmcell;
462b16650c8SDave May 
463521f74f9SMatthew G. Knepley   PetscFunctionBegin;
464*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank));
465b16650c8SDave May 
466*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetCellDM(dm,&dmcell));
4672c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!dmcell,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid if cell DM provided");
468b16650c8SDave May   isdmda = PETSC_FALSE;
469b16650c8SDave May   PetscObjectTypeCompare((PetscObject)dmcell,DMDA,&isdmda);
4702c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!isdmda,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only DMDA support for CollectBoundingBox");
471b16650c8SDave May 
472*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm,&dim));
473b16650c8SDave May   sizeof_bbox_ctx = sizeof(CollectBBox);
474b16650c8SDave May   PetscMalloc1(1,&bbox);
475b16650c8SDave May   bbox->owner_rank = rank;
476b16650c8SDave May 
477b16650c8SDave May   /* compute the bounding box based on the overlapping / stenctil size */
478b16650c8SDave May   {
479b16650c8SDave May     Vec lcoor;
480b16650c8SDave May 
481*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetCoordinatesLocal(dmcell,&lcoor));
482fe39f135SDave May     if (dim >= 1) {
483*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecStrideMin(lcoor,0,NULL,&bbox->min[0]));
484*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecStrideMax(lcoor,0,NULL,&bbox->max[0]));
485fe39f135SDave May     }
486fe39f135SDave May     if (dim >= 2) {
487*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecStrideMin(lcoor,1,NULL,&bbox->min[1]));
488*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecStrideMax(lcoor,1,NULL,&bbox->max[1]));
489b16650c8SDave May     }
490fe39f135SDave May     if (dim == 3) {
491*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecStrideMin(lcoor,2,NULL,&bbox->min[2]));
492*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecStrideMax(lcoor,2,NULL,&bbox->max[2]));
493fe39f135SDave May     }
494fe39f135SDave May   }
495*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL));
496b16650c8SDave May   *globalsize = npoints;
497*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval));
498*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm),0,&de));
499b16650c8SDave May   /* use DMDA neighbours */
500*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetNeighbors(dmcell,&dmneighborranks));
5018dbd68bcSDave May   if (dim == 1) {
5028dbd68bcSDave May     neighbour_cells = 3;
5038dbd68bcSDave May   } else if (dim == 2) {
5048dbd68bcSDave May     neighbour_cells = 9;
5058dbd68bcSDave May   } else {
5068dbd68bcSDave May     neighbour_cells = 27;
5078dbd68bcSDave May   }
508*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExTopologyInitialize(de));
509b16650c8SDave May   for (p=0; p<neighbour_cells; p++) {
510b16650c8SDave May     if ((dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank)) {
511*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSwarmDataExTopologyAddNeighbour(de,dmneighborranks[p]));
512b16650c8SDave May     }
513b16650c8SDave May   }
514*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExTopologyFinalize(de));
515*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExInitializeSendCount(de));
516b16650c8SDave May   for (p=0; p<neighbour_cells; p++) {
517b16650c8SDave May     if ((dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank)) {
518*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSwarmDataExAddToSendCount(de,dmneighborranks[p],1));
519b16650c8SDave May     }
520b16650c8SDave May   }
521*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExFinalizeSendCount(de));
522b16650c8SDave May   /* send bounding boxes */
523*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExPackInitialize(de,sizeof_bbox_ctx));
524b16650c8SDave May   for (p=0; p<neighbour_cells; p++) {
525b16650c8SDave May     nrank = dmneighborranks[p];
526b16650c8SDave May     if ((nrank >= 0) && (nrank != rank)) {
52777048351SPatrick Sanan       /* insert bbox buffer into DMSwarmDataExchanger */
528*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSwarmDataExPackData(de,nrank,1,bbox));
529b16650c8SDave May     }
530b16650c8SDave May   }
531*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExPackFinalize(de));
532b16650c8SDave May   /* recv bounding boxes */
533*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExBegin(de));
534*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExEnd(de));
535*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExGetRecvData(de,&n_bbox_recv,(void**)&recv_bbox));
536298827fbSBarry Smith   /*  Wrong, should not be using PETSC_COMM_WORLD */
537b16650c8SDave May   for (p=0; p<n_bbox_recv; p++) {
538298827fbSBarry 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,
539298827fbSBarry 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);
540b16650c8SDave May   }
541*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD,stdout));
542b16650c8SDave May   /* of course this is stupid as this "generic" function should have a better way to know what the coordinates are called */
543*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExInitializeSendCount(de));
544b16650c8SDave May   for (pk=0; pk<n_bbox_recv; pk++) {
545b16650c8SDave May     PetscReal *array_x,*array_y;
546b16650c8SDave May 
547*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x));
548*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y));
549b16650c8SDave May     for (p=0; p<npoints; p++) {
550b16650c8SDave May       if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0])) {
551b16650c8SDave May         if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1])) {
552*5f80ce2aSJacob Faibussowitsch           CHKERRQ(DMSwarmDataExAddToSendCount(de,recv_bbox[pk].owner_rank,1));
553b16650c8SDave May         }
554b16650c8SDave May       }
555b16650c8SDave May     }
556*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y));
557*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x));
558b16650c8SDave May   }
559*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExFinalizeSendCount(de));
560*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer));
561*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExPackInitialize(de,sizeof_dmswarm_point));
562b16650c8SDave May   for (pk=0; pk<n_bbox_recv; pk++) {
563b16650c8SDave May     PetscReal *array_x,*array_y;
564b16650c8SDave May 
565*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x));
566*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y));
567b16650c8SDave May     for (p=0; p<npoints; p++) {
568b16650c8SDave May       if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0])) {
569b16650c8SDave May         if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1])) {
570521f74f9SMatthew G. Knepley           /* copy point into buffer */
571*5f80ce2aSJacob Faibussowitsch           CHKERRQ(DMSwarmDataBucketFillPackedArray(swarm->db,p,point_buffer));
57277048351SPatrick Sanan           /* insert point buffer into DMSwarmDataExchanger */
573*5f80ce2aSJacob Faibussowitsch           CHKERRQ(DMSwarmDataExPackData(de,recv_bbox[pk].owner_rank,1,point_buffer));
574b16650c8SDave May         }
575b16650c8SDave May       }
576b16650c8SDave May     }
577*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y));
578*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x));
579b16650c8SDave May   }
580*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExPackFinalize(de));
581*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval));
582*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExBegin(de));
583*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExEnd(de));
584*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExGetRecvData(de,&n_points_recv,(void**)&recv_points));
585*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL));
586*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketSetSizes(swarm->db,npoints + n_points_recv,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
587b16650c8SDave May   for (p=0; p<n_points_recv; p++) {
588b16650c8SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point);
589b16650c8SDave May 
590*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataBucketInsertPackedArray(swarm->db,npoints+p,data_p));
591b16650c8SDave May   }
592*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketDestroyPackedArray(swarm->db,&point_buffer));
593b16650c8SDave May   PetscFree(bbox);
594*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExView(de));
595*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExDestroy(de));
596b16650c8SDave May   PetscFunctionReturn(0);
597b16650c8SDave May }
598a9fd7477SDave May 
599a9fd7477SDave May /* General collection when no order, or neighbour information is provided */
600a9fd7477SDave May /*
601a9fd7477SDave May  User provides context and collect() method
602a9fd7477SDave May  Broadcast user context
603a9fd7477SDave May 
604a9fd7477SDave May  for each context / rank {
605a9fd7477SDave May    collect(swarm,context,n,list)
606a9fd7477SDave May  }
607a9fd7477SDave May */
608a9fd7477SDave May PETSC_EXTERN PetscErrorCode DMSwarmCollect_General(DM dm,PetscErrorCode (*collect)(DM,void*,PetscInt*,PetscInt**),size_t ctx_size,void *ctx,PetscInt *globalsize)
609a9fd7477SDave May {
610a9fd7477SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
61177048351SPatrick Sanan   DMSwarmDataEx  de;
612a9fd7477SDave May   PetscInt       p,r,npoints,n_points_recv;
613e4fbd051SBarry Smith   PetscMPIInt    size,rank;
614a9fd7477SDave May   void           *point_buffer,*recv_points;
615a9fd7477SDave May   void           *ctxlist;
616a9fd7477SDave May   PetscInt       *n2collect,**collectlist;
617a9fd7477SDave May   size_t         sizeof_dmswarm_point;
618a9fd7477SDave May 
619521f74f9SMatthew G. Knepley   PetscFunctionBegin;
620*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size));
621*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank));
622*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL));
623a9fd7477SDave May   *globalsize = npoints;
624a9fd7477SDave May   /* Broadcast user context */
625e4fbd051SBarry Smith   PetscMalloc(ctx_size*size,&ctxlist);
626*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Allgather(ctx,ctx_size,MPI_CHAR,ctxlist,ctx_size,MPI_CHAR,PetscObjectComm((PetscObject)dm)));
627*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(size,&n2collect));
628*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(size,&collectlist));
629e4fbd051SBarry Smith   for (r=0; r<size; r++) {
630a9fd7477SDave May     PetscInt _n2collect;
631a9fd7477SDave May     PetscInt *_collectlist;
632a9fd7477SDave May     void     *_ctx_r;
633a9fd7477SDave May 
634a9fd7477SDave May     _n2collect   = 0;
635a9fd7477SDave May     _collectlist = NULL;
636a9fd7477SDave May     if (r != rank) { /* don't collect data from yourself */
637a9fd7477SDave May       _ctx_r = (void*)( (char*)ctxlist + r * ctx_size);
638*5f80ce2aSJacob Faibussowitsch       CHKERRQ(collect(dm,_ctx_r,&_n2collect,&_collectlist));
639a9fd7477SDave May     }
640a9fd7477SDave May     n2collect[r]   = _n2collect;
641a9fd7477SDave May     collectlist[r] = _collectlist;
642a9fd7477SDave May   }
643*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm),0,&de));
644a9fd7477SDave May   /* Define topology */
645*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExTopologyInitialize(de));
646e4fbd051SBarry Smith   for (r=0; r<size; r++) {
647a9fd7477SDave May     if (n2collect[r] > 0) {
648*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSwarmDataExTopologyAddNeighbour(de,(PetscMPIInt)r));
649a9fd7477SDave May     }
650a9fd7477SDave May   }
651*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExTopologyFinalize(de));
652a9fd7477SDave May   /* Define send counts */
653*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExInitializeSendCount(de));
654e4fbd051SBarry Smith   for (r=0; r<size; r++) {
655a9fd7477SDave May     if (n2collect[r] > 0) {
656*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSwarmDataExAddToSendCount(de,r,n2collect[r]));
657a9fd7477SDave May     }
658a9fd7477SDave May   }
659*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExFinalizeSendCount(de));
660a9fd7477SDave May   /* Pack data */
661*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer));
662*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExPackInitialize(de,sizeof_dmswarm_point));
663e4fbd051SBarry Smith   for (r=0; r<size; r++) {
664a9fd7477SDave May     for (p=0; p<n2collect[r]; p++) {
665*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSwarmDataBucketFillPackedArray(swarm->db,collectlist[r][p],point_buffer));
666a9fd7477SDave May       /* insert point buffer into the data exchanger */
667*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSwarmDataExPackData(de,r,1,point_buffer));
668a9fd7477SDave May     }
669a9fd7477SDave May   }
670*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExPackFinalize(de));
671a9fd7477SDave May   /* Scatter */
672*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExBegin(de));
673*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExEnd(de));
674a9fd7477SDave May   /* Collect data in DMSwarm container */
675*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExGetRecvData(de,&n_points_recv,(void**)&recv_points));
676*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL));
677*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketSetSizes(swarm->db,npoints + n_points_recv,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
678a9fd7477SDave May   for (p=0; p<n_points_recv; p++) {
679a9fd7477SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point);
680a9fd7477SDave May 
681*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataBucketInsertPackedArray(swarm->db,npoints+p,data_p));
682a9fd7477SDave May   }
683a9fd7477SDave May   /* Release memory */
684e4fbd051SBarry Smith   for (r=0; r<size; r++) {
685a9fd7477SDave May     if (collectlist[r]) PetscFree(collectlist[r]);
686a9fd7477SDave May   }
687*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(collectlist));
688*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(n2collect));
689*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(ctxlist));
690*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketDestroyPackedArray(swarm->db,&point_buffer));
691*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExView(de));
692*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataExDestroy(de));
693a9fd7477SDave May   PetscFunctionReturn(0);
694a9fd7477SDave May }
695