xref: /petsc/src/dm/impls/swarm/swarm_migrate.c (revision 953fc092c940ca8d0a1429c40027a832a7470e56)
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
106497c311SBarry Smith 
116497c311SBarry Smith  It should be storing the rank information as MPIInt not Int
12480eef7bSDave May */
13d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmMigrate_Push_Basic(DM dm, PetscBool remove_sent_points)
14d71ae5a4SJacob Faibussowitsch {
15df21e3a8SDave May   DM_Swarm     *swarm = (DM_Swarm *)dm->data;
1677048351SPatrick Sanan   DMSwarmDataEx de;
17df21e3a8SDave May   PetscInt      p, npoints, *rankval, n_points_recv;
18df21e3a8SDave May   PetscMPIInt   rank, nrank;
19df21e3a8SDave May   void         *point_buffer, *recv_points;
20df21e3a8SDave May   size_t        sizeof_dmswarm_point;
21356ed814SMatthew G. Knepley   PetscBool     debug = PETSC_FALSE;
22df21e3a8SDave May 
23521f74f9SMatthew G. Knepley   PetscFunctionBegin;
249566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
25df21e3a8SDave May 
269566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
279566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
289566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm), 0, &de));
299566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyInitialize(de));
30521f74f9SMatthew G. Knepley   for (p = 0; p < npoints; ++p) {
31835f2295SStefano Zampini     PetscCall(PetscMPIIntCast(rankval[p], &nrank));
3248a46eb9SPierre Jolivet     if (nrank != rank) PetscCall(DMSwarmDataExTopologyAddNeighbour(de, nrank));
33df21e3a8SDave May   }
349566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyFinalize(de));
359566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExInitializeSendCount(de));
36df21e3a8SDave May   for (p = 0; p < npoints; p++) {
37835f2295SStefano Zampini     PetscCall(PetscMPIIntCast(rankval[p], &nrank));
3848a46eb9SPierre Jolivet     if (nrank != rank) PetscCall(DMSwarmDataExAddToSendCount(de, nrank, 1));
39df21e3a8SDave May   }
409566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExFinalizeSendCount(de));
419566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCreatePackedArray(swarm->db, &sizeof_dmswarm_point, &point_buffer));
429566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackInitialize(de, sizeof_dmswarm_point));
43df21e3a8SDave May   for (p = 0; p < npoints; p++) {
44835f2295SStefano Zampini     PetscCall(PetscMPIIntCast(rankval[p], &nrank));
45df21e3a8SDave May     if (nrank != rank) {
46df21e3a8SDave May       /* copy point into buffer */
479566063dSJacob Faibussowitsch       PetscCall(DMSwarmDataBucketFillPackedArray(swarm->db, p, point_buffer));
4877048351SPatrick Sanan       /* insert point buffer into DMSwarmDataExchanger */
499566063dSJacob Faibussowitsch       PetscCall(DMSwarmDataExPackData(de, nrank, 1, point_buffer));
50df21e3a8SDave May     }
51df21e3a8SDave May   }
529566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackFinalize(de));
539566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
54df21e3a8SDave May 
55df21e3a8SDave May   if (remove_sent_points) {
5677048351SPatrick Sanan     DMSwarmDataField gfield;
5722a417f9SDave May 
589566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, DMSwarmField_rank, &gfield));
599566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataFieldGetAccess(gfield));
609566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataFieldGetEntries(gfield, (void **)&rankval));
6122a417f9SDave May 
62df21e3a8SDave May     /* remove points which left processor */
639566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
64df21e3a8SDave May     for (p = 0; p < npoints; p++) {
65835f2295SStefano Zampini       PetscCall(PetscMPIIntCast(rankval[p], &nrank));
66df21e3a8SDave May       if (nrank != rank) {
67df21e3a8SDave May         /* kill point */
689566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataFieldRestoreAccess(gfield));
6922a417f9SDave May 
709566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db, p));
7122a417f9SDave May 
729566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL)); /* you need to update npoints as the list size decreases! */
739566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataFieldGetAccess(gfield));
749566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataFieldGetEntries(gfield, (void **)&rankval));
75df21e3a8SDave May         p--; /* check replacement point */
76df21e3a8SDave May       }
77df21e3a8SDave May     }
789566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataFieldRestoreEntries(gfield, (void **)&rankval));
799566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataFieldRestoreAccess(gfield));
80df21e3a8SDave May   }
819566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExBegin(de));
829566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExEnd(de));
83835f2295SStefano Zampini   PetscCall(DMSwarmDataExGetRecvData(de, &n_points_recv, &recv_points));
849566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
859566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketSetSizes(swarm->db, npoints + n_points_recv, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
86df21e3a8SDave May   for (p = 0; p < n_points_recv; p++) {
87df21e3a8SDave May     void *data_p = (void *)((char *)recv_points + p * sizeof_dmswarm_point);
88df21e3a8SDave May 
899566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketInsertPackedArray(swarm->db, npoints + p, data_p));
90df21e3a8SDave May   }
91356ed814SMatthew G. Knepley   if (debug) PetscCall(DMSwarmDataExView(de));
929566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketDestroyPackedArray(swarm->db, &point_buffer));
939566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExDestroy(de));
943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
95df21e3a8SDave May }
962712d1f2SDave May 
9766976f2fSJacob Faibussowitsch static PetscErrorCode DMSwarmMigrate_DMNeighborScatter(DM dm, DM dmcell, PetscBool remove_sent_points, PetscInt *npoints_prior_migration)
98d71ae5a4SJacob Faibussowitsch {
9940c453e9SDave May   DM_Swarm          *swarm = (DM_Swarm *)dm->data;
10019307e5cSMatthew G. Knepley   DMSwarmCellDM      celldm;
10177048351SPatrick Sanan   DMSwarmDataEx      de;
1021f70fafbSZach Atkins   PetscInt           r, p, npoints, *p_cellid, 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;
10919307e5cSMatthew G. Knepley   const char        *cellid;
11040c453e9SDave May 
111521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1129566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
1139566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
11419307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetCellDMActive(dm, &celldm));
11519307e5cSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetCellID(celldm, &cellid));
11619307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetField(dm, cellid, NULL, NULL, (void **)&p_cellid));
1179566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm), 0, &de));
1189566063dSJacob Faibussowitsch   PetscCall(DMGetNeighbors(dmcell, &nneighbors, &neighbourranks));
1199566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyInitialize(de));
12040c453e9SDave May   for (r = 0; r < nneighbors; r++) {
12140c453e9SDave May     _rank = neighbourranks[r];
12248a46eb9SPierre Jolivet     if ((_rank != rank) && (_rank > 0)) PetscCall(DMSwarmDataExTopologyAddNeighbour(de, _rank));
12340c453e9SDave May   }
1249566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyFinalize(de));
1259566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyGetNeighbours(de, &mynneigh, &myneigh));
1269566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExInitializeSendCount(de));
12740c453e9SDave May   for (p = 0; p < npoints; p++) {
1281f70fafbSZach Atkins     if (p_cellid[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
1297c6d1d28SDave May       for (r = 0; r < mynneigh; r++) {
1307c6d1d28SDave May         _rank = myneigh[r];
1319566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataExAddToSendCount(de, _rank, 1));
13240c453e9SDave May       }
13340c453e9SDave May     }
13440c453e9SDave May   }
1359566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExFinalizeSendCount(de));
1369566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCreatePackedArray(swarm->db, &sizeof_dmswarm_point, &point_buffer));
1379566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackInitialize(de, sizeof_dmswarm_point));
13840c453e9SDave May   for (p = 0; p < npoints; p++) {
1391f70fafbSZach Atkins     if (p_cellid[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
1407c6d1d28SDave May       for (r = 0; r < mynneigh; r++) {
1417c6d1d28SDave May         _rank = myneigh[r];
14240c453e9SDave May         /* copy point into buffer */
1439566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataBucketFillPackedArray(swarm->db, p, point_buffer));
14477048351SPatrick Sanan         /* insert point buffer into DMSwarmDataExchanger */
1459566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataExPackData(de, _rank, 1, point_buffer));
14640c453e9SDave May       }
14740c453e9SDave May     }
14840c453e9SDave May   }
1499566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackFinalize(de));
15019307e5cSMatthew G. Knepley   PetscCall(DMSwarmRestoreField(dm, cellid, NULL, NULL, (void **)&p_cellid));
15140c453e9SDave May   if (remove_sent_points) {
15277048351SPatrick Sanan     DMSwarmDataField PField;
1537c6d1d28SDave May 
15419307e5cSMatthew G. Knepley     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, cellid, &PField));
1551f70fafbSZach Atkins     PetscCall(DMSwarmDataFieldGetEntries(PField, (void **)&p_cellid));
15640c453e9SDave May     /* remove points which left processor */
1579566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
15840c453e9SDave May     for (p = 0; p < npoints; p++) {
1591f70fafbSZach Atkins       if (p_cellid[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
16040c453e9SDave May         /* kill point */
1619566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db, p));
1629566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL)); /* you need to update npoints as the list size decreases! */
1631f70fafbSZach Atkins         PetscCall(DMSwarmDataFieldGetEntries(PField, (void **)&p_cellid));     /* update date point increase realloc performed */
16440c453e9SDave May         p--;                                                                   /* check replacement point */
16540c453e9SDave May       }
16640c453e9SDave May     }
16740c453e9SDave May   }
1689566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, npoints_prior_migration, NULL, NULL));
1699566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExBegin(de));
1709566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExEnd(de));
171835f2295SStefano Zampini   PetscCall(DMSwarmDataExGetRecvData(de, &n_points_recv, &recv_points));
1729566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
1739566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketSetSizes(swarm->db, npoints + n_points_recv, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
17440c453e9SDave May   for (p = 0; p < n_points_recv; p++) {
17540c453e9SDave May     void *data_p = (void *)((char *)recv_points + p * sizeof_dmswarm_point);
17640c453e9SDave May 
1779566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketInsertPackedArray(swarm->db, npoints + p, data_p));
17840c453e9SDave May   }
1799566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketDestroyPackedArray(swarm->db, &point_buffer));
1809566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExDestroy(de));
1813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18240c453e9SDave May }
183480eef7bSDave May 
184d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmMigrate_CellDMScatter(DM dm, PetscBool remove_sent_points)
185d71ae5a4SJacob Faibussowitsch {
186480eef7bSDave May   DM_Swarm          *swarm = (DM_Swarm *)dm->data;
18719307e5cSMatthew G. Knepley   DMSwarmCellDM      celldm;
18819307e5cSMatthew G. Knepley   PetscInt           p, npoints, npointsg = 0, npoints2, npoints2g, *rankval, *p_cellid, npoints_prior_migration, Nfc;
189bbe8250bSMatthew G. Knepley   PetscSF            sfcell = NULL;
190dfc14de9SMatthew G. Knepley   const PetscSFNode *LA_sfcell;
191480eef7bSDave May   DM                 dmcell;
192480eef7bSDave May   Vec                pos;
19340c453e9SDave May   PetscBool          error_check = swarm->migrate_error_on_missing_point;
19419307e5cSMatthew G. Knepley   const char       **coordFields;
195e4fbd051SBarry Smith   PetscMPIInt        size, rank;
19619307e5cSMatthew G. Knepley   const char        *cellid;
197480eef7bSDave May 
198521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1999566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(dm, &dmcell));
20028b400f6SJacob Faibussowitsch   PetscCheck(dmcell, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only valid if cell DM provided");
20119307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetCellDMActive(dm, &celldm));
20219307e5cSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetCellID(celldm, &cellid));
203480eef7bSDave May 
2049566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
2059566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
2067c6d1d28SDave May 
20743a82f2bSDave May #if 1
20843a82f2bSDave May   {
20943a82f2bSDave May     PetscInt     npoints_curr, range = 0;
21043a82f2bSDave May     PetscSFNode *sf_cells;
21143a82f2bSDave May 
2129566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints_curr, NULL, NULL));
2139566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(npoints_curr, &sf_cells));
21443a82f2bSDave May 
21519307e5cSMatthew G. Knepley     PetscCall(DMSwarmGetField(dm, cellid, NULL, NULL, (void **)&p_cellid));
21643a82f2bSDave May     for (p = 0; p < npoints_curr; p++) {
21743a82f2bSDave May       sf_cells[p].rank  = 0;
21843a82f2bSDave May       sf_cells[p].index = p_cellid[p];
219ad540459SPierre Jolivet       if (p_cellid[p] > range) range = p_cellid[p];
22043a82f2bSDave May     }
22119307e5cSMatthew G. Knepley     PetscCall(DMSwarmRestoreField(dm, cellid, NULL, NULL, (void **)&p_cellid));
22243a82f2bSDave May 
2239566063dSJacob Faibussowitsch     /* PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm),&sfcell)); */
2249566063dSJacob Faibussowitsch     PetscCall(PetscSFCreate(PETSC_COMM_SELF, &sfcell));
2259566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(sfcell, range, npoints_curr, NULL, PETSC_OWN_POINTER, sf_cells, PETSC_OWN_POINTER));
22643a82f2bSDave May   }
22743a82f2bSDave May #endif
22843a82f2bSDave May 
22919307e5cSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Nfc, &coordFields));
23019307e5cSMatthew G. Knepley   PetscCall(DMSwarmCreateLocalVectorFromFields(dm, Nfc, coordFields, &pos));
2319566063dSJacob Faibussowitsch   PetscCall(DMLocatePoints(dmcell, pos, DM_POINTLOCATION_NONE, &sfcell));
23219307e5cSMatthew G. Knepley   PetscCall(DMSwarmDestroyLocalVectorFromFields(dm, Nfc, coordFields, &pos));
233480eef7bSDave May 
23448a46eb9SPierre Jolivet   if (error_check) PetscCall(DMSwarmGetSize(dm, &npointsg));
2359566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
23619307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetField(dm, cellid, NULL, NULL, (void **)&p_cellid));
2379566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
2389566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell));
2391f70fafbSZach Atkins 
2401f70fafbSZach Atkins   for (p = 0; p < npoints; p++) {
2411f70fafbSZach Atkins     p_cellid[p] = LA_sfcell[p].index;
2421f70fafbSZach Atkins     rankval[p]  = rank;
2431f70fafbSZach Atkins   }
2449566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
24519307e5cSMatthew G. Knepley   PetscCall(DMSwarmRestoreField(dm, cellid, NULL, NULL, (void **)&p_cellid));
2469566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfcell));
247480eef7bSDave May 
248e4fbd051SBarry Smith   if (size > 1) {
2499566063dSJacob Faibussowitsch     PetscCall(DMSwarmMigrate_DMNeighborScatter(dm, dmcell, remove_sent_points, &npoints_prior_migration));
2506fbf25f8SDave May   } else {
25177048351SPatrick Sanan     DMSwarmDataField PField;
2520ed23c7fSDave May     PetscInt         npoints_curr;
2530ed23c7fSDave May 
2540ed23c7fSDave May     /* remove points which the domain */
25519307e5cSMatthew G. Knepley     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, cellid, &PField));
2561f70fafbSZach Atkins     PetscCall(DMSwarmDataFieldGetEntries(PField, (void **)&p_cellid));
2570ed23c7fSDave May 
2589566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints_curr, NULL, NULL));
2590ed23c7fSDave May     for (p = 0; p < npoints_curr; p++) {
2601f70fafbSZach Atkins       if (p_cellid[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
2610ed23c7fSDave May         /* kill point */
2629566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db, p));
2639566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints_curr, NULL, NULL)); /* you need to update npoints as the list size decreases! */
2641f70fafbSZach Atkins         PetscCall(DMSwarmDataFieldGetEntries(PField, (void **)&p_cellid));          /* update date point in case realloc performed */
2650ed23c7fSDave May         p--;                                                                        /* check replacement point */
2660ed23c7fSDave May       }
2670ed23c7fSDave May     }
2681f70fafbSZach Atkins     PetscCall(DMSwarmGetLocalSize(dm, &npoints_prior_migration));
2696fbf25f8SDave May   }
270480eef7bSDave May 
2712d4ee042Sprj-   /* locate points newly received */
2729566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints2, NULL, NULL));
273009b43efSDave May 
2745627991aSBarry Smith   { /* perform two point locations: (i) on the initial points set prior to communication; and (ii) on the new (received) points */
27519307e5cSMatthew G. Knepley     Vec              npos;
27619307e5cSMatthew G. Knepley     IS               nis;
277009b43efSDave May     PetscInt         npoints_from_neighbours, bs;
27877048351SPatrick Sanan     DMSwarmDataField PField;
279009b43efSDave May 
280009b43efSDave May     npoints_from_neighbours = npoints2 - npoints_prior_migration;
281009b43efSDave May 
28219307e5cSMatthew G. Knepley     PetscCall(DMSwarmCreateLocalVectorFromFields(dm, Nfc, coordFields, &pos));
28319307e5cSMatthew G. Knepley     PetscCall(VecGetBlockSize(pos, &bs));
28419307e5cSMatthew G. Knepley     PetscCall(ISCreateStride(PETSC_COMM_SELF, npoints_from_neighbours * bs, npoints_prior_migration * bs, 1, &nis));
285*953fc092SRylanor /*
286*953fc092SRylanor   Device VecGetSubVector to zero sized subvector triggers
287*953fc092SRylanor   debug error with mismatching memory types due to the device
288*953fc092SRylanor   pointer being host memtype without anything to point to in
289*953fc092SRylanor   Vec"Type"GetArrays(...), and we still need to pass something to
290*953fc092SRylanor   DMLocatePoints to avoid deadlock
291*953fc092SRylanor */
292*953fc092SRylanor #if defined(PETSC_HAVE_DEVICE)
293*953fc092SRylanor     if (npoints_from_neighbours > 0) {
29419307e5cSMatthew G. Knepley       PetscCall(VecGetSubVector(pos, nis, &npos));
29519307e5cSMatthew G. Knepley       PetscCall(DMLocatePoints(dmcell, npos, DM_POINTLOCATION_NONE, &sfcell));
29619307e5cSMatthew G. Knepley       PetscCall(VecRestoreSubVector(pos, nis, &npos));
297*953fc092SRylanor     } else {
298*953fc092SRylanor       PetscCall(VecCreate(PETSC_COMM_SELF, &npos));
299*953fc092SRylanor       PetscCall(VecSetSizes(npos, 0, PETSC_DETERMINE));
300*953fc092SRylanor       PetscCall(VecSetBlockSize(npos, bs));
301*953fc092SRylanor       PetscCall(VecSetType(npos, dm->vectype));
302*953fc092SRylanor       PetscCall(DMLocatePoints(dmcell, npos, DM_POINTLOCATION_NONE, &sfcell));
303*953fc092SRylanor       PetscCall(VecDestroy(&npos));
304*953fc092SRylanor     }
305*953fc092SRylanor #else
306*953fc092SRylanor     PetscCall(VecGetSubVector(pos, nis, &npos));
307*953fc092SRylanor     PetscCall(DMLocatePoints(dmcell, npos, DM_POINTLOCATION_NONE, &sfcell));
308*953fc092SRylanor     PetscCall(VecRestoreSubVector(pos, nis, &npos));
309*953fc092SRylanor #endif
31019307e5cSMatthew G. Knepley     PetscCall(ISDestroy(&nis));
31119307e5cSMatthew G. Knepley     PetscCall(DMSwarmDestroyLocalVectorFromFields(dm, Nfc, coordFields, &pos));
312009b43efSDave May 
3139566063dSJacob Faibussowitsch     PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell));
3149566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
31519307e5cSMatthew G. Knepley     PetscCall(DMSwarmGetField(dm, cellid, NULL, NULL, (void **)&p_cellid));
3161f70fafbSZach Atkins     for (p = 0; p < npoints_from_neighbours; p++) {
3171f70fafbSZach Atkins       rankval[npoints_prior_migration + p]  = rank;
3181f70fafbSZach Atkins       p_cellid[npoints_prior_migration + p] = LA_sfcell[p].index;
3191f70fafbSZach Atkins     }
3201f70fafbSZach Atkins 
32119307e5cSMatthew G. Knepley     PetscCall(DMSwarmRestoreField(dm, cellid, NULL, NULL, (void **)&p_cellid));
3229566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
323*953fc092SRylanor 
3249566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&sfcell));
325009b43efSDave May 
326009b43efSDave May     /* remove points which left processor */
32719307e5cSMatthew G. Knepley     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, cellid, &PField));
3281f70fafbSZach Atkins     PetscCall(DMSwarmDataFieldGetEntries(PField, (void **)&p_cellid));
329009b43efSDave May 
3309566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints2, NULL, NULL));
331009b43efSDave May     for (p = npoints_prior_migration; p < npoints2; p++) {
3321f70fafbSZach Atkins       if (p_cellid[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
333009b43efSDave May         /* kill point */
3349566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db, p));
3359566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints2, NULL, NULL)); /* you need to update npoints as the list size decreases! */
3361f70fafbSZach Atkins         PetscCall(DMSwarmDataFieldGetEntries(PField, (void **)&p_cellid));      /* update date point in case realloc performed */
337009b43efSDave May         p--;                                                                    /* check replacement point */
338009b43efSDave May       }
339009b43efSDave May     }
340009b43efSDave May   }
341009b43efSDave 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   }
3473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
348480eef7bSDave May }
349480eef7bSDave May 
350d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmMigrate_CellDMExact(DM dm, PetscBool remove_sent_points)
351d71ae5a4SJacob Faibussowitsch {
352521f74f9SMatthew G. Knepley   PetscFunctionBegin;
3533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
35408056efcSDave May }
35508056efcSDave May 
356480eef7bSDave May /*
357480eef7bSDave May  Redundant as this assumes points can only be sent to a single rank
358480eef7bSDave May */
359d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm, PetscInt *globalsize)
360d71ae5a4SJacob Faibussowitsch {
3612712d1f2SDave May   DM_Swarm     *swarm = (DM_Swarm *)dm->data;
36277048351SPatrick Sanan   DMSwarmDataEx de;
3632712d1f2SDave May   PetscInt      p, npoints, *rankval, n_points_recv;
3642712d1f2SDave May   PetscMPIInt   rank, nrank, negrank;
3652712d1f2SDave May   void         *point_buffer, *recv_points;
3662712d1f2SDave May   size_t        sizeof_dmswarm_point;
3672712d1f2SDave May 
368521f74f9SMatthew G. Knepley   PetscFunctionBegin;
3699566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
3709566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
3712712d1f2SDave May   *globalsize = npoints;
3729566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
3739566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm), 0, &de));
3749566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyInitialize(de));
3752712d1f2SDave May   for (p = 0; p < npoints; p++) {
376835f2295SStefano Zampini     PetscCall(PetscMPIIntCast(rankval[p], &negrank));
3772712d1f2SDave May     if (negrank < 0) {
3782712d1f2SDave May       nrank = -negrank - 1;
3799566063dSJacob Faibussowitsch       PetscCall(DMSwarmDataExTopologyAddNeighbour(de, nrank));
3802712d1f2SDave May     }
3812712d1f2SDave May   }
3829566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyFinalize(de));
3839566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExInitializeSendCount(de));
3842712d1f2SDave May   for (p = 0; p < npoints; p++) {
385835f2295SStefano Zampini     PetscCall(PetscMPIIntCast(rankval[p], &negrank));
3862712d1f2SDave May     if (negrank < 0) {
3872712d1f2SDave May       nrank = -negrank - 1;
3889566063dSJacob Faibussowitsch       PetscCall(DMSwarmDataExAddToSendCount(de, nrank, 1));
3892712d1f2SDave May     }
3902712d1f2SDave May   }
3919566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExFinalizeSendCount(de));
3929566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCreatePackedArray(swarm->db, &sizeof_dmswarm_point, &point_buffer));
3939566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackInitialize(de, sizeof_dmswarm_point));
3942712d1f2SDave May   for (p = 0; p < npoints; p++) {
395835f2295SStefano Zampini     PetscCall(PetscMPIIntCast(rankval[p], &negrank));
3962712d1f2SDave May     if (negrank < 0) {
3972712d1f2SDave May       nrank      = -negrank - 1;
3982712d1f2SDave May       rankval[p] = nrank;
3992712d1f2SDave May       /* copy point into buffer */
4009566063dSJacob Faibussowitsch       PetscCall(DMSwarmDataBucketFillPackedArray(swarm->db, p, point_buffer));
40177048351SPatrick Sanan       /* insert point buffer into DMSwarmDataExchanger */
4029566063dSJacob Faibussowitsch       PetscCall(DMSwarmDataExPackData(de, nrank, 1, point_buffer));
4032712d1f2SDave May       rankval[p] = negrank;
4042712d1f2SDave May     }
4052712d1f2SDave May   }
4069566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackFinalize(de));
4079566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
4089566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExBegin(de));
4099566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExEnd(de));
410835f2295SStefano Zampini   PetscCall(DMSwarmDataExGetRecvData(de, &n_points_recv, &recv_points));
4119566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
4129566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketSetSizes(swarm->db, npoints + n_points_recv, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
4132712d1f2SDave May   for (p = 0; p < n_points_recv; p++) {
4142712d1f2SDave May     void *data_p = (void *)((char *)recv_points + p * sizeof_dmswarm_point);
4152712d1f2SDave May 
4169566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketInsertPackedArray(swarm->db, npoints + p, data_p));
4172712d1f2SDave May   }
4189566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExView(de));
4199566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketDestroyPackedArray(swarm->db, &point_buffer));
4209566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExDestroy(de));
4213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4222712d1f2SDave May }
423b16650c8SDave May 
424b16650c8SDave May typedef struct {
425b16650c8SDave May   PetscMPIInt owner_rank;
426b16650c8SDave May   PetscReal   min[3], max[3];
427b16650c8SDave May } CollectBBox;
428b16650c8SDave May 
429d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMSwarmCollect_DMDABoundingBox(DM dm, PetscInt *globalsize)
430d71ae5a4SJacob Faibussowitsch {
431b16650c8SDave May   DM_Swarm          *swarm = (DM_Swarm *)dm->data;
43277048351SPatrick Sanan   DMSwarmDataEx      de;
433b16650c8SDave May   PetscInt           p, pk, npoints, *rankval, n_points_recv, n_bbox_recv, dim, neighbour_cells;
434b16650c8SDave May   PetscMPIInt        rank, nrank;
435b16650c8SDave May   void              *point_buffer, *recv_points;
436b16650c8SDave May   size_t             sizeof_dmswarm_point, sizeof_bbox_ctx;
437b16650c8SDave May   PetscBool          isdmda;
438b16650c8SDave May   CollectBBox       *bbox, *recv_bbox;
439b16650c8SDave May   const PetscMPIInt *dmneighborranks;
440b16650c8SDave May   DM                 dmcell;
441b16650c8SDave May 
442521f74f9SMatthew G. Knepley   PetscFunctionBegin;
4439566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
444b16650c8SDave May 
4459566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(dm, &dmcell));
44628b400f6SJacob Faibussowitsch   PetscCheck(dmcell, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only valid if cell DM provided");
447b16650c8SDave May   isdmda = PETSC_FALSE;
4483ba16761SJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dmcell, DMDA, &isdmda));
44928b400f6SJacob Faibussowitsch   PetscCheck(isdmda, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only DMDA support for CollectBoundingBox");
450b16650c8SDave May 
4519566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
452b16650c8SDave May   sizeof_bbox_ctx = sizeof(CollectBBox);
4533ba16761SJacob Faibussowitsch   PetscCall(PetscMalloc1(1, &bbox));
454b16650c8SDave May   bbox->owner_rank = rank;
455b16650c8SDave May 
456b16650c8SDave May   /* compute the bounding box based on the overlapping / stenctil size */
457b16650c8SDave May   {
458b16650c8SDave May     Vec lcoor;
459b16650c8SDave May 
4609566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(dmcell, &lcoor));
461fe39f135SDave May     if (dim >= 1) {
4629566063dSJacob Faibussowitsch       PetscCall(VecStrideMin(lcoor, 0, NULL, &bbox->min[0]));
4639566063dSJacob Faibussowitsch       PetscCall(VecStrideMax(lcoor, 0, NULL, &bbox->max[0]));
464fe39f135SDave May     }
465fe39f135SDave May     if (dim >= 2) {
4669566063dSJacob Faibussowitsch       PetscCall(VecStrideMin(lcoor, 1, NULL, &bbox->min[1]));
4679566063dSJacob Faibussowitsch       PetscCall(VecStrideMax(lcoor, 1, NULL, &bbox->max[1]));
468b16650c8SDave May     }
469fe39f135SDave May     if (dim == 3) {
4709566063dSJacob Faibussowitsch       PetscCall(VecStrideMin(lcoor, 2, NULL, &bbox->min[2]));
4719566063dSJacob Faibussowitsch       PetscCall(VecStrideMax(lcoor, 2, NULL, &bbox->max[2]));
472fe39f135SDave May     }
473fe39f135SDave May   }
4749566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
475b16650c8SDave May   *globalsize = npoints;
4769566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
4779566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm), 0, &de));
478b16650c8SDave May   /* use DMDA neighbours */
4799566063dSJacob Faibussowitsch   PetscCall(DMDAGetNeighbors(dmcell, &dmneighborranks));
4808dbd68bcSDave May   if (dim == 1) {
4818dbd68bcSDave May     neighbour_cells = 3;
4828dbd68bcSDave May   } else if (dim == 2) {
4838dbd68bcSDave May     neighbour_cells = 9;
4848dbd68bcSDave May   } else {
4858dbd68bcSDave May     neighbour_cells = 27;
4868dbd68bcSDave May   }
4879566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyInitialize(de));
488b16650c8SDave May   for (p = 0; p < neighbour_cells; p++) {
48948a46eb9SPierre Jolivet     if ((dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank)) PetscCall(DMSwarmDataExTopologyAddNeighbour(de, dmneighborranks[p]));
490b16650c8SDave May   }
4919566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyFinalize(de));
4929566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExInitializeSendCount(de));
493b16650c8SDave May   for (p = 0; p < neighbour_cells; p++) {
49448a46eb9SPierre Jolivet     if ((dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank)) PetscCall(DMSwarmDataExAddToSendCount(de, dmneighborranks[p], 1));
495b16650c8SDave May   }
4969566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExFinalizeSendCount(de));
497b16650c8SDave May   /* send bounding boxes */
4989566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackInitialize(de, sizeof_bbox_ctx));
499b16650c8SDave May   for (p = 0; p < neighbour_cells; p++) {
500b16650c8SDave May     nrank = dmneighborranks[p];
501b16650c8SDave May     if ((nrank >= 0) && (nrank != rank)) {
50277048351SPatrick Sanan       /* insert bbox buffer into DMSwarmDataExchanger */
5039566063dSJacob Faibussowitsch       PetscCall(DMSwarmDataExPackData(de, nrank, 1, bbox));
504b16650c8SDave May     }
505b16650c8SDave May   }
5069566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackFinalize(de));
507b16650c8SDave May   /* recv bounding boxes */
5089566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExBegin(de));
5099566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExEnd(de));
5109566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExGetRecvData(de, &n_bbox_recv, (void **)&recv_bbox));
511298827fbSBarry Smith   /*  Wrong, should not be using PETSC_COMM_WORLD */
512b16650c8SDave May   for (p = 0; p < n_bbox_recv; p++) {
5139371c9d4SSatish 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],
5149371c9d4SSatish Balay                                       (double)recv_bbox[p].max[1]));
515b16650c8SDave May   }
5169566063dSJacob Faibussowitsch   PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, stdout));
517b16650c8SDave May   /* of course this is stupid as this "generic" function should have a better way to know what the coordinates are called */
5189566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExInitializeSendCount(de));
519b16650c8SDave May   for (pk = 0; pk < n_bbox_recv; pk++) {
520b16650c8SDave May     PetscReal *array_x, *array_y;
521b16650c8SDave May 
5229566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dm, "coorx", NULL, NULL, (void **)&array_x));
5239566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dm, "coory", NULL, NULL, (void **)&array_y));
524b16650c8SDave May     for (p = 0; p < npoints; p++) {
525b16650c8SDave May       if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0])) {
52648a46eb9SPierre 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));
527b16650c8SDave May       }
528b16650c8SDave May     }
5299566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dm, "coory", NULL, NULL, (void **)&array_y));
5309566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dm, "coorx", NULL, NULL, (void **)&array_x));
531b16650c8SDave May   }
5329566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExFinalizeSendCount(de));
5339566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCreatePackedArray(swarm->db, &sizeof_dmswarm_point, &point_buffer));
5349566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackInitialize(de, sizeof_dmswarm_point));
535b16650c8SDave May   for (pk = 0; pk < n_bbox_recv; pk++) {
536b16650c8SDave May     PetscReal *array_x, *array_y;
537b16650c8SDave May 
5389566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dm, "coorx", NULL, NULL, (void **)&array_x));
5399566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dm, "coory", NULL, NULL, (void **)&array_y));
540b16650c8SDave May     for (p = 0; p < npoints; p++) {
541b16650c8SDave May       if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0])) {
542b16650c8SDave May         if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1])) {
543521f74f9SMatthew G. Knepley           /* copy point into buffer */
5449566063dSJacob Faibussowitsch           PetscCall(DMSwarmDataBucketFillPackedArray(swarm->db, p, point_buffer));
54577048351SPatrick Sanan           /* insert point buffer into DMSwarmDataExchanger */
5469566063dSJacob Faibussowitsch           PetscCall(DMSwarmDataExPackData(de, recv_bbox[pk].owner_rank, 1, point_buffer));
547b16650c8SDave May         }
548b16650c8SDave May       }
549b16650c8SDave May     }
5509566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dm, "coory", NULL, NULL, (void **)&array_y));
5519566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dm, "coorx", NULL, NULL, (void **)&array_x));
552b16650c8SDave May   }
5539566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackFinalize(de));
5549566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
5559566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExBegin(de));
5569566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExEnd(de));
557835f2295SStefano Zampini   PetscCall(DMSwarmDataExGetRecvData(de, &n_points_recv, &recv_points));
5589566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
5599566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketSetSizes(swarm->db, npoints + n_points_recv, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
560b16650c8SDave May   for (p = 0; p < n_points_recv; p++) {
561b16650c8SDave May     void *data_p = (void *)((char *)recv_points + p * sizeof_dmswarm_point);
562b16650c8SDave May 
5639566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketInsertPackedArray(swarm->db, npoints + p, data_p));
564b16650c8SDave May   }
5659566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketDestroyPackedArray(swarm->db, &point_buffer));
566e57ab8abSSatish Balay   PetscCall(PetscFree(bbox));
5679566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExView(de));
5689566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExDestroy(de));
5693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
570b16650c8SDave May }
571a9fd7477SDave May 
572a9fd7477SDave May /* General collection when no order, or neighbour information is provided */
573a9fd7477SDave May /*
574a9fd7477SDave May  User provides context and collect() method
575a9fd7477SDave May  Broadcast user context
576a9fd7477SDave May 
577a9fd7477SDave May  for each context / rank {
578a9fd7477SDave May    collect(swarm,context,n,list)
579a9fd7477SDave May  }
580a9fd7477SDave May */
581d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMSwarmCollect_General(DM dm, PetscErrorCode (*collect)(DM, void *, PetscInt *, PetscInt **), size_t ctx_size, void *ctx, PetscInt *globalsize)
582d71ae5a4SJacob Faibussowitsch {
583a9fd7477SDave May   DM_Swarm     *swarm = (DM_Swarm *)dm->data;
58477048351SPatrick Sanan   DMSwarmDataEx de;
5856497c311SBarry Smith   PetscInt      p, npoints, n_points_recv;
5866497c311SBarry Smith   PetscMPIInt   size, rank, len;
587a9fd7477SDave May   void         *point_buffer, *recv_points;
588a9fd7477SDave May   void         *ctxlist;
589a9fd7477SDave May   PetscInt     *n2collect, **collectlist;
590a9fd7477SDave May   size_t        sizeof_dmswarm_point;
591a9fd7477SDave May 
592521f74f9SMatthew G. Knepley   PetscFunctionBegin;
5939566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
5949566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
5959566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
596a9fd7477SDave May   *globalsize = npoints;
597a9fd7477SDave May   /* Broadcast user context */
5986497c311SBarry Smith   PetscCall(PetscMPIIntCast(ctx_size, &len));
5993ba16761SJacob Faibussowitsch   PetscCall(PetscMalloc(ctx_size * size, &ctxlist));
6006497c311SBarry Smith   PetscCallMPI(MPI_Allgather(ctx, len, MPI_CHAR, ctxlist, len, MPI_CHAR, PetscObjectComm((PetscObject)dm)));
6019566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size, &n2collect));
6029566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size, &collectlist));
6036497c311SBarry Smith   for (PetscMPIInt r = 0; r < size; r++) {
604a9fd7477SDave May     PetscInt  _n2collect;
605a9fd7477SDave May     PetscInt *_collectlist;
606a9fd7477SDave May     void     *_ctx_r;
607a9fd7477SDave May 
608a9fd7477SDave May     _n2collect   = 0;
609a9fd7477SDave May     _collectlist = NULL;
610a9fd7477SDave May     if (r != rank) { /* don't collect data from yourself */
611a9fd7477SDave May       _ctx_r = (void *)((char *)ctxlist + r * ctx_size);
6129566063dSJacob Faibussowitsch       PetscCall(collect(dm, _ctx_r, &_n2collect, &_collectlist));
613a9fd7477SDave May     }
614a9fd7477SDave May     n2collect[r]   = _n2collect;
615a9fd7477SDave May     collectlist[r] = _collectlist;
616a9fd7477SDave May   }
6179566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm), 0, &de));
618a9fd7477SDave May   /* Define topology */
6199566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyInitialize(de));
6206497c311SBarry Smith   for (PetscMPIInt r = 0; r < size; r++) {
621835f2295SStefano Zampini     if (n2collect[r] > 0) PetscCall(DMSwarmDataExTopologyAddNeighbour(de, r));
622a9fd7477SDave May   }
6239566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyFinalize(de));
624a9fd7477SDave May   /* Define send counts */
6259566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExInitializeSendCount(de));
6266497c311SBarry Smith   for (PetscMPIInt r = 0; r < size; r++) {
62748a46eb9SPierre Jolivet     if (n2collect[r] > 0) PetscCall(DMSwarmDataExAddToSendCount(de, r, n2collect[r]));
628a9fd7477SDave May   }
6299566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExFinalizeSendCount(de));
630a9fd7477SDave May   /* Pack data */
6319566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCreatePackedArray(swarm->db, &sizeof_dmswarm_point, &point_buffer));
6329566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackInitialize(de, sizeof_dmswarm_point));
6336497c311SBarry Smith   for (PetscMPIInt r = 0; r < size; r++) {
634a9fd7477SDave May     for (p = 0; p < n2collect[r]; p++) {
6359566063dSJacob Faibussowitsch       PetscCall(DMSwarmDataBucketFillPackedArray(swarm->db, collectlist[r][p], point_buffer));
636a9fd7477SDave May       /* insert point buffer into the data exchanger */
6379566063dSJacob Faibussowitsch       PetscCall(DMSwarmDataExPackData(de, r, 1, point_buffer));
638a9fd7477SDave May     }
639a9fd7477SDave May   }
6409566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackFinalize(de));
641a9fd7477SDave May   /* Scatter */
6429566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExBegin(de));
6439566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExEnd(de));
644a9fd7477SDave May   /* Collect data in DMSwarm container */
645835f2295SStefano Zampini   PetscCall(DMSwarmDataExGetRecvData(de, &n_points_recv, &recv_points));
6469566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
6479566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketSetSizes(swarm->db, npoints + n_points_recv, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
648a9fd7477SDave May   for (p = 0; p < n_points_recv; p++) {
649a9fd7477SDave May     void *data_p = (void *)((char *)recv_points + p * sizeof_dmswarm_point);
650a9fd7477SDave May 
6519566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketInsertPackedArray(swarm->db, npoints + p, data_p));
652a9fd7477SDave May   }
653a9fd7477SDave May   /* Release memory */
6546497c311SBarry Smith   for (PetscMPIInt r = 0; r < size; r++) {
655e57ab8abSSatish Balay     if (collectlist[r]) PetscCall(PetscFree(collectlist[r]));
656a9fd7477SDave May   }
6579566063dSJacob Faibussowitsch   PetscCall(PetscFree(collectlist));
6589566063dSJacob Faibussowitsch   PetscCall(PetscFree(n2collect));
6599566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctxlist));
6609566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketDestroyPackedArray(swarm->db, &point_buffer));
6619566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExView(de));
6629566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExDestroy(de));
6633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
664a9fd7477SDave May }
66547ce4f4bSMatthew G. Knepley 
66647ce4f4bSMatthew G. Knepley /*@
66747ce4f4bSMatthew G. Knepley   DMSwarmGetMigrateType - Get the style of point migration
66847ce4f4bSMatthew G. Knepley 
66920f4b53cSBarry Smith   Logically Collective
67047ce4f4bSMatthew G. Knepley 
67160225df5SJacob Faibussowitsch   Input Parameter:
67220f4b53cSBarry Smith . dm - the `DMSWARM`
67347ce4f4bSMatthew G. Knepley 
67460225df5SJacob Faibussowitsch   Output Parameter:
67520f4b53cSBarry Smith . mtype - The migration type, see `DMSwarmMigrateType`
67647ce4f4bSMatthew G. Knepley 
67747ce4f4bSMatthew G. Knepley   Level: intermediate
67847ce4f4bSMatthew G. Knepley 
67942747ad1SJacob Faibussowitsch .seealso: `DM`, `DMSWARM`, `DMSwarmMigrateType`, `DMSwarmMigrate()`
68047ce4f4bSMatthew G. Knepley @*/
68147ce4f4bSMatthew G. Knepley PetscErrorCode DMSwarmGetMigrateType(DM dm, DMSwarmMigrateType *mtype)
68247ce4f4bSMatthew G. Knepley {
68347ce4f4bSMatthew G. Knepley   PetscFunctionBegin;
68447ce4f4bSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
6854f572ea9SToby Isaac   PetscAssertPointer(mtype, 2);
68647ce4f4bSMatthew G. Knepley   *mtype = ((DM_Swarm *)dm->data)->migrate_type;
6873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
68847ce4f4bSMatthew G. Knepley }
68947ce4f4bSMatthew G. Knepley 
69047ce4f4bSMatthew G. Knepley /*@
69147ce4f4bSMatthew G. Knepley   DMSwarmSetMigrateType - Set the style of point migration
69247ce4f4bSMatthew G. Knepley 
69320f4b53cSBarry Smith   Logically Collective
69447ce4f4bSMatthew G. Knepley 
69560225df5SJacob Faibussowitsch   Input Parameters:
69620f4b53cSBarry Smith + dm    - the `DMSWARM`
69720f4b53cSBarry Smith - mtype - The migration type, see `DMSwarmMigrateType`
69847ce4f4bSMatthew G. Knepley 
69947ce4f4bSMatthew G. Knepley   Level: intermediate
70047ce4f4bSMatthew G. Knepley 
70160225df5SJacob Faibussowitsch .seealso: `DM`, `DMSWARM`, `DMSwarmMigrateType`, `DMSwarmGetMigrateType()`, `DMSwarmMigrate()`
70247ce4f4bSMatthew G. Knepley @*/
70347ce4f4bSMatthew G. Knepley PetscErrorCode DMSwarmSetMigrateType(DM dm, DMSwarmMigrateType mtype)
70447ce4f4bSMatthew G. Knepley {
70547ce4f4bSMatthew G. Knepley   PetscFunctionBegin;
70647ce4f4bSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
70747ce4f4bSMatthew G. Knepley   PetscValidLogicalCollectiveInt(dm, mtype, 2);
70847ce4f4bSMatthew G. Knepley   ((DM_Swarm *)dm->data)->migrate_type = mtype;
7093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
71047ce4f4bSMatthew G. Knepley }
711