xref: /petsc/src/dm/impls/swarm/swarm_migrate.c (revision 6497c311e7b976d467be1503c1effce92a60525c)
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
10*6497c311SBarry Smith 
11*6497c311SBarry 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) {
31*6497c311SBarry Smith     nrank = (PetscMPIInt)rankval[p];
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++) {
37*6497c311SBarry Smith     nrank = (PetscMPIInt)rankval[p];
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++) {
44*6497c311SBarry Smith     nrank = (PetscMPIInt)rankval[p];
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++) {
65*6497c311SBarry Smith       nrank = (PetscMPIInt)rankval[p];
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));
839566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExGetRecvData(de, &n_points_recv, (void **)&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;
10077048351SPatrick Sanan   DMSwarmDataEx      de;
1011f70fafbSZach Atkins   PetscInt           r, p, npoints, *p_cellid, n_points_recv;
10240c453e9SDave May   PetscMPIInt        rank, _rank;
10340c453e9SDave May   const PetscMPIInt *neighbourranks;
10440c453e9SDave May   void              *point_buffer, *recv_points;
10540c453e9SDave May   size_t             sizeof_dmswarm_point;
10640c453e9SDave May   PetscInt           nneighbors;
1077c6d1d28SDave May   PetscMPIInt        mynneigh, *myneigh;
10840c453e9SDave May 
109521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1109566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
1119566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
1121f70fafbSZach Atkins   PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&p_cellid));
1139566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm), 0, &de));
1149566063dSJacob Faibussowitsch   PetscCall(DMGetNeighbors(dmcell, &nneighbors, &neighbourranks));
1159566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyInitialize(de));
11640c453e9SDave May   for (r = 0; r < nneighbors; r++) {
11740c453e9SDave May     _rank = neighbourranks[r];
11848a46eb9SPierre Jolivet     if ((_rank != rank) && (_rank > 0)) PetscCall(DMSwarmDataExTopologyAddNeighbour(de, _rank));
11940c453e9SDave May   }
1209566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyFinalize(de));
1219566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyGetNeighbours(de, &mynneigh, &myneigh));
1229566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExInitializeSendCount(de));
12340c453e9SDave May   for (p = 0; p < npoints; p++) {
1241f70fafbSZach Atkins     if (p_cellid[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
1257c6d1d28SDave May       for (r = 0; r < mynneigh; r++) {
1267c6d1d28SDave May         _rank = myneigh[r];
1279566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataExAddToSendCount(de, _rank, 1));
12840c453e9SDave May       }
12940c453e9SDave May     }
13040c453e9SDave May   }
1319566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExFinalizeSendCount(de));
1329566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCreatePackedArray(swarm->db, &sizeof_dmswarm_point, &point_buffer));
1339566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackInitialize(de, sizeof_dmswarm_point));
13440c453e9SDave May   for (p = 0; p < npoints; p++) {
1351f70fafbSZach Atkins     if (p_cellid[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
1367c6d1d28SDave May       for (r = 0; r < mynneigh; r++) {
1377c6d1d28SDave May         _rank = myneigh[r];
13840c453e9SDave May         /* copy point into buffer */
1399566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataBucketFillPackedArray(swarm->db, p, point_buffer));
14077048351SPatrick Sanan         /* insert point buffer into DMSwarmDataExchanger */
1419566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataExPackData(de, _rank, 1, point_buffer));
14240c453e9SDave May       }
14340c453e9SDave May     }
14440c453e9SDave May   }
1459566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackFinalize(de));
1461f70fafbSZach Atkins   PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&p_cellid));
14740c453e9SDave May   if (remove_sent_points) {
14877048351SPatrick Sanan     DMSwarmDataField PField;
1497c6d1d28SDave May 
1501f70fafbSZach Atkins     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, DMSwarmPICField_cellid, &PField));
1511f70fafbSZach Atkins     PetscCall(DMSwarmDataFieldGetEntries(PField, (void **)&p_cellid));
15240c453e9SDave May     /* remove points which left processor */
1539566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
15440c453e9SDave May     for (p = 0; p < npoints; p++) {
1551f70fafbSZach Atkins       if (p_cellid[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
15640c453e9SDave May         /* kill point */
1579566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db, p));
1589566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL)); /* you need to update npoints as the list size decreases! */
1591f70fafbSZach Atkins         PetscCall(DMSwarmDataFieldGetEntries(PField, (void **)&p_cellid));     /* update date point increase realloc performed */
16040c453e9SDave May         p--;                                                                   /* check replacement point */
16140c453e9SDave May       }
16240c453e9SDave May     }
16340c453e9SDave May   }
1649566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, npoints_prior_migration, NULL, NULL));
1659566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExBegin(de));
1669566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExEnd(de));
1679566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExGetRecvData(de, &n_points_recv, (void **)&recv_points));
1689566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
1699566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketSetSizes(swarm->db, npoints + n_points_recv, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
17040c453e9SDave May   for (p = 0; p < n_points_recv; p++) {
17140c453e9SDave May     void *data_p = (void *)((char *)recv_points + p * sizeof_dmswarm_point);
17240c453e9SDave May 
1739566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketInsertPackedArray(swarm->db, npoints + p, data_p));
17440c453e9SDave May   }
1759566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketDestroyPackedArray(swarm->db, &point_buffer));
1769566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExDestroy(de));
1773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17840c453e9SDave May }
179480eef7bSDave May 
180d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmMigrate_CellDMScatter(DM dm, PetscBool remove_sent_points)
181d71ae5a4SJacob Faibussowitsch {
182480eef7bSDave May   DM_Swarm          *swarm = (DM_Swarm *)dm->data;
1831f70fafbSZach Atkins   PetscInt           p, npoints, npointsg = 0, npoints2, npoints2g, *rankval, *p_cellid, npoints_prior_migration;
184bbe8250bSMatthew G. Knepley   PetscSF            sfcell = NULL;
185dfc14de9SMatthew G. Knepley   const PetscSFNode *LA_sfcell;
186480eef7bSDave May   DM                 dmcell;
187480eef7bSDave May   Vec                pos;
18840c453e9SDave May   PetscBool          error_check = swarm->migrate_error_on_missing_point;
189e4fbd051SBarry Smith   PetscMPIInt        size, rank;
190480eef7bSDave May 
191521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1929566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(dm, &dmcell));
19328b400f6SJacob Faibussowitsch   PetscCheck(dmcell, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only valid if cell DM provided");
194480eef7bSDave May 
1959566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
1969566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
1977c6d1d28SDave May 
19843a82f2bSDave May #if 1
19943a82f2bSDave May   {
20043a82f2bSDave May     PetscInt     npoints_curr, range = 0;
20143a82f2bSDave May     PetscSFNode *sf_cells;
20243a82f2bSDave May 
2039566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints_curr, NULL, NULL));
2049566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(npoints_curr, &sf_cells));
20543a82f2bSDave May 
2069566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&p_cellid));
20743a82f2bSDave May     for (p = 0; p < npoints_curr; p++) {
20843a82f2bSDave May       sf_cells[p].rank  = 0;
20943a82f2bSDave May       sf_cells[p].index = p_cellid[p];
210ad540459SPierre Jolivet       if (p_cellid[p] > range) range = p_cellid[p];
21143a82f2bSDave May     }
2129566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&p_cellid));
21343a82f2bSDave May 
2149566063dSJacob Faibussowitsch     /* PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm),&sfcell)); */
2159566063dSJacob Faibussowitsch     PetscCall(PetscSFCreate(PETSC_COMM_SELF, &sfcell));
2169566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(sfcell, range, npoints_curr, NULL, PETSC_OWN_POINTER, sf_cells, PETSC_OWN_POINTER));
21743a82f2bSDave May   }
21843a82f2bSDave May #endif
21943a82f2bSDave May 
2209566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateLocalVectorFromField(dm, DMSwarmPICField_coor, &pos));
2219566063dSJacob Faibussowitsch   PetscCall(DMLocatePoints(dmcell, pos, DM_POINTLOCATION_NONE, &sfcell));
2229566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyLocalVectorFromField(dm, DMSwarmPICField_coor, &pos));
223480eef7bSDave May 
22448a46eb9SPierre Jolivet   if (error_check) PetscCall(DMSwarmGetSize(dm, &npointsg));
2259566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
2261f70fafbSZach Atkins   PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&p_cellid));
2279566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
2289566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell));
2291f70fafbSZach Atkins 
2301f70fafbSZach Atkins   for (p = 0; p < npoints; p++) {
2311f70fafbSZach Atkins     p_cellid[p] = LA_sfcell[p].index;
2321f70fafbSZach Atkins     rankval[p]  = rank;
2331f70fafbSZach Atkins   }
2349566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
2351f70fafbSZach Atkins   PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&p_cellid));
2369566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfcell));
237480eef7bSDave May 
238e4fbd051SBarry Smith   if (size > 1) {
2399566063dSJacob Faibussowitsch     PetscCall(DMSwarmMigrate_DMNeighborScatter(dm, dmcell, remove_sent_points, &npoints_prior_migration));
2406fbf25f8SDave May   } else {
24177048351SPatrick Sanan     DMSwarmDataField PField;
2420ed23c7fSDave May     PetscInt         npoints_curr;
2430ed23c7fSDave May 
2440ed23c7fSDave May     /* remove points which the domain */
2451f70fafbSZach Atkins     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, DMSwarmPICField_cellid, &PField));
2461f70fafbSZach Atkins     PetscCall(DMSwarmDataFieldGetEntries(PField, (void **)&p_cellid));
2470ed23c7fSDave May 
2489566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints_curr, NULL, NULL));
2490ed23c7fSDave May     for (p = 0; p < npoints_curr; p++) {
2501f70fafbSZach Atkins       if (p_cellid[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
2510ed23c7fSDave May         /* kill point */
2529566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db, p));
2539566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints_curr, NULL, NULL)); /* you need to update npoints as the list size decreases! */
2541f70fafbSZach Atkins         PetscCall(DMSwarmDataFieldGetEntries(PField, (void **)&p_cellid));          /* update date point in case realloc performed */
2550ed23c7fSDave May         p--;                                                                        /* check replacement point */
2560ed23c7fSDave May       }
2570ed23c7fSDave May     }
2581f70fafbSZach Atkins     PetscCall(DMSwarmGetLocalSize(dm, &npoints_prior_migration));
2596fbf25f8SDave May   }
260480eef7bSDave May 
2612d4ee042Sprj-   /* locate points newly received */
2629566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints2, NULL, NULL));
263009b43efSDave May 
2647c6d1d28SDave May #if 0
2652d4ee042Sprj-   { /* safe alternative - however this performs two point locations on: (i) the initial points set and; (ii) the (initial + received) point set */
2667c6d1d28SDave May     PetscScalar *LA_coor;
2677c6d1d28SDave May     PetscInt bs;
26877048351SPatrick Sanan     DMSwarmDataField PField;
2697c6d1d28SDave May 
2709566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor));
2719566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,bs,bs*npoints2,(const PetscScalar*)LA_coor,&pos));
2729566063dSJacob Faibussowitsch     PetscCall(DMLocatePoints(dmcell,pos,DM_POINTLOCATION_NONE,&sfcell));
2737c6d1d28SDave May 
2749566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&pos));
2759566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor));
2767c6d1d28SDave May 
2779566063dSJacob Faibussowitsch     PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell));
2789566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval));
2797c6d1d28SDave May     for (p=0; p<npoints2; p++) {
280dfc14de9SMatthew G. Knepley       rankval[p] = LA_sfcell[p].index;
2817c6d1d28SDave May     }
2829566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&sfcell));
2839566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval));
28440c453e9SDave May 
2857c6d1d28SDave May     /* remove points which left processor */
2869566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,DMSwarmField_rank,&PField));
2879566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataFieldGetEntries(PField,(void**)&rankval));
2887c6d1d28SDave May 
2899566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints2,NULL,NULL));
2907c6d1d28SDave May     for (p=0; p<npoints2; p++) {
291f954cb40SDave May       if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
2927c6d1d28SDave May         /* kill point */
2939566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db,p));
2949566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataBucketGetSizes(swarm->db,&npoints2,NULL,NULL)); /* you need to update npoints as the list size decreases! */
2959566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataFieldGetEntries(PField,(void**)&rankval)); /* update date point increase realloc performed */
2967c6d1d28SDave May         p--; /* check replacement point */
2977c6d1d28SDave May       }
2987c6d1d28SDave May     }
29940c453e9SDave May   }
300009b43efSDave May #endif
301009b43efSDave May 
3025627991aSBarry Smith   { /* perform two point locations: (i) on the initial points set prior to communication; and (ii) on the new (received) points */
303009b43efSDave May     PetscScalar     *LA_coor;
304009b43efSDave May     PetscInt         npoints_from_neighbours, bs;
30577048351SPatrick Sanan     DMSwarmDataField PField;
306009b43efSDave May 
307009b43efSDave May     npoints_from_neighbours = npoints2 - npoints_prior_migration;
308009b43efSDave May 
3099566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, &bs, NULL, (void **)&LA_coor));
3109566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, bs, bs * npoints_from_neighbours, (const PetscScalar *)&LA_coor[bs * npoints_prior_migration], &pos));
311009b43efSDave May 
3129566063dSJacob Faibussowitsch     PetscCall(DMLocatePoints(dmcell, pos, DM_POINTLOCATION_NONE, &sfcell));
313009b43efSDave May 
3149566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&pos));
3159566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, &bs, NULL, (void **)&LA_coor));
316009b43efSDave May 
3179566063dSJacob Faibussowitsch     PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell));
3189566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
3191f70fafbSZach Atkins     PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&p_cellid));
3201f70fafbSZach Atkins     for (p = 0; p < npoints_from_neighbours; p++) {
3211f70fafbSZach Atkins       rankval[npoints_prior_migration + p]  = rank;
3221f70fafbSZach Atkins       p_cellid[npoints_prior_migration + p] = LA_sfcell[p].index;
3231f70fafbSZach Atkins     }
3241f70fafbSZach Atkins 
3251f70fafbSZach Atkins     PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&p_cellid));
3269566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
3279566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&sfcell));
328009b43efSDave May 
329009b43efSDave May     /* remove points which left processor */
3301f70fafbSZach Atkins     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, DMSwarmPICField_cellid, &PField));
3311f70fafbSZach Atkins     PetscCall(DMSwarmDataFieldGetEntries(PField, (void **)&p_cellid));
332009b43efSDave May 
3339566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints2, NULL, NULL));
334009b43efSDave May     for (p = npoints_prior_migration; p < npoints2; p++) {
3351f70fafbSZach Atkins       if (p_cellid[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
336009b43efSDave May         /* kill point */
3379566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db, p));
3389566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints2, NULL, NULL)); /* you need to update npoints as the list size decreases! */
3391f70fafbSZach Atkins         PetscCall(DMSwarmDataFieldGetEntries(PField, (void **)&p_cellid));      /* update date point in case realloc performed */
340009b43efSDave May         p--;                                                                    /* check replacement point */
341009b43efSDave May       }
342009b43efSDave May     }
343009b43efSDave May   }
344009b43efSDave May 
34540c453e9SDave May   /* check for error on removed points */
34640c453e9SDave May   if (error_check) {
3479566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dm, &npoints2g));
34863a3b9bcSJacob 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);
34940c453e9SDave May   }
3503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
351480eef7bSDave May }
352480eef7bSDave May 
353d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmMigrate_CellDMExact(DM dm, PetscBool remove_sent_points)
354d71ae5a4SJacob Faibussowitsch {
355521f74f9SMatthew G. Knepley   PetscFunctionBegin;
3563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
35708056efcSDave May }
35808056efcSDave May 
359480eef7bSDave May /*
360480eef7bSDave May  Redundant as this assumes points can only be sent to a single rank
361480eef7bSDave May */
362d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm, PetscInt *globalsize)
363d71ae5a4SJacob Faibussowitsch {
3642712d1f2SDave May   DM_Swarm     *swarm = (DM_Swarm *)dm->data;
36577048351SPatrick Sanan   DMSwarmDataEx de;
3662712d1f2SDave May   PetscInt      p, npoints, *rankval, n_points_recv;
3672712d1f2SDave May   PetscMPIInt   rank, nrank, negrank;
3682712d1f2SDave May   void         *point_buffer, *recv_points;
3692712d1f2SDave May   size_t        sizeof_dmswarm_point;
3702712d1f2SDave May 
371521f74f9SMatthew G. Knepley   PetscFunctionBegin;
3729566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
3739566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
3742712d1f2SDave May   *globalsize = npoints;
3759566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
3769566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm), 0, &de));
3779566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyInitialize(de));
3782712d1f2SDave May   for (p = 0; p < npoints; p++) {
379*6497c311SBarry Smith     negrank = (PetscMPIInt)rankval[p];
3802712d1f2SDave May     if (negrank < 0) {
3812712d1f2SDave May       nrank = -negrank - 1;
3829566063dSJacob Faibussowitsch       PetscCall(DMSwarmDataExTopologyAddNeighbour(de, nrank));
3832712d1f2SDave May     }
3842712d1f2SDave May   }
3859566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyFinalize(de));
3869566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExInitializeSendCount(de));
3872712d1f2SDave May   for (p = 0; p < npoints; p++) {
388*6497c311SBarry Smith     negrank = (PetscMPIInt)rankval[p];
3892712d1f2SDave May     if (negrank < 0) {
3902712d1f2SDave May       nrank = -negrank - 1;
3919566063dSJacob Faibussowitsch       PetscCall(DMSwarmDataExAddToSendCount(de, nrank, 1));
3922712d1f2SDave May     }
3932712d1f2SDave May   }
3949566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExFinalizeSendCount(de));
3959566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCreatePackedArray(swarm->db, &sizeof_dmswarm_point, &point_buffer));
3969566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackInitialize(de, sizeof_dmswarm_point));
3972712d1f2SDave May   for (p = 0; p < npoints; p++) {
398*6497c311SBarry Smith     negrank = (PetscMPIInt)rankval[p];
3992712d1f2SDave May     if (negrank < 0) {
4002712d1f2SDave May       nrank      = -negrank - 1;
4012712d1f2SDave May       rankval[p] = nrank;
4022712d1f2SDave May       /* copy point into buffer */
4039566063dSJacob Faibussowitsch       PetscCall(DMSwarmDataBucketFillPackedArray(swarm->db, p, point_buffer));
40477048351SPatrick Sanan       /* insert point buffer into DMSwarmDataExchanger */
4059566063dSJacob Faibussowitsch       PetscCall(DMSwarmDataExPackData(de, nrank, 1, point_buffer));
4062712d1f2SDave May       rankval[p] = negrank;
4072712d1f2SDave May     }
4082712d1f2SDave May   }
4099566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackFinalize(de));
4109566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
4119566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExBegin(de));
4129566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExEnd(de));
4139566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExGetRecvData(de, &n_points_recv, (void **)&recv_points));
4149566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
4159566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketSetSizes(swarm->db, npoints + n_points_recv, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
4162712d1f2SDave May   for (p = 0; p < n_points_recv; p++) {
4172712d1f2SDave May     void *data_p = (void *)((char *)recv_points + p * sizeof_dmswarm_point);
4182712d1f2SDave May 
4199566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketInsertPackedArray(swarm->db, npoints + p, data_p));
4202712d1f2SDave May   }
4219566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExView(de));
4229566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketDestroyPackedArray(swarm->db, &point_buffer));
4239566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExDestroy(de));
4243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4252712d1f2SDave May }
426b16650c8SDave May 
427b16650c8SDave May typedef struct {
428b16650c8SDave May   PetscMPIInt owner_rank;
429b16650c8SDave May   PetscReal   min[3], max[3];
430b16650c8SDave May } CollectBBox;
431b16650c8SDave May 
432d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMSwarmCollect_DMDABoundingBox(DM dm, PetscInt *globalsize)
433d71ae5a4SJacob Faibussowitsch {
434b16650c8SDave May   DM_Swarm          *swarm = (DM_Swarm *)dm->data;
43577048351SPatrick Sanan   DMSwarmDataEx      de;
436b16650c8SDave May   PetscInt           p, pk, npoints, *rankval, n_points_recv, n_bbox_recv, dim, neighbour_cells;
437b16650c8SDave May   PetscMPIInt        rank, nrank;
438b16650c8SDave May   void              *point_buffer, *recv_points;
439b16650c8SDave May   size_t             sizeof_dmswarm_point, sizeof_bbox_ctx;
440b16650c8SDave May   PetscBool          isdmda;
441b16650c8SDave May   CollectBBox       *bbox, *recv_bbox;
442b16650c8SDave May   const PetscMPIInt *dmneighborranks;
443b16650c8SDave May   DM                 dmcell;
444b16650c8SDave May 
445521f74f9SMatthew G. Knepley   PetscFunctionBegin;
4469566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
447b16650c8SDave May 
4489566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(dm, &dmcell));
44928b400f6SJacob Faibussowitsch   PetscCheck(dmcell, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only valid if cell DM provided");
450b16650c8SDave May   isdmda = PETSC_FALSE;
4513ba16761SJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dmcell, DMDA, &isdmda));
45228b400f6SJacob Faibussowitsch   PetscCheck(isdmda, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only DMDA support for CollectBoundingBox");
453b16650c8SDave May 
4549566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
455b16650c8SDave May   sizeof_bbox_ctx = sizeof(CollectBBox);
4563ba16761SJacob Faibussowitsch   PetscCall(PetscMalloc1(1, &bbox));
457b16650c8SDave May   bbox->owner_rank = rank;
458b16650c8SDave May 
459b16650c8SDave May   /* compute the bounding box based on the overlapping / stenctil size */
460b16650c8SDave May   {
461b16650c8SDave May     Vec lcoor;
462b16650c8SDave May 
4639566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(dmcell, &lcoor));
464fe39f135SDave May     if (dim >= 1) {
4659566063dSJacob Faibussowitsch       PetscCall(VecStrideMin(lcoor, 0, NULL, &bbox->min[0]));
4669566063dSJacob Faibussowitsch       PetscCall(VecStrideMax(lcoor, 0, NULL, &bbox->max[0]));
467fe39f135SDave May     }
468fe39f135SDave May     if (dim >= 2) {
4699566063dSJacob Faibussowitsch       PetscCall(VecStrideMin(lcoor, 1, NULL, &bbox->min[1]));
4709566063dSJacob Faibussowitsch       PetscCall(VecStrideMax(lcoor, 1, NULL, &bbox->max[1]));
471b16650c8SDave May     }
472fe39f135SDave May     if (dim == 3) {
4739566063dSJacob Faibussowitsch       PetscCall(VecStrideMin(lcoor, 2, NULL, &bbox->min[2]));
4749566063dSJacob Faibussowitsch       PetscCall(VecStrideMax(lcoor, 2, NULL, &bbox->max[2]));
475fe39f135SDave May     }
476fe39f135SDave May   }
4779566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
478b16650c8SDave May   *globalsize = npoints;
4799566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
4809566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm), 0, &de));
481b16650c8SDave May   /* use DMDA neighbours */
4829566063dSJacob Faibussowitsch   PetscCall(DMDAGetNeighbors(dmcell, &dmneighborranks));
4838dbd68bcSDave May   if (dim == 1) {
4848dbd68bcSDave May     neighbour_cells = 3;
4858dbd68bcSDave May   } else if (dim == 2) {
4868dbd68bcSDave May     neighbour_cells = 9;
4878dbd68bcSDave May   } else {
4888dbd68bcSDave May     neighbour_cells = 27;
4898dbd68bcSDave May   }
4909566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyInitialize(de));
491b16650c8SDave May   for (p = 0; p < neighbour_cells; p++) {
49248a46eb9SPierre Jolivet     if ((dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank)) PetscCall(DMSwarmDataExTopologyAddNeighbour(de, dmneighborranks[p]));
493b16650c8SDave May   }
4949566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyFinalize(de));
4959566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExInitializeSendCount(de));
496b16650c8SDave May   for (p = 0; p < neighbour_cells; p++) {
49748a46eb9SPierre Jolivet     if ((dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank)) PetscCall(DMSwarmDataExAddToSendCount(de, dmneighborranks[p], 1));
498b16650c8SDave May   }
4999566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExFinalizeSendCount(de));
500b16650c8SDave May   /* send bounding boxes */
5019566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackInitialize(de, sizeof_bbox_ctx));
502b16650c8SDave May   for (p = 0; p < neighbour_cells; p++) {
503b16650c8SDave May     nrank = dmneighborranks[p];
504b16650c8SDave May     if ((nrank >= 0) && (nrank != rank)) {
50577048351SPatrick Sanan       /* insert bbox buffer into DMSwarmDataExchanger */
5069566063dSJacob Faibussowitsch       PetscCall(DMSwarmDataExPackData(de, nrank, 1, bbox));
507b16650c8SDave May     }
508b16650c8SDave May   }
5099566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackFinalize(de));
510b16650c8SDave May   /* recv bounding boxes */
5119566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExBegin(de));
5129566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExEnd(de));
5139566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExGetRecvData(de, &n_bbox_recv, (void **)&recv_bbox));
514298827fbSBarry Smith   /*  Wrong, should not be using PETSC_COMM_WORLD */
515b16650c8SDave May   for (p = 0; p < n_bbox_recv; p++) {
5169371c9d4SSatish 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],
5179371c9d4SSatish Balay                                       (double)recv_bbox[p].max[1]));
518b16650c8SDave May   }
5199566063dSJacob Faibussowitsch   PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, stdout));
520b16650c8SDave May   /* of course this is stupid as this "generic" function should have a better way to know what the coordinates are called */
5219566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExInitializeSendCount(de));
522b16650c8SDave May   for (pk = 0; pk < n_bbox_recv; pk++) {
523b16650c8SDave May     PetscReal *array_x, *array_y;
524b16650c8SDave May 
5259566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dm, "coorx", NULL, NULL, (void **)&array_x));
5269566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dm, "coory", NULL, NULL, (void **)&array_y));
527b16650c8SDave May     for (p = 0; p < npoints; p++) {
528b16650c8SDave May       if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0])) {
52948a46eb9SPierre 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));
530b16650c8SDave May       }
531b16650c8SDave May     }
5329566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dm, "coory", NULL, NULL, (void **)&array_y));
5339566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dm, "coorx", NULL, NULL, (void **)&array_x));
534b16650c8SDave May   }
5359566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExFinalizeSendCount(de));
5369566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCreatePackedArray(swarm->db, &sizeof_dmswarm_point, &point_buffer));
5379566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackInitialize(de, sizeof_dmswarm_point));
538b16650c8SDave May   for (pk = 0; pk < n_bbox_recv; pk++) {
539b16650c8SDave May     PetscReal *array_x, *array_y;
540b16650c8SDave May 
5419566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dm, "coorx", NULL, NULL, (void **)&array_x));
5429566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dm, "coory", NULL, NULL, (void **)&array_y));
543b16650c8SDave May     for (p = 0; p < npoints; p++) {
544b16650c8SDave May       if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0])) {
545b16650c8SDave May         if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1])) {
546521f74f9SMatthew G. Knepley           /* copy point into buffer */
5479566063dSJacob Faibussowitsch           PetscCall(DMSwarmDataBucketFillPackedArray(swarm->db, p, point_buffer));
54877048351SPatrick Sanan           /* insert point buffer into DMSwarmDataExchanger */
5499566063dSJacob Faibussowitsch           PetscCall(DMSwarmDataExPackData(de, recv_bbox[pk].owner_rank, 1, point_buffer));
550b16650c8SDave May         }
551b16650c8SDave May       }
552b16650c8SDave May     }
5539566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dm, "coory", NULL, NULL, (void **)&array_y));
5549566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dm, "coorx", NULL, NULL, (void **)&array_x));
555b16650c8SDave May   }
5569566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackFinalize(de));
5579566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
5589566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExBegin(de));
5599566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExEnd(de));
5609566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExGetRecvData(de, &n_points_recv, (void **)&recv_points));
5619566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
5629566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketSetSizes(swarm->db, npoints + n_points_recv, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
563b16650c8SDave May   for (p = 0; p < n_points_recv; p++) {
564b16650c8SDave May     void *data_p = (void *)((char *)recv_points + p * sizeof_dmswarm_point);
565b16650c8SDave May 
5669566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketInsertPackedArray(swarm->db, npoints + p, data_p));
567b16650c8SDave May   }
5689566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketDestroyPackedArray(swarm->db, &point_buffer));
569e57ab8abSSatish Balay   PetscCall(PetscFree(bbox));
5709566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExView(de));
5719566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExDestroy(de));
5723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
573b16650c8SDave May }
574a9fd7477SDave May 
575a9fd7477SDave May /* General collection when no order, or neighbour information is provided */
576a9fd7477SDave May /*
577a9fd7477SDave May  User provides context and collect() method
578a9fd7477SDave May  Broadcast user context
579a9fd7477SDave May 
580a9fd7477SDave May  for each context / rank {
581a9fd7477SDave May    collect(swarm,context,n,list)
582a9fd7477SDave May  }
583a9fd7477SDave May */
584d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMSwarmCollect_General(DM dm, PetscErrorCode (*collect)(DM, void *, PetscInt *, PetscInt **), size_t ctx_size, void *ctx, PetscInt *globalsize)
585d71ae5a4SJacob Faibussowitsch {
586a9fd7477SDave May   DM_Swarm     *swarm = (DM_Swarm *)dm->data;
58777048351SPatrick Sanan   DMSwarmDataEx de;
588*6497c311SBarry Smith   PetscInt      p, npoints, n_points_recv;
589*6497c311SBarry Smith   PetscMPIInt   size, rank, len;
590a9fd7477SDave May   void         *point_buffer, *recv_points;
591a9fd7477SDave May   void         *ctxlist;
592a9fd7477SDave May   PetscInt     *n2collect, **collectlist;
593a9fd7477SDave May   size_t        sizeof_dmswarm_point;
594a9fd7477SDave May 
595521f74f9SMatthew G. Knepley   PetscFunctionBegin;
5969566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
5979566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
5989566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
599a9fd7477SDave May   *globalsize = npoints;
600a9fd7477SDave May   /* Broadcast user context */
601*6497c311SBarry Smith   PetscCall(PetscMPIIntCast(ctx_size, &len));
6023ba16761SJacob Faibussowitsch   PetscCall(PetscMalloc(ctx_size * size, &ctxlist));
603*6497c311SBarry Smith   PetscCallMPI(MPI_Allgather(ctx, len, MPI_CHAR, ctxlist, len, MPI_CHAR, PetscObjectComm((PetscObject)dm)));
6049566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size, &n2collect));
6059566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size, &collectlist));
606*6497c311SBarry Smith   for (PetscMPIInt r = 0; r < size; r++) {
607a9fd7477SDave May     PetscInt  _n2collect;
608a9fd7477SDave May     PetscInt *_collectlist;
609a9fd7477SDave May     void     *_ctx_r;
610a9fd7477SDave May 
611a9fd7477SDave May     _n2collect   = 0;
612a9fd7477SDave May     _collectlist = NULL;
613a9fd7477SDave May     if (r != rank) { /* don't collect data from yourself */
614a9fd7477SDave May       _ctx_r = (void *)((char *)ctxlist + r * ctx_size);
6159566063dSJacob Faibussowitsch       PetscCall(collect(dm, _ctx_r, &_n2collect, &_collectlist));
616a9fd7477SDave May     }
617a9fd7477SDave May     n2collect[r]   = _n2collect;
618a9fd7477SDave May     collectlist[r] = _collectlist;
619a9fd7477SDave May   }
6209566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm), 0, &de));
621a9fd7477SDave May   /* Define topology */
6229566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyInitialize(de));
623*6497c311SBarry Smith   for (PetscMPIInt r = 0; r < size; r++) {
62448a46eb9SPierre Jolivet     if (n2collect[r] > 0) PetscCall(DMSwarmDataExTopologyAddNeighbour(de, (PetscMPIInt)r));
625a9fd7477SDave May   }
6269566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExTopologyFinalize(de));
627a9fd7477SDave May   /* Define send counts */
6289566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExInitializeSendCount(de));
629*6497c311SBarry Smith   for (PetscMPIInt r = 0; r < size; r++) {
63048a46eb9SPierre Jolivet     if (n2collect[r] > 0) PetscCall(DMSwarmDataExAddToSendCount(de, r, n2collect[r]));
631a9fd7477SDave May   }
6329566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExFinalizeSendCount(de));
633a9fd7477SDave May   /* Pack data */
6349566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCreatePackedArray(swarm->db, &sizeof_dmswarm_point, &point_buffer));
6359566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackInitialize(de, sizeof_dmswarm_point));
636*6497c311SBarry Smith   for (PetscMPIInt r = 0; r < size; r++) {
637a9fd7477SDave May     for (p = 0; p < n2collect[r]; p++) {
6389566063dSJacob Faibussowitsch       PetscCall(DMSwarmDataBucketFillPackedArray(swarm->db, collectlist[r][p], point_buffer));
639a9fd7477SDave May       /* insert point buffer into the data exchanger */
6409566063dSJacob Faibussowitsch       PetscCall(DMSwarmDataExPackData(de, r, 1, point_buffer));
641a9fd7477SDave May     }
642a9fd7477SDave May   }
6439566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExPackFinalize(de));
644a9fd7477SDave May   /* Scatter */
6459566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExBegin(de));
6469566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExEnd(de));
647a9fd7477SDave May   /* Collect data in DMSwarm container */
6489566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExGetRecvData(de, &n_points_recv, (void **)&recv_points));
6499566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
6509566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketSetSizes(swarm->db, npoints + n_points_recv, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
651a9fd7477SDave May   for (p = 0; p < n_points_recv; p++) {
652a9fd7477SDave May     void *data_p = (void *)((char *)recv_points + p * sizeof_dmswarm_point);
653a9fd7477SDave May 
6549566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketInsertPackedArray(swarm->db, npoints + p, data_p));
655a9fd7477SDave May   }
656a9fd7477SDave May   /* Release memory */
657*6497c311SBarry Smith   for (PetscMPIInt r = 0; r < size; r++) {
658e57ab8abSSatish Balay     if (collectlist[r]) PetscCall(PetscFree(collectlist[r]));
659a9fd7477SDave May   }
6609566063dSJacob Faibussowitsch   PetscCall(PetscFree(collectlist));
6619566063dSJacob Faibussowitsch   PetscCall(PetscFree(n2collect));
6629566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctxlist));
6639566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketDestroyPackedArray(swarm->db, &point_buffer));
6649566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExView(de));
6659566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataExDestroy(de));
6663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
667a9fd7477SDave May }
66847ce4f4bSMatthew G. Knepley 
66947ce4f4bSMatthew G. Knepley /*@
67047ce4f4bSMatthew G. Knepley   DMSwarmGetMigrateType - Get the style of point migration
67147ce4f4bSMatthew G. Knepley 
67220f4b53cSBarry Smith   Logically Collective
67347ce4f4bSMatthew G. Knepley 
67460225df5SJacob Faibussowitsch   Input Parameter:
67520f4b53cSBarry Smith . dm - the `DMSWARM`
67647ce4f4bSMatthew G. Knepley 
67760225df5SJacob Faibussowitsch   Output Parameter:
67820f4b53cSBarry Smith . mtype - The migration type, see `DMSwarmMigrateType`
67947ce4f4bSMatthew G. Knepley 
68047ce4f4bSMatthew G. Knepley   Level: intermediate
68147ce4f4bSMatthew G. Knepley 
68242747ad1SJacob Faibussowitsch .seealso: `DM`, `DMSWARM`, `DMSwarmMigrateType`, `DMSwarmMigrate()`
68347ce4f4bSMatthew G. Knepley @*/
68447ce4f4bSMatthew G. Knepley PetscErrorCode DMSwarmGetMigrateType(DM dm, DMSwarmMigrateType *mtype)
68547ce4f4bSMatthew G. Knepley {
68647ce4f4bSMatthew G. Knepley   PetscFunctionBegin;
68747ce4f4bSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
6884f572ea9SToby Isaac   PetscAssertPointer(mtype, 2);
68947ce4f4bSMatthew G. Knepley   *mtype = ((DM_Swarm *)dm->data)->migrate_type;
6903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
69147ce4f4bSMatthew G. Knepley }
69247ce4f4bSMatthew G. Knepley 
69347ce4f4bSMatthew G. Knepley /*@
69447ce4f4bSMatthew G. Knepley   DMSwarmSetMigrateType - Set the style of point migration
69547ce4f4bSMatthew G. Knepley 
69620f4b53cSBarry Smith   Logically Collective
69747ce4f4bSMatthew G. Knepley 
69860225df5SJacob Faibussowitsch   Input Parameters:
69920f4b53cSBarry Smith + dm    - the `DMSWARM`
70020f4b53cSBarry Smith - mtype - The migration type, see `DMSwarmMigrateType`
70147ce4f4bSMatthew G. Knepley 
70247ce4f4bSMatthew G. Knepley   Level: intermediate
70347ce4f4bSMatthew G. Knepley 
70460225df5SJacob Faibussowitsch .seealso: `DM`, `DMSWARM`, `DMSwarmMigrateType`, `DMSwarmGetMigrateType()`, `DMSwarmMigrate()`
70547ce4f4bSMatthew G. Knepley @*/
70647ce4f4bSMatthew G. Knepley PetscErrorCode DMSwarmSetMigrateType(DM dm, DMSwarmMigrateType mtype)
70747ce4f4bSMatthew G. Knepley {
70847ce4f4bSMatthew G. Knepley   PetscFunctionBegin;
70947ce4f4bSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
71047ce4f4bSMatthew G. Knepley   PetscValidLogicalCollectiveInt(dm, mtype, 2);
71147ce4f4bSMatthew G. Knepley   ((DM_Swarm *)dm->data)->migrate_type = mtype;
7123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
71347ce4f4bSMatthew G. Knepley }
714