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 */
DMSwarmMigrate_Push_Basic(DM dm,PetscBool remove_sent_points)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
DMSwarmMigrate_DMNeighborScatter(DM dm,DM dmcell,PetscBool remove_sent_points,PetscInt * npoints_prior_migration)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
DMSwarmMigrate_CellDMScatter(DM dm,PetscBool remove_sent_points)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));
259131df328Sjosephpu if (remove_sent_points) {
2600ed23c7fSDave May for (p = 0; p < npoints_curr; p++) {
2611f70fafbSZach Atkins if (p_cellid[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
2620ed23c7fSDave May /* kill point */
2639566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db, p));
2649566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints_curr, NULL, NULL)); /* you need to update npoints as the list size decreases! */
2651f70fafbSZach Atkins PetscCall(DMSwarmDataFieldGetEntries(PField, (void **)&p_cellid)); /* update date point in case realloc performed */
2660ed23c7fSDave May p--; /* check replacement point */
2670ed23c7fSDave May }
2680ed23c7fSDave May }
269131df328Sjosephpu }
2701f70fafbSZach Atkins PetscCall(DMSwarmGetLocalSize(dm, &npoints_prior_migration));
2716fbf25f8SDave May }
272480eef7bSDave May
2732d4ee042Sprj- /* locate points newly received */
2749566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints2, NULL, NULL));
275009b43efSDave May
2765627991aSBarry Smith { /* perform two point locations: (i) on the initial points set prior to communication; and (ii) on the new (received) points */
27719307e5cSMatthew G. Knepley Vec npos;
27819307e5cSMatthew G. Knepley IS nis;
279009b43efSDave May PetscInt npoints_from_neighbours, bs;
28077048351SPatrick Sanan DMSwarmDataField PField;
281009b43efSDave May
282009b43efSDave May npoints_from_neighbours = npoints2 - npoints_prior_migration;
283009b43efSDave May
28419307e5cSMatthew G. Knepley PetscCall(DMSwarmCreateLocalVectorFromFields(dm, Nfc, coordFields, &pos));
28519307e5cSMatthew G. Knepley PetscCall(VecGetBlockSize(pos, &bs));
28619307e5cSMatthew G. Knepley PetscCall(ISCreateStride(PETSC_COMM_SELF, npoints_from_neighbours * bs, npoints_prior_migration * bs, 1, &nis));
287953fc092SRylanor /*
288953fc092SRylanor Device VecGetSubVector to zero sized subvector triggers
289953fc092SRylanor debug error with mismatching memory types due to the device
290953fc092SRylanor pointer being host memtype without anything to point to in
291953fc092SRylanor Vec"Type"GetArrays(...), and we still need to pass something to
292953fc092SRylanor DMLocatePoints to avoid deadlock
293953fc092SRylanor */
294953fc092SRylanor #if defined(PETSC_HAVE_DEVICE)
295953fc092SRylanor if (npoints_from_neighbours > 0) {
29619307e5cSMatthew G. Knepley PetscCall(VecGetSubVector(pos, nis, &npos));
29719307e5cSMatthew G. Knepley PetscCall(DMLocatePoints(dmcell, npos, DM_POINTLOCATION_NONE, &sfcell));
29819307e5cSMatthew G. Knepley PetscCall(VecRestoreSubVector(pos, nis, &npos));
299953fc092SRylanor } else {
300953fc092SRylanor PetscCall(VecCreate(PETSC_COMM_SELF, &npos));
301953fc092SRylanor PetscCall(VecSetSizes(npos, 0, PETSC_DETERMINE));
302953fc092SRylanor PetscCall(VecSetBlockSize(npos, bs));
303953fc092SRylanor PetscCall(VecSetType(npos, dm->vectype));
304953fc092SRylanor PetscCall(DMLocatePoints(dmcell, npos, DM_POINTLOCATION_NONE, &sfcell));
305953fc092SRylanor PetscCall(VecDestroy(&npos));
306953fc092SRylanor }
307953fc092SRylanor #else
308953fc092SRylanor PetscCall(VecGetSubVector(pos, nis, &npos));
309953fc092SRylanor PetscCall(DMLocatePoints(dmcell, npos, DM_POINTLOCATION_NONE, &sfcell));
310953fc092SRylanor PetscCall(VecRestoreSubVector(pos, nis, &npos));
311953fc092SRylanor #endif
31219307e5cSMatthew G. Knepley PetscCall(ISDestroy(&nis));
31319307e5cSMatthew G. Knepley PetscCall(DMSwarmDestroyLocalVectorFromFields(dm, Nfc, coordFields, &pos));
314009b43efSDave May
3159566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell));
3169566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
31719307e5cSMatthew G. Knepley PetscCall(DMSwarmGetField(dm, cellid, NULL, NULL, (void **)&p_cellid));
3181f70fafbSZach Atkins for (p = 0; p < npoints_from_neighbours; p++) {
3191f70fafbSZach Atkins rankval[npoints_prior_migration + p] = rank;
3201f70fafbSZach Atkins p_cellid[npoints_prior_migration + p] = LA_sfcell[p].index;
3211f70fafbSZach Atkins }
3221f70fafbSZach Atkins
32319307e5cSMatthew G. Knepley PetscCall(DMSwarmRestoreField(dm, cellid, NULL, NULL, (void **)&p_cellid));
3249566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
325953fc092SRylanor
3269566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfcell));
327009b43efSDave May
328009b43efSDave May /* remove points which left processor */
32919307e5cSMatthew G. Knepley PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, cellid, &PField));
3301f70fafbSZach Atkins PetscCall(DMSwarmDataFieldGetEntries(PField, (void **)&p_cellid));
331009b43efSDave May
3329566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints2, NULL, NULL));
333131df328Sjosephpu if (remove_sent_points) {
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 }
344131df328Sjosephpu }
345009b43efSDave May
34640c453e9SDave May /* check for error on removed points */
34740c453e9SDave May if (error_check) {
3489566063dSJacob Faibussowitsch PetscCall(DMSwarmGetSize(dm, &npoints2g));
34963a3b9bcSJacob 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);
35040c453e9SDave May }
3513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
352480eef7bSDave May }
353480eef7bSDave May
DMSwarmMigrate_CellDMExact(DM dm,PetscBool remove_sent_points)354d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmMigrate_CellDMExact(DM dm, PetscBool remove_sent_points)
355d71ae5a4SJacob Faibussowitsch {
356521f74f9SMatthew G. Knepley PetscFunctionBegin;
3573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
35808056efcSDave May }
35908056efcSDave May
360480eef7bSDave May /*
361480eef7bSDave May Redundant as this assumes points can only be sent to a single rank
362480eef7bSDave May */
DMSwarmMigrate_GlobalToLocal_Basic(DM dm,PetscInt * globalsize)363d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm, PetscInt *globalsize)
364d71ae5a4SJacob Faibussowitsch {
3652712d1f2SDave May DM_Swarm *swarm = (DM_Swarm *)dm->data;
36677048351SPatrick Sanan DMSwarmDataEx de;
3672712d1f2SDave May PetscInt p, npoints, *rankval, n_points_recv;
3682712d1f2SDave May PetscMPIInt rank, nrank, negrank;
3692712d1f2SDave May void *point_buffer, *recv_points;
3702712d1f2SDave May size_t sizeof_dmswarm_point;
3712712d1f2SDave May
372521f74f9SMatthew G. Knepley PetscFunctionBegin;
3739566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
3749566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
3752712d1f2SDave May *globalsize = npoints;
3769566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
3779566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm), 0, &de));
3789566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExTopologyInitialize(de));
3792712d1f2SDave May for (p = 0; p < npoints; p++) {
380835f2295SStefano Zampini PetscCall(PetscMPIIntCast(rankval[p], &negrank));
3812712d1f2SDave May if (negrank < 0) {
3822712d1f2SDave May nrank = -negrank - 1;
3839566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExTopologyAddNeighbour(de, nrank));
3842712d1f2SDave May }
3852712d1f2SDave May }
3869566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExTopologyFinalize(de));
3879566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExInitializeSendCount(de));
3882712d1f2SDave May for (p = 0; p < npoints; p++) {
389835f2295SStefano Zampini PetscCall(PetscMPIIntCast(rankval[p], &negrank));
3902712d1f2SDave May if (negrank < 0) {
3912712d1f2SDave May nrank = -negrank - 1;
3929566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExAddToSendCount(de, nrank, 1));
3932712d1f2SDave May }
3942712d1f2SDave May }
3959566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExFinalizeSendCount(de));
3969566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketCreatePackedArray(swarm->db, &sizeof_dmswarm_point, &point_buffer));
3979566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExPackInitialize(de, sizeof_dmswarm_point));
3982712d1f2SDave May for (p = 0; p < npoints; p++) {
399835f2295SStefano Zampini PetscCall(PetscMPIIntCast(rankval[p], &negrank));
4002712d1f2SDave May if (negrank < 0) {
4012712d1f2SDave May nrank = -negrank - 1;
4022712d1f2SDave May rankval[p] = nrank;
4032712d1f2SDave May /* copy point into buffer */
4049566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketFillPackedArray(swarm->db, p, point_buffer));
40577048351SPatrick Sanan /* insert point buffer into DMSwarmDataExchanger */
4069566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExPackData(de, nrank, 1, point_buffer));
4072712d1f2SDave May rankval[p] = negrank;
4082712d1f2SDave May }
4092712d1f2SDave May }
4109566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExPackFinalize(de));
4119566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
4129566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExBegin(de));
4139566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExEnd(de));
414835f2295SStefano Zampini PetscCall(DMSwarmDataExGetRecvData(de, &n_points_recv, &recv_points));
4159566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
4169566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketSetSizes(swarm->db, npoints + n_points_recv, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
4172712d1f2SDave May for (p = 0; p < n_points_recv; p++) {
4182712d1f2SDave May void *data_p = (void *)((char *)recv_points + p * sizeof_dmswarm_point);
4192712d1f2SDave May
4209566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketInsertPackedArray(swarm->db, npoints + p, data_p));
4212712d1f2SDave May }
4229566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExView(de));
4239566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketDestroyPackedArray(swarm->db, &point_buffer));
4249566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExDestroy(de));
4253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4262712d1f2SDave May }
427b16650c8SDave May
428b16650c8SDave May typedef struct {
429b16650c8SDave May PetscMPIInt owner_rank;
430b16650c8SDave May PetscReal min[3], max[3];
431b16650c8SDave May } CollectBBox;
432b16650c8SDave May
DMSwarmCollect_DMDABoundingBox(DM dm,PetscInt * globalsize)433d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMSwarmCollect_DMDABoundingBox(DM dm, PetscInt *globalsize)
434d71ae5a4SJacob Faibussowitsch {
435b16650c8SDave May DM_Swarm *swarm = (DM_Swarm *)dm->data;
43677048351SPatrick Sanan DMSwarmDataEx de;
437b16650c8SDave May PetscInt p, pk, npoints, *rankval, n_points_recv, n_bbox_recv, dim, neighbour_cells;
438b16650c8SDave May PetscMPIInt rank, nrank;
439b16650c8SDave May void *point_buffer, *recv_points;
440b16650c8SDave May size_t sizeof_dmswarm_point, sizeof_bbox_ctx;
441b16650c8SDave May PetscBool isdmda;
442b16650c8SDave May CollectBBox *bbox, *recv_bbox;
443b16650c8SDave May const PetscMPIInt *dmneighborranks;
444b16650c8SDave May DM dmcell;
445b16650c8SDave May
446521f74f9SMatthew G. Knepley PetscFunctionBegin;
4479566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
448b16650c8SDave May
4499566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &dmcell));
45028b400f6SJacob Faibussowitsch PetscCheck(dmcell, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only valid if cell DM provided");
451b16650c8SDave May isdmda = PETSC_FALSE;
4523ba16761SJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dmcell, DMDA, &isdmda));
45328b400f6SJacob Faibussowitsch PetscCheck(isdmda, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only DMDA support for CollectBoundingBox");
454b16650c8SDave May
4559566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim));
456b16650c8SDave May sizeof_bbox_ctx = sizeof(CollectBBox);
4573ba16761SJacob Faibussowitsch PetscCall(PetscMalloc1(1, &bbox));
458b16650c8SDave May bbox->owner_rank = rank;
459b16650c8SDave May
460b16650c8SDave May /* compute the bounding box based on the overlapping / stenctil size */
461b16650c8SDave May {
462b16650c8SDave May Vec lcoor;
463b16650c8SDave May
4649566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dmcell, &lcoor));
465fe39f135SDave May if (dim >= 1) {
4669566063dSJacob Faibussowitsch PetscCall(VecStrideMin(lcoor, 0, NULL, &bbox->min[0]));
4679566063dSJacob Faibussowitsch PetscCall(VecStrideMax(lcoor, 0, NULL, &bbox->max[0]));
468fe39f135SDave May }
469fe39f135SDave May if (dim >= 2) {
4709566063dSJacob Faibussowitsch PetscCall(VecStrideMin(lcoor, 1, NULL, &bbox->min[1]));
4719566063dSJacob Faibussowitsch PetscCall(VecStrideMax(lcoor, 1, NULL, &bbox->max[1]));
472b16650c8SDave May }
473fe39f135SDave May if (dim == 3) {
4749566063dSJacob Faibussowitsch PetscCall(VecStrideMin(lcoor, 2, NULL, &bbox->min[2]));
4759566063dSJacob Faibussowitsch PetscCall(VecStrideMax(lcoor, 2, NULL, &bbox->max[2]));
476fe39f135SDave May }
477fe39f135SDave May }
4789566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
479b16650c8SDave May *globalsize = npoints;
4809566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
4819566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm), 0, &de));
482b16650c8SDave May /* use DMDA neighbours */
4839566063dSJacob Faibussowitsch PetscCall(DMDAGetNeighbors(dmcell, &dmneighborranks));
4848dbd68bcSDave May if (dim == 1) {
4858dbd68bcSDave May neighbour_cells = 3;
4868dbd68bcSDave May } else if (dim == 2) {
4878dbd68bcSDave May neighbour_cells = 9;
4888dbd68bcSDave May } else {
4898dbd68bcSDave May neighbour_cells = 27;
4908dbd68bcSDave May }
4919566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExTopologyInitialize(de));
492b16650c8SDave May for (p = 0; p < neighbour_cells; p++) {
49348a46eb9SPierre Jolivet if ((dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank)) PetscCall(DMSwarmDataExTopologyAddNeighbour(de, dmneighborranks[p]));
494b16650c8SDave May }
4959566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExTopologyFinalize(de));
4969566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExInitializeSendCount(de));
497b16650c8SDave May for (p = 0; p < neighbour_cells; p++) {
49848a46eb9SPierre Jolivet if ((dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank)) PetscCall(DMSwarmDataExAddToSendCount(de, dmneighborranks[p], 1));
499b16650c8SDave May }
5009566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExFinalizeSendCount(de));
501b16650c8SDave May /* send bounding boxes */
5029566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExPackInitialize(de, sizeof_bbox_ctx));
503b16650c8SDave May for (p = 0; p < neighbour_cells; p++) {
504b16650c8SDave May nrank = dmneighborranks[p];
505b16650c8SDave May if ((nrank >= 0) && (nrank != rank)) {
50677048351SPatrick Sanan /* insert bbox buffer into DMSwarmDataExchanger */
5079566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExPackData(de, nrank, 1, bbox));
508b16650c8SDave May }
509b16650c8SDave May }
5109566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExPackFinalize(de));
511b16650c8SDave May /* recv bounding boxes */
5129566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExBegin(de));
5139566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExEnd(de));
5149566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExGetRecvData(de, &n_bbox_recv, (void **)&recv_bbox));
515298827fbSBarry Smith /* Wrong, should not be using PETSC_COMM_WORLD */
516b16650c8SDave May for (p = 0; p < n_bbox_recv; p++) {
5179371c9d4SSatish 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],
5189371c9d4SSatish Balay (double)recv_bbox[p].max[1]));
519b16650c8SDave May }
5209566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, stdout));
521b16650c8SDave May /* of course this is stupid as this "generic" function should have a better way to know what the coordinates are called */
5229566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExInitializeSendCount(de));
523b16650c8SDave May for (pk = 0; pk < n_bbox_recv; pk++) {
524b16650c8SDave May PetscReal *array_x, *array_y;
525b16650c8SDave May
5269566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, "coorx", NULL, NULL, (void **)&array_x));
5279566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, "coory", NULL, NULL, (void **)&array_y));
528b16650c8SDave May for (p = 0; p < npoints; p++) {
529b16650c8SDave May if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0])) {
53048a46eb9SPierre 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));
531b16650c8SDave May }
532b16650c8SDave May }
5339566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, "coory", NULL, NULL, (void **)&array_y));
5349566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, "coorx", NULL, NULL, (void **)&array_x));
535b16650c8SDave May }
5369566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExFinalizeSendCount(de));
5379566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketCreatePackedArray(swarm->db, &sizeof_dmswarm_point, &point_buffer));
5389566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExPackInitialize(de, sizeof_dmswarm_point));
539b16650c8SDave May for (pk = 0; pk < n_bbox_recv; pk++) {
540b16650c8SDave May PetscReal *array_x, *array_y;
541b16650c8SDave May
5429566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, "coorx", NULL, NULL, (void **)&array_x));
5439566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, "coory", NULL, NULL, (void **)&array_y));
544b16650c8SDave May for (p = 0; p < npoints; p++) {
545b16650c8SDave May if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0])) {
546b16650c8SDave May if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1])) {
547521f74f9SMatthew G. Knepley /* copy point into buffer */
5489566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketFillPackedArray(swarm->db, p, point_buffer));
54977048351SPatrick Sanan /* insert point buffer into DMSwarmDataExchanger */
5509566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExPackData(de, recv_bbox[pk].owner_rank, 1, point_buffer));
551b16650c8SDave May }
552b16650c8SDave May }
553b16650c8SDave May }
5549566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, "coory", NULL, NULL, (void **)&array_y));
5559566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, "coorx", NULL, NULL, (void **)&array_x));
556b16650c8SDave May }
5579566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExPackFinalize(de));
5589566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
5599566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExBegin(de));
5609566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExEnd(de));
561835f2295SStefano Zampini PetscCall(DMSwarmDataExGetRecvData(de, &n_points_recv, &recv_points));
5629566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
5639566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketSetSizes(swarm->db, npoints + n_points_recv, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
564b16650c8SDave May for (p = 0; p < n_points_recv; p++) {
565b16650c8SDave May void *data_p = (void *)((char *)recv_points + p * sizeof_dmswarm_point);
566b16650c8SDave May
5679566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketInsertPackedArray(swarm->db, npoints + p, data_p));
568b16650c8SDave May }
5699566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketDestroyPackedArray(swarm->db, &point_buffer));
570e57ab8abSSatish Balay PetscCall(PetscFree(bbox));
5719566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExView(de));
5729566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExDestroy(de));
5733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
574b16650c8SDave May }
575a9fd7477SDave May
576a9fd7477SDave May /* General collection when no order, or neighbour information is provided */
577a9fd7477SDave May /*
578a9fd7477SDave May User provides context and collect() method
579a9fd7477SDave May Broadcast user context
580a9fd7477SDave May
581a9fd7477SDave May for each context / rank {
582a9fd7477SDave May collect(swarm,context,n,list)
583a9fd7477SDave May }
584a9fd7477SDave May */
DMSwarmCollect_General(DM dm,PetscErrorCode (* collect)(DM,void *,PetscInt *,PetscInt **),size_t ctx_size,PetscCtx ctx,PetscInt * globalsize)585*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMSwarmCollect_General(DM dm, PetscErrorCode (*collect)(DM, void *, PetscInt *, PetscInt **), size_t ctx_size, PetscCtx ctx, PetscInt *globalsize)
586d71ae5a4SJacob Faibussowitsch {
587a9fd7477SDave May DM_Swarm *swarm = (DM_Swarm *)dm->data;
58877048351SPatrick Sanan DMSwarmDataEx de;
5896497c311SBarry Smith PetscInt p, npoints, n_points_recv;
5906497c311SBarry Smith PetscMPIInt size, rank, len;
591a9fd7477SDave May void *point_buffer, *recv_points;
592a9fd7477SDave May void *ctxlist;
593a9fd7477SDave May PetscInt *n2collect, **collectlist;
594a9fd7477SDave May size_t sizeof_dmswarm_point;
595a9fd7477SDave May
596521f74f9SMatthew G. Knepley PetscFunctionBegin;
5979566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
5989566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
5999566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
600a9fd7477SDave May *globalsize = npoints;
601a9fd7477SDave May /* Broadcast user context */
6026497c311SBarry Smith PetscCall(PetscMPIIntCast(ctx_size, &len));
6033ba16761SJacob Faibussowitsch PetscCall(PetscMalloc(ctx_size * size, &ctxlist));
6046497c311SBarry Smith PetscCallMPI(MPI_Allgather(ctx, len, MPI_CHAR, ctxlist, len, MPI_CHAR, PetscObjectComm((PetscObject)dm)));
6059566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &n2collect));
6069566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &collectlist));
6076497c311SBarry Smith for (PetscMPIInt r = 0; r < size; r++) {
608a9fd7477SDave May PetscInt _n2collect;
609a9fd7477SDave May PetscInt *_collectlist;
610a9fd7477SDave May void *_ctx_r;
611a9fd7477SDave May
612a9fd7477SDave May _n2collect = 0;
613a9fd7477SDave May _collectlist = NULL;
614a9fd7477SDave May if (r != rank) { /* don't collect data from yourself */
615a9fd7477SDave May _ctx_r = (void *)((char *)ctxlist + r * ctx_size);
6169566063dSJacob Faibussowitsch PetscCall(collect(dm, _ctx_r, &_n2collect, &_collectlist));
617a9fd7477SDave May }
618a9fd7477SDave May n2collect[r] = _n2collect;
619a9fd7477SDave May collectlist[r] = _collectlist;
620a9fd7477SDave May }
6219566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExCreate(PetscObjectComm((PetscObject)dm), 0, &de));
622a9fd7477SDave May /* Define topology */
6239566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExTopologyInitialize(de));
6246497c311SBarry Smith for (PetscMPIInt r = 0; r < size; r++) {
625835f2295SStefano Zampini if (n2collect[r] > 0) PetscCall(DMSwarmDataExTopologyAddNeighbour(de, r));
626a9fd7477SDave May }
6279566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExTopologyFinalize(de));
628a9fd7477SDave May /* Define send counts */
6299566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExInitializeSendCount(de));
6306497c311SBarry Smith for (PetscMPIInt r = 0; r < size; r++) {
63148a46eb9SPierre Jolivet if (n2collect[r] > 0) PetscCall(DMSwarmDataExAddToSendCount(de, r, n2collect[r]));
632a9fd7477SDave May }
6339566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExFinalizeSendCount(de));
634a9fd7477SDave May /* Pack data */
6359566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketCreatePackedArray(swarm->db, &sizeof_dmswarm_point, &point_buffer));
6369566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExPackInitialize(de, sizeof_dmswarm_point));
6376497c311SBarry Smith for (PetscMPIInt r = 0; r < size; r++) {
638a9fd7477SDave May for (p = 0; p < n2collect[r]; p++) {
6399566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketFillPackedArray(swarm->db, collectlist[r][p], point_buffer));
640a9fd7477SDave May /* insert point buffer into the data exchanger */
6419566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExPackData(de, r, 1, point_buffer));
642a9fd7477SDave May }
643a9fd7477SDave May }
6449566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExPackFinalize(de));
645a9fd7477SDave May /* Scatter */
6469566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExBegin(de));
6479566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExEnd(de));
648a9fd7477SDave May /* Collect data in DMSwarm container */
649835f2295SStefano Zampini PetscCall(DMSwarmDataExGetRecvData(de, &n_points_recv, &recv_points));
6509566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
6519566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketSetSizes(swarm->db, npoints + n_points_recv, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
652a9fd7477SDave May for (p = 0; p < n_points_recv; p++) {
653a9fd7477SDave May void *data_p = (void *)((char *)recv_points + p * sizeof_dmswarm_point);
654a9fd7477SDave May
6559566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketInsertPackedArray(swarm->db, npoints + p, data_p));
656a9fd7477SDave May }
657a9fd7477SDave May /* Release memory */
6586497c311SBarry Smith for (PetscMPIInt r = 0; r < size; r++) {
659e57ab8abSSatish Balay if (collectlist[r]) PetscCall(PetscFree(collectlist[r]));
660a9fd7477SDave May }
6619566063dSJacob Faibussowitsch PetscCall(PetscFree(collectlist));
6629566063dSJacob Faibussowitsch PetscCall(PetscFree(n2collect));
6639566063dSJacob Faibussowitsch PetscCall(PetscFree(ctxlist));
6649566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketDestroyPackedArray(swarm->db, &point_buffer));
6659566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExView(de));
6669566063dSJacob Faibussowitsch PetscCall(DMSwarmDataExDestroy(de));
6673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
668a9fd7477SDave May }
66947ce4f4bSMatthew G. Knepley
67047ce4f4bSMatthew G. Knepley /*@
67147ce4f4bSMatthew G. Knepley DMSwarmGetMigrateType - Get the style of point migration
67247ce4f4bSMatthew G. Knepley
67320f4b53cSBarry Smith Logically Collective
67447ce4f4bSMatthew G. Knepley
67560225df5SJacob Faibussowitsch Input Parameter:
67620f4b53cSBarry Smith . dm - the `DMSWARM`
67747ce4f4bSMatthew G. Knepley
67860225df5SJacob Faibussowitsch Output Parameter:
67920f4b53cSBarry Smith . mtype - The migration type, see `DMSwarmMigrateType`
68047ce4f4bSMatthew G. Knepley
68147ce4f4bSMatthew G. Knepley Level: intermediate
68247ce4f4bSMatthew G. Knepley
68342747ad1SJacob Faibussowitsch .seealso: `DM`, `DMSWARM`, `DMSwarmMigrateType`, `DMSwarmMigrate()`
68447ce4f4bSMatthew G. Knepley @*/
DMSwarmGetMigrateType(DM dm,DMSwarmMigrateType * mtype)68547ce4f4bSMatthew G. Knepley PetscErrorCode DMSwarmGetMigrateType(DM dm, DMSwarmMigrateType *mtype)
68647ce4f4bSMatthew G. Knepley {
68747ce4f4bSMatthew G. Knepley PetscFunctionBegin;
68847ce4f4bSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
6894f572ea9SToby Isaac PetscAssertPointer(mtype, 2);
69047ce4f4bSMatthew G. Knepley *mtype = ((DM_Swarm *)dm->data)->migrate_type;
6913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
69247ce4f4bSMatthew G. Knepley }
69347ce4f4bSMatthew G. Knepley
69447ce4f4bSMatthew G. Knepley /*@
69547ce4f4bSMatthew G. Knepley DMSwarmSetMigrateType - Set the style of point migration
69647ce4f4bSMatthew G. Knepley
69720f4b53cSBarry Smith Logically Collective
69847ce4f4bSMatthew G. Knepley
69960225df5SJacob Faibussowitsch Input Parameters:
70020f4b53cSBarry Smith + dm - the `DMSWARM`
70120f4b53cSBarry Smith - mtype - The migration type, see `DMSwarmMigrateType`
70247ce4f4bSMatthew G. Knepley
70347ce4f4bSMatthew G. Knepley Level: intermediate
70447ce4f4bSMatthew G. Knepley
70560225df5SJacob Faibussowitsch .seealso: `DM`, `DMSWARM`, `DMSwarmMigrateType`, `DMSwarmGetMigrateType()`, `DMSwarmMigrate()`
70647ce4f4bSMatthew G. Knepley @*/
DMSwarmSetMigrateType(DM dm,DMSwarmMigrateType mtype)70747ce4f4bSMatthew G. Knepley PetscErrorCode DMSwarmSetMigrateType(DM dm, DMSwarmMigrateType mtype)
70847ce4f4bSMatthew G. Knepley {
70947ce4f4bSMatthew G. Knepley PetscFunctionBegin;
71047ce4f4bSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
71147ce4f4bSMatthew G. Knepley PetscValidLogicalCollectiveInt(dm, mtype, 2);
71247ce4f4bSMatthew G. Knepley ((DM_Swarm *)dm->data)->migrate_type = mtype;
7133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
71447ce4f4bSMatthew G. Knepley }
715