xref: /petsc/src/dm/impls/swarm/swarm_migrate.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 */
119371c9d4SSatish Balay PetscErrorCode DMSwarmMigrate_Push_Basic(DM dm, PetscBool remove_sent_points) {
12df21e3a8SDave May   DM_Swarm     *swarm = (DM_Swarm *)dm->data;
1377048351SPatrick Sanan   DMSwarmDataEx de;
14df21e3a8SDave May   PetscInt      p, npoints, *rankval, n_points_recv;
15df21e3a8SDave May   PetscMPIInt   rank, nrank;
16df21e3a8SDave May   void         *point_buffer, *recv_points;
17df21e3a8SDave May   size_t        sizeof_dmswarm_point;
18356ed814SMatthew G. Knepley   PetscBool     debug = PETSC_FALSE;
19df21e3a8SDave May 
20521f74f9SMatthew G. Knepley   PetscFunctionBegin;
219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
22df21e3a8SDave May 
239566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
249566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
259566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm), 0, &de));
269566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyInitialize(de));
27521f74f9SMatthew G. Knepley   for (p = 0; p < npoints; ++p) {
28df21e3a8SDave May     nrank = rankval[p];
29*48a46eb9SPierre Jolivet     if (nrank != rank) PetscCall(DMSwarmDataExTopologyAddNeighbour(de, nrank));
30df21e3a8SDave May   }
319566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyFinalize(de));
329566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExInitializeSendCount(de));
33df21e3a8SDave May   for (p = 0; p < npoints; p++) {
34df21e3a8SDave May     nrank = rankval[p];
35*48a46eb9SPierre Jolivet     if (nrank != rank) PetscCall(DMSwarmDataExAddToSendCount(de, nrank, 1));
36df21e3a8SDave May   }
379566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExFinalizeSendCount(de));
389566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCreatePackedArray(swarm->db, &sizeof_dmswarm_point, &point_buffer));
399566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackInitialize(de, sizeof_dmswarm_point));
40df21e3a8SDave May   for (p = 0; p < npoints; p++) {
41df21e3a8SDave May     nrank = rankval[p];
42df21e3a8SDave May     if (nrank != rank) {
43df21e3a8SDave May       /* copy point into buffer */
449566063dSJacob Faibussowitsch       PetscCall(DMSwarmDataBucketFillPackedArray(swarm->db, p, point_buffer));
4577048351SPatrick Sanan       /* insert point buffer into DMSwarmDataExchanger */
469566063dSJacob Faibussowitsch       PetscCall(DMSwarmDataExPackData(de, nrank, 1, point_buffer));
47df21e3a8SDave May     }
48df21e3a8SDave May   }
499566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackFinalize(de));
509566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
51df21e3a8SDave May 
52df21e3a8SDave May   if (remove_sent_points) {
5377048351SPatrick Sanan     DMSwarmDataField gfield;
5422a417f9SDave May 
559566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, DMSwarmField_rank, &gfield));
569566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataFieldGetAccess(gfield));
579566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataFieldGetEntries(gfield, (void **)&rankval));
5822a417f9SDave May 
59df21e3a8SDave May     /* remove points which left processor */
609566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
61df21e3a8SDave May     for (p = 0; p < npoints; p++) {
62df21e3a8SDave May       nrank = rankval[p];
63df21e3a8SDave May       if (nrank != rank) {
64df21e3a8SDave May         /* kill point */
659566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataFieldRestoreAccess(gfield));
6622a417f9SDave May 
679566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db, p));
6822a417f9SDave May 
699566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL)); /* you need to update npoints as the list size decreases! */
709566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataFieldGetAccess(gfield));
719566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataFieldGetEntries(gfield, (void **)&rankval));
72df21e3a8SDave May         p--; /* check replacement point */
73df21e3a8SDave May       }
74df21e3a8SDave May     }
759566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataFieldRestoreEntries(gfield, (void **)&rankval));
769566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataFieldRestoreAccess(gfield));
77df21e3a8SDave May   }
789566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExBegin(de));
799566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExEnd(de));
809566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExGetRecvData(de, &n_points_recv, (void **)&recv_points));
819566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
829566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketSetSizes(swarm->db, npoints + n_points_recv, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
83df21e3a8SDave May   for (p = 0; p < n_points_recv; p++) {
84df21e3a8SDave May     void *data_p = (void *)((char *)recv_points + p * sizeof_dmswarm_point);
85df21e3a8SDave May 
869566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketInsertPackedArray(swarm->db, npoints + p, data_p));
87df21e3a8SDave May   }
88356ed814SMatthew G. Knepley   if (debug) PetscCall(DMSwarmDataExView(de));
899566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketDestroyPackedArray(swarm->db, &point_buffer));
909566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExDestroy(de));
91df21e3a8SDave May   PetscFunctionReturn(0);
92df21e3a8SDave May }
932712d1f2SDave May 
949371c9d4SSatish Balay PetscErrorCode DMSwarmMigrate_DMNeighborScatter(DM dm, DM dmcell, PetscBool remove_sent_points, PetscInt *npoints_prior_migration) {
9540c453e9SDave May   DM_Swarm          *swarm = (DM_Swarm *)dm->data;
9677048351SPatrick Sanan   DMSwarmDataEx      de;
9740c453e9SDave May   PetscInt           r, p, npoints, *rankval, n_points_recv;
9840c453e9SDave May   PetscMPIInt        rank, _rank;
9940c453e9SDave May   const PetscMPIInt *neighbourranks;
10040c453e9SDave May   void              *point_buffer, *recv_points;
10140c453e9SDave May   size_t             sizeof_dmswarm_point;
10240c453e9SDave May   PetscInt           nneighbors;
1037c6d1d28SDave May   PetscMPIInt        mynneigh, *myneigh;
10440c453e9SDave May 
105521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1069566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
1079566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
1089566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
1099566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm), 0, &de));
1109566063dSJacob Faibussowitsch   PetscCall(DMGetNeighbors(dmcell, &nneighbors, &neighbourranks));
1119566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyInitialize(de));
11240c453e9SDave May   for (r = 0; r < nneighbors; r++) {
11340c453e9SDave May     _rank = neighbourranks[r];
114*48a46eb9SPierre Jolivet     if ((_rank != rank) && (_rank > 0)) PetscCall(DMSwarmDataExTopologyAddNeighbour(de, _rank));
11540c453e9SDave May   }
1169566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyFinalize(de));
1179566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyGetNeighbours(de, &mynneigh, &myneigh));
1189566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExInitializeSendCount(de));
11940c453e9SDave May   for (p = 0; p < npoints; p++) {
120f954cb40SDave May     if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
1217c6d1d28SDave May       for (r = 0; r < mynneigh; r++) {
1227c6d1d28SDave May         _rank = myneigh[r];
1239566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataExAddToSendCount(de, _rank, 1));
12440c453e9SDave May       }
12540c453e9SDave May     }
12640c453e9SDave May   }
1279566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExFinalizeSendCount(de));
1289566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCreatePackedArray(swarm->db, &sizeof_dmswarm_point, &point_buffer));
1299566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackInitialize(de, sizeof_dmswarm_point));
13040c453e9SDave May   for (p = 0; p < npoints; p++) {
131f954cb40SDave May     if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
1327c6d1d28SDave May       for (r = 0; r < mynneigh; r++) {
1337c6d1d28SDave May         _rank = myneigh[r];
13440c453e9SDave May         /* copy point into buffer */
1359566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataBucketFillPackedArray(swarm->db, p, point_buffer));
13677048351SPatrick Sanan         /* insert point buffer into DMSwarmDataExchanger */
1379566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataExPackData(de, _rank, 1, point_buffer));
13840c453e9SDave May       }
13940c453e9SDave May     }
14040c453e9SDave May   }
1419566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackFinalize(de));
1429566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
14340c453e9SDave May   if (remove_sent_points) {
14477048351SPatrick Sanan     DMSwarmDataField PField;
1457c6d1d28SDave May 
1469566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, DMSwarmField_rank, &PField));
1479566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataFieldGetEntries(PField, (void **)&rankval));
14840c453e9SDave May     /* remove points which left processor */
1499566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
15040c453e9SDave May     for (p = 0; p < npoints; p++) {
151f954cb40SDave May       if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
15240c453e9SDave May         /* kill point */
1539566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db, p));
1549566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL)); /* you need to update npoints as the list size decreases! */
1559566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataFieldGetEntries(PField, (void **)&rankval));      /* update date point increase realloc performed */
15640c453e9SDave May         p--;                                                                   /* check replacement point */
15740c453e9SDave May       }
15840c453e9SDave May     }
15940c453e9SDave May   }
1609566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, npoints_prior_migration, NULL, NULL));
1619566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExBegin(de));
1629566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExEnd(de));
1639566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExGetRecvData(de, &n_points_recv, (void **)&recv_points));
1649566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
1659566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketSetSizes(swarm->db, npoints + n_points_recv, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
16640c453e9SDave May   for (p = 0; p < n_points_recv; p++) {
16740c453e9SDave May     void *data_p = (void *)((char *)recv_points + p * sizeof_dmswarm_point);
16840c453e9SDave May 
1699566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketInsertPackedArray(swarm->db, npoints + p, data_p));
17040c453e9SDave May   }
1719566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketDestroyPackedArray(swarm->db, &point_buffer));
1729566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExDestroy(de));
17340c453e9SDave May   PetscFunctionReturn(0);
17440c453e9SDave May }
175480eef7bSDave May 
1769371c9d4SSatish Balay PetscErrorCode DMSwarmMigrate_CellDMScatter(DM dm, PetscBool remove_sent_points) {
177480eef7bSDave May   DM_Swarm          *swarm = (DM_Swarm *)dm->data;
1786fbf25f8SDave May   PetscInt           p, npoints, npointsg = 0, npoints2, npoints2g, *rankval, npoints_prior_migration;
179bbe8250bSMatthew G. Knepley   PetscSF            sfcell = NULL;
180dfc14de9SMatthew G. Knepley   const PetscSFNode *LA_sfcell;
181480eef7bSDave May   DM                 dmcell;
182480eef7bSDave May   Vec                pos;
18340c453e9SDave May   PetscBool          error_check = swarm->migrate_error_on_missing_point;
184e4fbd051SBarry Smith   PetscMPIInt        size, rank;
185480eef7bSDave May 
186521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1879566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(dm, &dmcell));
18828b400f6SJacob Faibussowitsch   PetscCheck(dmcell, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only valid if cell DM provided");
189480eef7bSDave May 
1909566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
1919566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
1927c6d1d28SDave May 
19343a82f2bSDave May #if 1
19443a82f2bSDave May   {
19543a82f2bSDave May     PetscInt    *p_cellid;
19643a82f2bSDave May     PetscInt     npoints_curr, range = 0;
19743a82f2bSDave May     PetscSFNode *sf_cells;
19843a82f2bSDave May 
1999566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints_curr, NULL, NULL));
2009566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(npoints_curr, &sf_cells));
20143a82f2bSDave May 
2029566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
2039566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&p_cellid));
20443a82f2bSDave May     for (p = 0; p < npoints_curr; p++) {
20543a82f2bSDave May       sf_cells[p].rank  = 0;
20643a82f2bSDave May       sf_cells[p].index = p_cellid[p];
2079371c9d4SSatish Balay       if (p_cellid[p] > range) { range = p_cellid[p]; }
20843a82f2bSDave May     }
2099566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&p_cellid));
2109566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
21143a82f2bSDave May 
2129566063dSJacob Faibussowitsch     /* PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm),&sfcell)); */
2139566063dSJacob Faibussowitsch     PetscCall(PetscSFCreate(PETSC_COMM_SELF, &sfcell));
2149566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(sfcell, range, npoints_curr, NULL, PETSC_OWN_POINTER, sf_cells, PETSC_OWN_POINTER));
21543a82f2bSDave May   }
21643a82f2bSDave May #endif
21743a82f2bSDave May 
2189566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateLocalVectorFromField(dm, DMSwarmPICField_coor, &pos));
2199566063dSJacob Faibussowitsch   PetscCall(DMLocatePoints(dmcell, pos, DM_POINTLOCATION_NONE, &sfcell));
2209566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyLocalVectorFromField(dm, DMSwarmPICField_coor, &pos));
221480eef7bSDave May 
222*48a46eb9SPierre Jolivet   if (error_check) PetscCall(DMSwarmGetSize(dm, &npointsg));
2239566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
2249566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
2259566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell));
2269371c9d4SSatish Balay   for (p = 0; p < npoints; p++) { rankval[p] = LA_sfcell[p].index; }
2279566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
2289566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfcell));
229480eef7bSDave May 
230e4fbd051SBarry Smith   if (size > 1) {
2319566063dSJacob Faibussowitsch     PetscCall(DMSwarmMigrate_DMNeighborScatter(dm, dmcell, remove_sent_points, &npoints_prior_migration));
2326fbf25f8SDave May   } else {
23377048351SPatrick Sanan     DMSwarmDataField PField;
2340ed23c7fSDave May     PetscInt         npoints_curr;
2350ed23c7fSDave May 
2360ed23c7fSDave May     /* remove points which the domain */
2379566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, DMSwarmField_rank, &PField));
2389566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataFieldGetEntries(PField, (void **)&rankval));
2390ed23c7fSDave May 
2409566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints_curr, NULL, NULL));
2410ed23c7fSDave May     for (p = 0; p < npoints_curr; p++) {
2420ed23c7fSDave May       if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
2430ed23c7fSDave May         /* kill point */
2449566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db, p));
2459566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints_curr, NULL, NULL)); /* you need to update npoints as the list size decreases! */
2469566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataFieldGetEntries(PField, (void **)&rankval));           /* update date point increase realloc performed */
2470ed23c7fSDave May         p--;                                                                        /* check replacement point */
2480ed23c7fSDave May       }
2490ed23c7fSDave May     }
2509566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dm, &npoints_prior_migration));
2516fbf25f8SDave May   }
252480eef7bSDave May 
2532d4ee042Sprj-   /* locate points newly received */
2549566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints2, NULL, NULL));
255009b43efSDave May 
2567c6d1d28SDave May #if 0
2572d4ee042Sprj-   { /* safe alternative - however this performs two point locations on: (i) the initial points set and; (ii) the (initial + received) point set */
2587c6d1d28SDave May     PetscScalar *LA_coor;
2597c6d1d28SDave May     PetscInt bs;
26077048351SPatrick Sanan     DMSwarmDataField PField;
2617c6d1d28SDave May 
2629566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor));
2639566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,bs,bs*npoints2,(const PetscScalar*)LA_coor,&pos));
2649566063dSJacob Faibussowitsch     PetscCall(DMLocatePoints(dmcell,pos,DM_POINTLOCATION_NONE,&sfcell));
2657c6d1d28SDave May 
2669566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&pos));
2679566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor));
2687c6d1d28SDave May 
2699566063dSJacob Faibussowitsch     PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell));
2709566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval));
2717c6d1d28SDave May     for (p=0; p<npoints2; p++) {
272dfc14de9SMatthew G. Knepley       rankval[p] = LA_sfcell[p].index;
2737c6d1d28SDave May     }
2749566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&sfcell));
2759566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval));
27640c453e9SDave May 
2777c6d1d28SDave May     /* remove points which left processor */
2789566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,DMSwarmField_rank,&PField));
2799566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataFieldGetEntries(PField,(void**)&rankval));
2807c6d1d28SDave May 
2819566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints2,NULL,NULL));
2827c6d1d28SDave May     for (p=0; p<npoints2; p++) {
283f954cb40SDave May       if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
2847c6d1d28SDave May         /* kill point */
2859566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db,p));
2869566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints2,NULL,NULL)); /* you need to update npoints as the list size decreases! */
2879566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataFieldGetEntries(PField,(void**)&rankval)); /* update date point increase realloc performed */
2887c6d1d28SDave May         p--; /* check replacement point */
2897c6d1d28SDave May       }
2907c6d1d28SDave May     }
29140c453e9SDave May   }
292009b43efSDave May #endif
293009b43efSDave May 
2945627991aSBarry Smith   { /* perform two point locations: (i) on the initial points set prior to communication; and (ii) on the new (received) points */
295009b43efSDave May     PetscScalar     *LA_coor;
296009b43efSDave May     PetscInt         npoints_from_neighbours, bs;
29777048351SPatrick Sanan     DMSwarmDataField PField;
298009b43efSDave May 
299009b43efSDave May     npoints_from_neighbours = npoints2 - npoints_prior_migration;
300009b43efSDave May 
3019566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, &bs, NULL, (void **)&LA_coor));
3029566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, bs, bs * npoints_from_neighbours, (const PetscScalar *)&LA_coor[bs * npoints_prior_migration], &pos));
303009b43efSDave May 
3049566063dSJacob Faibussowitsch     PetscCall(DMLocatePoints(dmcell, pos, DM_POINTLOCATION_NONE, &sfcell));
305009b43efSDave May 
3069566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&pos));
3079566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, &bs, NULL, (void **)&LA_coor));
308009b43efSDave May 
3099566063dSJacob Faibussowitsch     PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell));
3109566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
3119371c9d4SSatish Balay     for (p = 0; p < npoints_from_neighbours; p++) { rankval[npoints_prior_migration + p] = LA_sfcell[p].index; }
3129566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
3139566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&sfcell));
314009b43efSDave May 
315009b43efSDave May     /* remove points which left processor */
3169566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, DMSwarmField_rank, &PField));
3179566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataFieldGetEntries(PField, (void **)&rankval));
318009b43efSDave May 
3199566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints2, NULL, NULL));
320009b43efSDave May     for (p = npoints_prior_migration; p < npoints2; p++) {
321009b43efSDave May       if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
322009b43efSDave May         /* kill point */
3239566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db, p));
3249566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints2, NULL, NULL)); /* you need to update npoints as the list size decreases! */
3259566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataFieldGetEntries(PField, (void **)&rankval));       /* update date point increase realloc performed */
326009b43efSDave May         p--;                                                                    /* check replacement point */
327009b43efSDave May       }
328009b43efSDave May     }
329009b43efSDave May   }
330009b43efSDave May 
3313151f1d5SDave May   {
3323151f1d5SDave May     PetscInt *p_cellid;
3333151f1d5SDave May 
3349566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints2, NULL, NULL));
3359566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
3369566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&p_cellid));
3379371c9d4SSatish Balay     for (p = 0; p < npoints2; p++) { p_cellid[p] = rankval[p]; }
3389566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&p_cellid));
3399566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
3403151f1d5SDave May   }
3413151f1d5SDave May 
34240c453e9SDave May   /* check for error on removed points */
34340c453e9SDave May   if (error_check) {
3449566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dm, &npoints2g));
34563a3b9bcSJacob Faibussowitsch     PetscCheck(npointsg == npoints2g, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Points from the DMSwarm must remain constant during migration (initial %" PetscInt_FMT " - final %" PetscInt_FMT ")", npointsg, npoints2g);
34640c453e9SDave May   }
347480eef7bSDave May   PetscFunctionReturn(0);
348480eef7bSDave May }
349480eef7bSDave May 
3509371c9d4SSatish Balay PetscErrorCode DMSwarmMigrate_CellDMExact(DM dm, PetscBool remove_sent_points) {
351521f74f9SMatthew G. Knepley   PetscFunctionBegin;
35208056efcSDave May   PetscFunctionReturn(0);
35308056efcSDave May }
35408056efcSDave May 
355480eef7bSDave May /*
356480eef7bSDave May  Redundant as this assumes points can only be sent to a single rank
357480eef7bSDave May */
3589371c9d4SSatish Balay PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm, PetscInt *globalsize) {
3592712d1f2SDave May   DM_Swarm     *swarm = (DM_Swarm *)dm->data;
36077048351SPatrick Sanan   DMSwarmDataEx de;
3612712d1f2SDave May   PetscInt      p, npoints, *rankval, n_points_recv;
3622712d1f2SDave May   PetscMPIInt   rank, nrank, negrank;
3632712d1f2SDave May   void         *point_buffer, *recv_points;
3642712d1f2SDave May   size_t        sizeof_dmswarm_point;
3652712d1f2SDave May 
366521f74f9SMatthew G. Knepley   PetscFunctionBegin;
3679566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
3689566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
3692712d1f2SDave May   *globalsize = npoints;
3709566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
3719566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm), 0, &de));
3729566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyInitialize(de));
3732712d1f2SDave May   for (p = 0; p < npoints; p++) {
3742712d1f2SDave May     negrank = rankval[p];
3752712d1f2SDave May     if (negrank < 0) {
3762712d1f2SDave May       nrank = -negrank - 1;
3779566063dSJacob Faibussowitsch       PetscCall(DMSwarmDataExTopologyAddNeighbour(de, nrank));
3782712d1f2SDave May     }
3792712d1f2SDave May   }
3809566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyFinalize(de));
3819566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExInitializeSendCount(de));
3822712d1f2SDave May   for (p = 0; p < npoints; p++) {
3832712d1f2SDave May     negrank = rankval[p];
3842712d1f2SDave May     if (negrank < 0) {
3852712d1f2SDave May       nrank = -negrank - 1;
3869566063dSJacob Faibussowitsch       PetscCall(DMSwarmDataExAddToSendCount(de, nrank, 1));
3872712d1f2SDave May     }
3882712d1f2SDave May   }
3899566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExFinalizeSendCount(de));
3909566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCreatePackedArray(swarm->db, &sizeof_dmswarm_point, &point_buffer));
3919566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackInitialize(de, sizeof_dmswarm_point));
3922712d1f2SDave May   for (p = 0; p < npoints; p++) {
3932712d1f2SDave May     negrank = rankval[p];
3942712d1f2SDave May     if (negrank < 0) {
3952712d1f2SDave May       nrank      = -negrank - 1;
3962712d1f2SDave May       rankval[p] = nrank;
3972712d1f2SDave May       /* copy point into buffer */
3989566063dSJacob Faibussowitsch       PetscCall(DMSwarmDataBucketFillPackedArray(swarm->db, p, point_buffer));
39977048351SPatrick Sanan       /* insert point buffer into DMSwarmDataExchanger */
4009566063dSJacob Faibussowitsch       PetscCall(DMSwarmDataExPackData(de, nrank, 1, point_buffer));
4012712d1f2SDave May       rankval[p] = negrank;
4022712d1f2SDave May     }
4032712d1f2SDave May   }
4049566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackFinalize(de));
4059566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
4069566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExBegin(de));
4079566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExEnd(de));
4089566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExGetRecvData(de, &n_points_recv, (void **)&recv_points));
4099566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
4109566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketSetSizes(swarm->db, npoints + n_points_recv, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
4112712d1f2SDave May   for (p = 0; p < n_points_recv; p++) {
4122712d1f2SDave May     void *data_p = (void *)((char *)recv_points + p * sizeof_dmswarm_point);
4132712d1f2SDave May 
4149566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketInsertPackedArray(swarm->db, npoints + p, data_p));
4152712d1f2SDave May   }
4169566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExView(de));
4179566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketDestroyPackedArray(swarm->db, &point_buffer));
4189566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExDestroy(de));
4192712d1f2SDave May   PetscFunctionReturn(0);
4202712d1f2SDave May }
421b16650c8SDave May 
422b16650c8SDave May typedef struct {
423b16650c8SDave May   PetscMPIInt owner_rank;
424b16650c8SDave May   PetscReal   min[3], max[3];
425b16650c8SDave May } CollectBBox;
426b16650c8SDave May 
4279371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMSwarmCollect_DMDABoundingBox(DM dm, PetscInt *globalsize) {
428b16650c8SDave May   DM_Swarm          *swarm = (DM_Swarm *)dm->data;
42977048351SPatrick Sanan   DMSwarmDataEx      de;
430b16650c8SDave May   PetscInt           p, pk, npoints, *rankval, n_points_recv, n_bbox_recv, dim, neighbour_cells;
431b16650c8SDave May   PetscMPIInt        rank, nrank;
432b16650c8SDave May   void              *point_buffer, *recv_points;
433b16650c8SDave May   size_t             sizeof_dmswarm_point, sizeof_bbox_ctx;
434b16650c8SDave May   PetscBool          isdmda;
435b16650c8SDave May   CollectBBox       *bbox, *recv_bbox;
436b16650c8SDave May   const PetscMPIInt *dmneighborranks;
437b16650c8SDave May   DM                 dmcell;
438b16650c8SDave May 
439521f74f9SMatthew G. Knepley   PetscFunctionBegin;
4409566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
441b16650c8SDave May 
4429566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(dm, &dmcell));
44328b400f6SJacob Faibussowitsch   PetscCheck(dmcell, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only valid if cell DM provided");
444b16650c8SDave May   isdmda = PETSC_FALSE;
445b16650c8SDave May   PetscObjectTypeCompare((PetscObject)dmcell, DMDA, &isdmda);
44628b400f6SJacob Faibussowitsch   PetscCheck(isdmda, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only DMDA support for CollectBoundingBox");
447b16650c8SDave May 
4489566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
449b16650c8SDave May   sizeof_bbox_ctx = sizeof(CollectBBox);
450b16650c8SDave May   PetscMalloc1(1, &bbox);
451b16650c8SDave May   bbox->owner_rank = rank;
452b16650c8SDave May 
453b16650c8SDave May   /* compute the bounding box based on the overlapping / stenctil size */
454b16650c8SDave May   {
455b16650c8SDave May     Vec lcoor;
456b16650c8SDave May 
4579566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(dmcell, &lcoor));
458fe39f135SDave May     if (dim >= 1) {
4599566063dSJacob Faibussowitsch       PetscCall(VecStrideMin(lcoor, 0, NULL, &bbox->min[0]));
4609566063dSJacob Faibussowitsch       PetscCall(VecStrideMax(lcoor, 0, NULL, &bbox->max[0]));
461fe39f135SDave May     }
462fe39f135SDave May     if (dim >= 2) {
4639566063dSJacob Faibussowitsch       PetscCall(VecStrideMin(lcoor, 1, NULL, &bbox->min[1]));
4649566063dSJacob Faibussowitsch       PetscCall(VecStrideMax(lcoor, 1, NULL, &bbox->max[1]));
465b16650c8SDave May     }
466fe39f135SDave May     if (dim == 3) {
4679566063dSJacob Faibussowitsch       PetscCall(VecStrideMin(lcoor, 2, NULL, &bbox->min[2]));
4689566063dSJacob Faibussowitsch       PetscCall(VecStrideMax(lcoor, 2, NULL, &bbox->max[2]));
469fe39f135SDave May     }
470fe39f135SDave May   }
4719566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
472b16650c8SDave May   *globalsize = npoints;
4739566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
4749566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm), 0, &de));
475b16650c8SDave May   /* use DMDA neighbours */
4769566063dSJacob Faibussowitsch   PetscCall(DMDAGetNeighbors(dmcell, &dmneighborranks));
4778dbd68bcSDave May   if (dim == 1) {
4788dbd68bcSDave May     neighbour_cells = 3;
4798dbd68bcSDave May   } else if (dim == 2) {
4808dbd68bcSDave May     neighbour_cells = 9;
4818dbd68bcSDave May   } else {
4828dbd68bcSDave May     neighbour_cells = 27;
4838dbd68bcSDave May   }
4849566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyInitialize(de));
485b16650c8SDave May   for (p = 0; p < neighbour_cells; p++) {
486*48a46eb9SPierre Jolivet     if ((dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank)) PetscCall(DMSwarmDataExTopologyAddNeighbour(de, dmneighborranks[p]));
487b16650c8SDave May   }
4889566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyFinalize(de));
4899566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExInitializeSendCount(de));
490b16650c8SDave May   for (p = 0; p < neighbour_cells; p++) {
491*48a46eb9SPierre Jolivet     if ((dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank)) PetscCall(DMSwarmDataExAddToSendCount(de, dmneighborranks[p], 1));
492b16650c8SDave May   }
4939566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExFinalizeSendCount(de));
494b16650c8SDave May   /* send bounding boxes */
4959566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackInitialize(de, sizeof_bbox_ctx));
496b16650c8SDave May   for (p = 0; p < neighbour_cells; p++) {
497b16650c8SDave May     nrank = dmneighborranks[p];
498b16650c8SDave May     if ((nrank >= 0) && (nrank != rank)) {
49977048351SPatrick Sanan       /* insert bbox buffer into DMSwarmDataExchanger */
5009566063dSJacob Faibussowitsch       PetscCall(DMSwarmDataExPackData(de, nrank, 1, bbox));
501b16650c8SDave May     }
502b16650c8SDave May   }
5039566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackFinalize(de));
504b16650c8SDave May   /* recv bounding boxes */
5059566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExBegin(de));
5069566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExEnd(de));
5079566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExGetRecvData(de, &n_bbox_recv, (void **)&recv_bbox));
508298827fbSBarry Smith   /*  Wrong, should not be using PETSC_COMM_WORLD */
509b16650c8SDave May   for (p = 0; p < n_bbox_recv; p++) {
5109371c9d4SSatish Balay     PetscCall(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, (double)recv_bbox[p].min[0], (double)recv_bbox[p].max[0], (double)recv_bbox[p].min[1],
5119371c9d4SSatish Balay                                       (double)recv_bbox[p].max[1]));
512b16650c8SDave May   }
5139566063dSJacob Faibussowitsch   PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, stdout));
514b16650c8SDave May   /* of course this is stupid as this "generic" function should have a better way to know what the coordinates are called */
5159566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExInitializeSendCount(de));
516b16650c8SDave May   for (pk = 0; pk < n_bbox_recv; pk++) {
517b16650c8SDave May     PetscReal *array_x, *array_y;
518b16650c8SDave May 
5199566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dm, "coorx", NULL, NULL, (void **)&array_x));
5209566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dm, "coory", NULL, NULL, (void **)&array_y));
521b16650c8SDave May     for (p = 0; p < npoints; p++) {
522b16650c8SDave May       if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0])) {
523*48a46eb9SPierre Jolivet         if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1])) PetscCall(DMSwarmDataExAddToSendCount(de, recv_bbox[pk].owner_rank, 1));
524b16650c8SDave May       }
525b16650c8SDave May     }
5269566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dm, "coory", NULL, NULL, (void **)&array_y));
5279566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dm, "coorx", NULL, NULL, (void **)&array_x));
528b16650c8SDave May   }
5299566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExFinalizeSendCount(de));
5309566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCreatePackedArray(swarm->db, &sizeof_dmswarm_point, &point_buffer));
5319566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackInitialize(de, sizeof_dmswarm_point));
532b16650c8SDave May   for (pk = 0; pk < n_bbox_recv; pk++) {
533b16650c8SDave May     PetscReal *array_x, *array_y;
534b16650c8SDave May 
5359566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dm, "coorx", NULL, NULL, (void **)&array_x));
5369566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dm, "coory", NULL, NULL, (void **)&array_y));
537b16650c8SDave May     for (p = 0; p < npoints; p++) {
538b16650c8SDave May       if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0])) {
539b16650c8SDave May         if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1])) {
540521f74f9SMatthew G. Knepley           /* copy point into buffer */
5419566063dSJacob Faibussowitsch           PetscCall(DMSwarmDataBucketFillPackedArray(swarm->db, p, point_buffer));
54277048351SPatrick Sanan           /* insert point buffer into DMSwarmDataExchanger */
5439566063dSJacob Faibussowitsch           PetscCall(DMSwarmDataExPackData(de, recv_bbox[pk].owner_rank, 1, point_buffer));
544b16650c8SDave May         }
545b16650c8SDave May       }
546b16650c8SDave May     }
5479566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dm, "coory", NULL, NULL, (void **)&array_y));
5489566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dm, "coorx", NULL, NULL, (void **)&array_x));
549b16650c8SDave May   }
5509566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackFinalize(de));
5519566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
5529566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExBegin(de));
5539566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExEnd(de));
5549566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExGetRecvData(de, &n_points_recv, (void **)&recv_points));
5559566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
5569566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketSetSizes(swarm->db, npoints + n_points_recv, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
557b16650c8SDave May   for (p = 0; p < n_points_recv; p++) {
558b16650c8SDave May     void *data_p = (void *)((char *)recv_points + p * sizeof_dmswarm_point);
559b16650c8SDave May 
5609566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketInsertPackedArray(swarm->db, npoints + p, data_p));
561b16650c8SDave May   }
5629566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketDestroyPackedArray(swarm->db, &point_buffer));
563b16650c8SDave May   PetscFree(bbox);
5649566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExView(de));
5659566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExDestroy(de));
566b16650c8SDave May   PetscFunctionReturn(0);
567b16650c8SDave May }
568a9fd7477SDave May 
569a9fd7477SDave May /* General collection when no order, or neighbour information is provided */
570a9fd7477SDave May /*
571a9fd7477SDave May  User provides context and collect() method
572a9fd7477SDave May  Broadcast user context
573a9fd7477SDave May 
574a9fd7477SDave May  for each context / rank {
575a9fd7477SDave May    collect(swarm,context,n,list)
576a9fd7477SDave May  }
577a9fd7477SDave May */
5789371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMSwarmCollect_General(DM dm, PetscErrorCode (*collect)(DM, void *, PetscInt *, PetscInt **), size_t ctx_size, void *ctx, PetscInt *globalsize) {
579a9fd7477SDave May   DM_Swarm     *swarm = (DM_Swarm *)dm->data;
58077048351SPatrick Sanan   DMSwarmDataEx de;
581a9fd7477SDave May   PetscInt      p, r, npoints, n_points_recv;
582e4fbd051SBarry Smith   PetscMPIInt   size, rank;
583a9fd7477SDave May   void         *point_buffer, *recv_points;
584a9fd7477SDave May   void         *ctxlist;
585a9fd7477SDave May   PetscInt     *n2collect, **collectlist;
586a9fd7477SDave May   size_t        sizeof_dmswarm_point;
587a9fd7477SDave May 
588521f74f9SMatthew G. Knepley   PetscFunctionBegin;
5899566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
5909566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
5919566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
592a9fd7477SDave May   *globalsize = npoints;
593a9fd7477SDave May   /* Broadcast user context */
594e4fbd051SBarry Smith   PetscMalloc(ctx_size * size, &ctxlist);
5959566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(ctx, ctx_size, MPI_CHAR, ctxlist, ctx_size, MPI_CHAR, PetscObjectComm((PetscObject)dm)));
5969566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size, &n2collect));
5979566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size, &collectlist));
598e4fbd051SBarry Smith   for (r = 0; r < size; r++) {
599a9fd7477SDave May     PetscInt  _n2collect;
600a9fd7477SDave May     PetscInt *_collectlist;
601a9fd7477SDave May     void     *_ctx_r;
602a9fd7477SDave May 
603a9fd7477SDave May     _n2collect   = 0;
604a9fd7477SDave May     _collectlist = NULL;
605a9fd7477SDave May     if (r != rank) { /* don't collect data from yourself */
606a9fd7477SDave May       _ctx_r = (void *)((char *)ctxlist + r * ctx_size);
6079566063dSJacob Faibussowitsch       PetscCall(collect(dm, _ctx_r, &_n2collect, &_collectlist));
608a9fd7477SDave May     }
609a9fd7477SDave May     n2collect[r]   = _n2collect;
610a9fd7477SDave May     collectlist[r] = _collectlist;
611a9fd7477SDave May   }
6129566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm), 0, &de));
613a9fd7477SDave May   /* Define topology */
6149566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyInitialize(de));
615e4fbd051SBarry Smith   for (r = 0; r < size; r++) {
616*48a46eb9SPierre Jolivet     if (n2collect[r] > 0) PetscCall(DMSwarmDataExTopologyAddNeighbour(de, (PetscMPIInt)r));
617a9fd7477SDave May   }
6189566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyFinalize(de));
619a9fd7477SDave May   /* Define send counts */
6209566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExInitializeSendCount(de));
621e4fbd051SBarry Smith   for (r = 0; r < size; r++) {
622*48a46eb9SPierre Jolivet     if (n2collect[r] > 0) PetscCall(DMSwarmDataExAddToSendCount(de, r, n2collect[r]));
623a9fd7477SDave May   }
6249566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExFinalizeSendCount(de));
625a9fd7477SDave May   /* Pack data */
6269566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCreatePackedArray(swarm->db, &sizeof_dmswarm_point, &point_buffer));
6279566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackInitialize(de, sizeof_dmswarm_point));
628e4fbd051SBarry Smith   for (r = 0; r < size; r++) {
629a9fd7477SDave May     for (p = 0; p < n2collect[r]; p++) {
6309566063dSJacob Faibussowitsch       PetscCall(DMSwarmDataBucketFillPackedArray(swarm->db, collectlist[r][p], point_buffer));
631a9fd7477SDave May       /* insert point buffer into the data exchanger */
6329566063dSJacob Faibussowitsch       PetscCall(DMSwarmDataExPackData(de, r, 1, point_buffer));
633a9fd7477SDave May     }
634a9fd7477SDave May   }
6359566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackFinalize(de));
636a9fd7477SDave May   /* Scatter */
6379566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExBegin(de));
6389566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExEnd(de));
639a9fd7477SDave May   /* Collect data in DMSwarm container */
6409566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExGetRecvData(de, &n_points_recv, (void **)&recv_points));
6419566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
6429566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketSetSizes(swarm->db, npoints + n_points_recv, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
643a9fd7477SDave May   for (p = 0; p < n_points_recv; p++) {
644a9fd7477SDave May     void *data_p = (void *)((char *)recv_points + p * sizeof_dmswarm_point);
645a9fd7477SDave May 
6469566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketInsertPackedArray(swarm->db, npoints + p, data_p));
647a9fd7477SDave May   }
648a9fd7477SDave May   /* Release memory */
649e4fbd051SBarry Smith   for (r = 0; r < size; r++) {
650a9fd7477SDave May     if (collectlist[r]) PetscFree(collectlist[r]);
651a9fd7477SDave May   }
6529566063dSJacob Faibussowitsch   PetscCall(PetscFree(collectlist));
6539566063dSJacob Faibussowitsch   PetscCall(PetscFree(n2collect));
6549566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctxlist));
6559566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketDestroyPackedArray(swarm->db, &point_buffer));
6569566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExView(de));
6579566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExDestroy(de));
658a9fd7477SDave May   PetscFunctionReturn(0);
659a9fd7477SDave May }
660