xref: /petsc/src/dm/impls/swarm/swarm_migrate.c (revision 3151f1d57dd52c070a1341b62bca41e71b9c6066)
1df21e3a8SDave May 
2dfc14de9SMatthew G. Knepley #include <petscsf.h>
308056efcSDave May #include <petscdmswarm.h>
4df21e3a8SDave May #include <petsc/private/dmswarmimpl.h>    /*I   "petscdmswarm.h"   I*/
5df21e3a8SDave May #include "data_bucket.h"
6df21e3a8SDave May #include "data_ex.h"
7df21e3a8SDave May 
8df21e3a8SDave May 
9480eef7bSDave May /*
10480eef7bSDave May  User loads desired location (MPI rank) into field DMSwarm_rank
11480eef7bSDave May */
12df21e3a8SDave May PetscErrorCode DMSwarmMigrate_Push_Basic(DM dm,PetscBool remove_sent_points)
13df21e3a8SDave May {
14df21e3a8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
15df21e3a8SDave May   PetscErrorCode ierr;
16df21e3a8SDave May   DataEx 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;
21df21e3a8SDave May 
22521f74f9SMatthew G. Knepley   PetscFunctionBegin;
23df21e3a8SDave May   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
24df21e3a8SDave May 
25df21e3a8SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
2608056efcSDave May   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
27521f74f9SMatthew G. Knepley   ierr = DataExCreate(PetscObjectComm((PetscObject)dm),0, &de);CHKERRQ(ierr);
28df21e3a8SDave May   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
29521f74f9SMatthew G. Knepley   for (p = 0; p < npoints; ++p) {
30df21e3a8SDave May     nrank = rankval[p];
31df21e3a8SDave May     if (nrank != rank) {
32df21e3a8SDave May       ierr = DataExTopologyAddNeighbour(de,nrank);CHKERRQ(ierr);
33df21e3a8SDave May     }
34df21e3a8SDave May   }
35df21e3a8SDave May   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
36df21e3a8SDave May   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
37df21e3a8SDave May   for (p=0; p<npoints; p++) {
38df21e3a8SDave May     nrank = rankval[p];
39df21e3a8SDave May     if (nrank != rank) {
40df21e3a8SDave May       ierr = DataExAddToSendCount(de,nrank,1);CHKERRQ(ierr);
41df21e3a8SDave May     }
42df21e3a8SDave May   }
43df21e3a8SDave May   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
44df21e3a8SDave May   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
45df21e3a8SDave May   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
46df21e3a8SDave May   for (p=0; p<npoints; p++) {
47df21e3a8SDave May     nrank = rankval[p];
48df21e3a8SDave May     if (nrank != rank) {
49df21e3a8SDave May       /* copy point into buffer */
50df21e3a8SDave May       ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr);
51df21e3a8SDave May       /* insert point buffer into DataExchanger */
52df21e3a8SDave May       ierr = DataExPackData(de,nrank,1,point_buffer);CHKERRQ(ierr);
53df21e3a8SDave May     }
54df21e3a8SDave May   }
55df21e3a8SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
5622a417f9SDave May   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
57df21e3a8SDave May 
58df21e3a8SDave May   if (remove_sent_points) {
5922a417f9SDave May     DataField gfield;
6022a417f9SDave May 
6122a417f9SDave May     ierr = DataBucketGetDataFieldByName(swarm->db,DMSwarmField_rank,&gfield);CHKERRQ(ierr);
6222a417f9SDave May     ierr = DataFieldGetAccess(gfield);CHKERRQ(ierr);
6322a417f9SDave May     ierr = DataFieldGetEntries(gfield,(void**)&rankval);CHKERRQ(ierr);
6422a417f9SDave May 
65df21e3a8SDave May     /* remove points which left processor */
66df21e3a8SDave May     ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
67df21e3a8SDave May     for (p=0; p<npoints; p++) {
68df21e3a8SDave May       nrank = rankval[p];
69df21e3a8SDave May       if (nrank != rank) {
70df21e3a8SDave May         /* kill point */
7122a417f9SDave May         ierr = DataFieldRestoreAccess(gfield);CHKERRQ(ierr);
7222a417f9SDave May 
73df21e3a8SDave May         ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr);
7422a417f9SDave May 
7522a417f9SDave May         ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */
7622a417f9SDave May         ierr = DataFieldGetAccess(gfield);CHKERRQ(ierr);
7722a417f9SDave May         ierr = DataFieldGetEntries(gfield,(void**)&rankval);CHKERRQ(ierr);
78df21e3a8SDave May         p--; /* check replacement point */
79df21e3a8SDave May       }
80df21e3a8SDave May     }
8122a417f9SDave May     ierr = DataFieldRestoreEntries(gfield,(void**)&rankval);CHKERRQ(ierr);
8222a417f9SDave May     ierr = DataFieldRestoreAccess(gfield);CHKERRQ(ierr);
83df21e3a8SDave May   }
84df21e3a8SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
85df21e3a8SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
86df21e3a8SDave May   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
87df21e3a8SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
8865141ba8SDave May   ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr);
89df21e3a8SDave May   for (p=0; p<n_points_recv; p++) {
90df21e3a8SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
91df21e3a8SDave May 
92df21e3a8SDave May     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
93df21e3a8SDave May   }
94df21e3a8SDave May   ierr = DataExView(de);CHKERRQ(ierr);
95df21e3a8SDave May   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
96df21e3a8SDave May   ierr = DataExDestroy(de);CHKERRQ(ierr);
97df21e3a8SDave May   PetscFunctionReturn(0);
98df21e3a8SDave May }
992712d1f2SDave May 
100889dbfe5SDave May PetscErrorCode DMSwarmMigrate_DMNeighborScatter(DM dm,DM dmcell,PetscBool remove_sent_points,PetscInt *npoints_prior_migration)
10140c453e9SDave May {
10240c453e9SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
10340c453e9SDave May   PetscErrorCode ierr;
10440c453e9SDave May   DataEx de;
10540c453e9SDave May   PetscInt r,p,npoints,*rankval,n_points_recv;
10640c453e9SDave May   PetscMPIInt rank,_rank;
10740c453e9SDave May   const PetscMPIInt *neighbourranks;
10840c453e9SDave May   void *point_buffer,*recv_points;
10940c453e9SDave May   size_t sizeof_dmswarm_point;
11040c453e9SDave May   PetscInt nneighbors;
1117c6d1d28SDave May   PetscMPIInt mynneigh,*myneigh;
11240c453e9SDave May 
113521f74f9SMatthew G. Knepley   PetscFunctionBegin;
11440c453e9SDave May   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
11540c453e9SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
11608056efcSDave May   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
117521f74f9SMatthew G. Knepley   ierr = DataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr);
11840c453e9SDave May   ierr = DMGetNeighbors(dmcell,&nneighbors,&neighbourranks);CHKERRQ(ierr);
11940c453e9SDave May   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
12040c453e9SDave May   for (r=0; r<nneighbors; r++) {
12140c453e9SDave May     _rank = neighbourranks[r];
12240c453e9SDave May     if ((_rank != rank) && (_rank > 0)) {
12340c453e9SDave May       ierr = DataExTopologyAddNeighbour(de,_rank);CHKERRQ(ierr);
12440c453e9SDave May     }
12540c453e9SDave May   }
12640c453e9SDave May   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
1277c6d1d28SDave May   ierr = DataExTopologyGetNeighbours(de,&mynneigh,&myneigh);CHKERRQ(ierr);
12840c453e9SDave May   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
12940c453e9SDave May   for (p=0; p<npoints; p++) {
130f954cb40SDave May     if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
1317c6d1d28SDave May       for (r=0; r<mynneigh; r++) {
1327c6d1d28SDave May         _rank = myneigh[r];
13340c453e9SDave May         ierr = DataExAddToSendCount(de,_rank,1);CHKERRQ(ierr);
13440c453e9SDave May       }
13540c453e9SDave May     }
13640c453e9SDave May   }
13740c453e9SDave May   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
13840c453e9SDave May   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
13940c453e9SDave May   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
14040c453e9SDave May   for (p=0; p<npoints; p++) {
141f954cb40SDave May     if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
1427c6d1d28SDave May       for (r=0; r<mynneigh; r++) {
1437c6d1d28SDave May         _rank = myneigh[r];
14440c453e9SDave May         /* copy point into buffer */
14540c453e9SDave May         ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr);
14640c453e9SDave May         /* insert point buffer into DataExchanger */
14740c453e9SDave May         ierr = DataExPackData(de,_rank,1,point_buffer);CHKERRQ(ierr);
14840c453e9SDave May       }
14940c453e9SDave May     }
15040c453e9SDave May   }
15140c453e9SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
1527c6d1d28SDave May   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
15340c453e9SDave May   if (remove_sent_points) {
1547c6d1d28SDave May     DataField PField;
1557c6d1d28SDave May 
1567c6d1d28SDave May     ierr = DataBucketGetDataFieldByName(swarm->db,DMSwarmField_rank,&PField);CHKERRQ(ierr);
1577c6d1d28SDave May     ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr);
15840c453e9SDave May     /* remove points which left processor */
15940c453e9SDave May     ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
16040c453e9SDave May     for (p=0; p<npoints; p++) {
161f954cb40SDave May       if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
16240c453e9SDave May         /* kill point */
16340c453e9SDave May         ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr);
1647c6d1d28SDave May         ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */
1657c6d1d28SDave May         ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); /* update date point increase realloc performed */
16640c453e9SDave May         p--; /* check replacement point */
16740c453e9SDave May       }
16840c453e9SDave May     }
16940c453e9SDave May   }
17040c453e9SDave May   ierr = DataBucketGetSizes(swarm->db,npoints_prior_migration,NULL,NULL);CHKERRQ(ierr);
17140c453e9SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
17240c453e9SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
17340c453e9SDave May   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
17440c453e9SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
17565141ba8SDave May   ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr);
17640c453e9SDave May   for (p=0; p<n_points_recv; p++) {
17740c453e9SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
17840c453e9SDave May 
17940c453e9SDave May     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
18040c453e9SDave May   }
18140c453e9SDave May   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
18240c453e9SDave May   ierr = DataExDestroy(de);CHKERRQ(ierr);
18340c453e9SDave May   PetscFunctionReturn(0);
18440c453e9SDave May }
185480eef7bSDave May 
18608056efcSDave May PetscErrorCode DMSwarmMigrate_CellDMScatter(DM dm,PetscBool remove_sent_points)
187480eef7bSDave May {
188480eef7bSDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
189480eef7bSDave May   PetscErrorCode ierr;
1906fbf25f8SDave May   PetscInt p,npoints,npointsg=0,npoints2,npoints2g,*rankval,npoints_prior_migration;
191bbe8250bSMatthew G. Knepley   PetscSF sfcell = NULL;
192dfc14de9SMatthew G. Knepley   const PetscSFNode *LA_sfcell;
193480eef7bSDave May   DM dmcell;
194480eef7bSDave May   Vec pos;
19540c453e9SDave May   PetscBool error_check = swarm->migrate_error_on_missing_point;
1967c6d1d28SDave May   PetscMPIInt commsize,rank;
197480eef7bSDave May 
198521f74f9SMatthew G. Knepley   PetscFunctionBegin;
199480eef7bSDave May   ierr = DMSwarmGetCellDM(dm,&dmcell);CHKERRQ(ierr);
200480eef7bSDave May   if (!dmcell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid if cell DM provided");
201480eef7bSDave May 
2027c6d1d28SDave May   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&commsize);CHKERRQ(ierr);
2037c6d1d28SDave May   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
2047c6d1d28SDave May 
205fb1bcc12SMatthew G. Knepley   ierr = DMSwarmCreateLocalVectorFromField(dm, DMSwarmPICField_coor, &pos);CHKERRQ(ierr);
206dfc14de9SMatthew G. Knepley   ierr = DMLocatePoints(dmcell, pos, DM_POINTLOCATION_NONE, &sfcell);CHKERRQ(ierr);
207fb1bcc12SMatthew G. Knepley   ierr = DMSwarmDestroyLocalVectorFromField(dm, DMSwarmPICField_coor, &pos);CHKERRQ(ierr);
208480eef7bSDave May 
20940c453e9SDave May   if (error_check) {
21040c453e9SDave May     ierr = DMSwarmGetSize(dm,&npointsg);CHKERRQ(ierr);
21140c453e9SDave May   }
212480eef7bSDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
21308056efcSDave May   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
214dfc14de9SMatthew G. Knepley   ierr = PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell);CHKERRQ(ierr);
215480eef7bSDave May   for (p=0; p<npoints; p++) {
216dfc14de9SMatthew G. Knepley     rankval[p] = LA_sfcell[p].index;
217480eef7bSDave May   }
21808056efcSDave May   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
219dfc14de9SMatthew G. Knepley   ierr = PetscSFDestroy(&sfcell);CHKERRQ(ierr);
220480eef7bSDave May 
2216fbf25f8SDave May   if (commsize > 1) {
222889dbfe5SDave May     ierr = DMSwarmMigrate_DMNeighborScatter(dm,dmcell,remove_sent_points,&npoints_prior_migration);CHKERRQ(ierr);
2236fbf25f8SDave May   } else {
2240ed23c7fSDave May     DataField PField;
2250ed23c7fSDave May     PetscInt npoints_curr;
2260ed23c7fSDave May 
2270ed23c7fSDave May     /* remove points which the domain */
2280ed23c7fSDave May     ierr = DataBucketGetDataFieldByName(swarm->db,DMSwarmField_rank,&PField);CHKERRQ(ierr);
2290ed23c7fSDave May     ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr);
2300ed23c7fSDave May 
2310ed23c7fSDave May     ierr = DataBucketGetSizes(swarm->db,&npoints_curr,NULL,NULL);CHKERRQ(ierr);
2320ed23c7fSDave May     for (p=0; p<npoints_curr; p++) {
2330ed23c7fSDave May       if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
2340ed23c7fSDave May         /* kill point */
2350ed23c7fSDave May         ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr);
2360ed23c7fSDave May         ierr = DataBucketGetSizes(swarm->db,&npoints_curr,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */
2370ed23c7fSDave May         ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); /* update date point increase realloc performed */
2380ed23c7fSDave May         p--; /* check replacement point */
2390ed23c7fSDave May       }
2400ed23c7fSDave May     }
2416fbf25f8SDave May     ierr = DMSwarmGetSize(dm,&npoints_prior_migration);CHKERRQ(ierr);
2420ed23c7fSDave May 
2436fbf25f8SDave May   }
244480eef7bSDave May 
24540c453e9SDave May   /* locate points newly recevied */
24640c453e9SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr);
247009b43efSDave May 
2487c6d1d28SDave May #if 0
249009b43efSDave May   { // safe alternative - however this performs two point locations on: (i) the intial points set and; (ii) the (intial + recieved) point set
2507c6d1d28SDave May     PetscScalar *LA_coor;
2517c6d1d28SDave May     PetscInt bs;
2527c6d1d28SDave May     DataField PField;
2537c6d1d28SDave May 
2547c6d1d28SDave May     ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr);
2557c6d1d28SDave May     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,bs*npoints2,(const PetscScalar*)LA_coor,&pos);CHKERRQ(ierr);
256dfc14de9SMatthew G. Knepley     ierr = DMLocatePoints(dmcell,pos,DM_POINTLOCATION_NONE,&sfcell);CHKERRQ(ierr);
2577c6d1d28SDave May 
2587c6d1d28SDave May     ierr = VecDestroy(&pos);CHKERRQ(ierr);
2597c6d1d28SDave May     ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr);
2607c6d1d28SDave May 
261dfc14de9SMatthew G. Knepley     ierr = PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell);CHKERRQ(ierr);
2627c6d1d28SDave May     ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
2637c6d1d28SDave May     for (p=0; p<npoints2; p++) {
264dfc14de9SMatthew G. Knepley       rankval[p] = LA_sfcell[p].index;
2657c6d1d28SDave May     }
266dfc14de9SMatthew G. Knepley     ierr = PetscSFDestroy(&sfcell);CHKERRQ(ierr);
26708056efcSDave May     ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
26840c453e9SDave May 
2697c6d1d28SDave May     /* remove points which left processor */
2707c6d1d28SDave May     ierr = DataBucketGetDataFieldByName(swarm->db,DMSwarmField_rank,&PField);CHKERRQ(ierr);
2717c6d1d28SDave May     ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr);
2727c6d1d28SDave May 
2737c6d1d28SDave May     ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr);
2747c6d1d28SDave May     for (p=0; p<npoints2; p++) {
275f954cb40SDave May       if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
2767c6d1d28SDave May         /* kill point */
2777c6d1d28SDave May         ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr);
2787c6d1d28SDave May         ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */
2797c6d1d28SDave May         ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); /* update date point increase realloc performed */
2807c6d1d28SDave May         p--; /* check replacement point */
2817c6d1d28SDave May       }
2827c6d1d28SDave May     }
28340c453e9SDave May   }
284009b43efSDave May #endif
285009b43efSDave May 
286009b43efSDave May   { // this performs two point locations: (i) on the intial points set prior to communication; and (ii) on the new (recieved) points
287009b43efSDave May     PetscScalar *LA_coor;
288009b43efSDave May     PetscInt npoints_from_neighbours,bs;
289009b43efSDave May     DataField PField;
290009b43efSDave May 
291009b43efSDave May     npoints_from_neighbours = npoints2 - npoints_prior_migration;
292009b43efSDave May 
293009b43efSDave May     ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr);
294009b43efSDave May     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,bs*npoints_from_neighbours,(const PetscScalar*)&LA_coor[bs*npoints_prior_migration],&pos);CHKERRQ(ierr);
295009b43efSDave May 
296009b43efSDave May     ierr = DMLocatePoints(dmcell,pos,DM_POINTLOCATION_NONE,&sfcell);CHKERRQ(ierr);
297009b43efSDave May 
298009b43efSDave May     ierr = VecDestroy(&pos);CHKERRQ(ierr);
299009b43efSDave May     ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr);
300009b43efSDave May 
301009b43efSDave May     ierr = PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell);CHKERRQ(ierr);
302009b43efSDave May     ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
303009b43efSDave May     for (p=0; p<npoints_from_neighbours; p++) {
304009b43efSDave May       rankval[npoints_prior_migration + p] = LA_sfcell[p].index;
305009b43efSDave May     }
306009b43efSDave May     ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
307009b43efSDave May     ierr = PetscSFDestroy(&sfcell);CHKERRQ(ierr);
308009b43efSDave May 
309009b43efSDave May     /* remove points which left processor */
310009b43efSDave May     ierr = DataBucketGetDataFieldByName(swarm->db,DMSwarmField_rank,&PField);CHKERRQ(ierr);
311009b43efSDave May     ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr);
312009b43efSDave May 
313009b43efSDave May     ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr);
314009b43efSDave May     for (p=npoints_prior_migration; p<npoints2; p++) {
315009b43efSDave May       if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
316009b43efSDave May         /* kill point */
317009b43efSDave May         ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr);
318009b43efSDave May         ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */
319009b43efSDave May         ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); /* update date point increase realloc performed */
320009b43efSDave May         p--; /* check replacement point */
321009b43efSDave May       }
322009b43efSDave May     }
323009b43efSDave May   }
324009b43efSDave May 
325*3151f1d5SDave May   {
326*3151f1d5SDave May     PetscInt *p_cellid;
327*3151f1d5SDave May 
328*3151f1d5SDave May     ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr);
329*3151f1d5SDave May     ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
330*3151f1d5SDave May     ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&p_cellid);CHKERRQ(ierr);
331*3151f1d5SDave May     for (p=0; p<npoints2; p++) {
332*3151f1d5SDave May       p_cellid[p] = rankval[p];
333*3151f1d5SDave May     }
334*3151f1d5SDave May     ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&p_cellid);CHKERRQ(ierr);
335*3151f1d5SDave May     ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
336*3151f1d5SDave May   }
337*3151f1d5SDave May 
33840c453e9SDave May   /* check for error on removed points */
33940c453e9SDave May   if (error_check) {
34040c453e9SDave May     ierr = DMSwarmGetSize(dm,&npoints2g);CHKERRQ(ierr);
34140c453e9SDave May     if (npointsg != npoints2g) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Points from the DMSwarm must remain constant during migration (initial %D - final %D)",npointsg,npoints2g);
34240c453e9SDave May   }
343480eef7bSDave May   PetscFunctionReturn(0);
344480eef7bSDave May }
345480eef7bSDave May 
34608056efcSDave May PetscErrorCode DMSwarmMigrate_CellDMExact(DM dm,PetscBool remove_sent_points)
34708056efcSDave May {
348521f74f9SMatthew G. Knepley   PetscFunctionBegin;
34908056efcSDave May   PetscFunctionReturn(0);
35008056efcSDave May }
35108056efcSDave May 
352480eef7bSDave May /*
353480eef7bSDave May  Redundant as this assumes points can only be sent to a single rank
354480eef7bSDave May */
3552712d1f2SDave May PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm,PetscInt *globalsize)
3562712d1f2SDave May {
3572712d1f2SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
3582712d1f2SDave May   PetscErrorCode ierr;
3592712d1f2SDave May   DataEx de;
3602712d1f2SDave May   PetscInt p,npoints,*rankval,n_points_recv;
3612712d1f2SDave May   PetscMPIInt rank,nrank,negrank;
3622712d1f2SDave May   void *point_buffer,*recv_points;
3632712d1f2SDave May   size_t sizeof_dmswarm_point;
3642712d1f2SDave May 
365521f74f9SMatthew G. Knepley   PetscFunctionBegin;
3662712d1f2SDave May   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
3672712d1f2SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
3682712d1f2SDave May   *globalsize = npoints;
36908056efcSDave May   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
370521f74f9SMatthew G. Knepley   ierr = DataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr);
3712712d1f2SDave May   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
3722712d1f2SDave May   for (p=0; p<npoints; p++) {
3732712d1f2SDave May     negrank = rankval[p];
3742712d1f2SDave May     if (negrank < 0) {
3752712d1f2SDave May       nrank = -negrank - 1;
3762712d1f2SDave May       ierr = DataExTopologyAddNeighbour(de,nrank);CHKERRQ(ierr);
3772712d1f2SDave May     }
3782712d1f2SDave May   }
3792712d1f2SDave May   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
3802712d1f2SDave May   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
3812712d1f2SDave May   for (p=0; p<npoints; p++) {
3822712d1f2SDave May     negrank = rankval[p];
3832712d1f2SDave May     if (negrank < 0) {
3842712d1f2SDave May       nrank = -negrank - 1;
3852712d1f2SDave May       ierr = DataExAddToSendCount(de,nrank,1);CHKERRQ(ierr);
3862712d1f2SDave May     }
3872712d1f2SDave May   }
3882712d1f2SDave May   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
3892712d1f2SDave May   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
3902712d1f2SDave May   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
3912712d1f2SDave May   for (p=0; p<npoints; p++) {
3922712d1f2SDave May     negrank = rankval[p];
3932712d1f2SDave May     if (negrank < 0) {
3942712d1f2SDave May       nrank = -negrank - 1;
3952712d1f2SDave May       rankval[p] = nrank;
3962712d1f2SDave May       /* copy point into buffer */
3972712d1f2SDave May       ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr);
3982712d1f2SDave May       /* insert point buffer into DataExchanger */
3992712d1f2SDave May       ierr = DataExPackData(de,nrank,1,point_buffer);CHKERRQ(ierr);
4002712d1f2SDave May       rankval[p] = negrank;
4012712d1f2SDave May     }
4022712d1f2SDave May   }
4032712d1f2SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
40408056efcSDave May   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
4052712d1f2SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
4062712d1f2SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
4072712d1f2SDave May   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
4082712d1f2SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
40965141ba8SDave May   ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr);
4102712d1f2SDave May   for (p=0; p<n_points_recv; p++) {
4112712d1f2SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
4122712d1f2SDave May 
4132712d1f2SDave May     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
4142712d1f2SDave May   }
4152712d1f2SDave May   ierr = DataExView(de);CHKERRQ(ierr);
4162712d1f2SDave May   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
4172712d1f2SDave May   ierr = DataExDestroy(de);CHKERRQ(ierr);
4182712d1f2SDave May   PetscFunctionReturn(0);
4192712d1f2SDave May }
420b16650c8SDave May 
421b16650c8SDave May typedef struct {
422b16650c8SDave May   PetscMPIInt owner_rank;
423b16650c8SDave May   PetscReal min[3],max[3];
424b16650c8SDave May } CollectBBox;
425b16650c8SDave May 
426fe39f135SDave May PETSC_EXTERN PetscErrorCode DMSwarmCollect_DMDABoundingBox(DM dm,PetscInt *globalsize)
427b16650c8SDave May {
428b16650c8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
429b16650c8SDave May   PetscErrorCode ierr;
430b16650c8SDave May   DataEx de;
431b16650c8SDave May   PetscInt p,pk,npoints,*rankval,n_points_recv,n_bbox_recv,dim,neighbour_cells;
432b16650c8SDave May   PetscMPIInt rank,nrank;
433b16650c8SDave May   void *point_buffer,*recv_points;
434b16650c8SDave May   size_t sizeof_dmswarm_point,sizeof_bbox_ctx;
435b16650c8SDave May   PetscBool isdmda;
436b16650c8SDave May   CollectBBox *bbox,*recv_bbox;
437b16650c8SDave May   const PetscMPIInt *dmneighborranks;
438b16650c8SDave May   DM dmcell;
439b16650c8SDave May 
440521f74f9SMatthew G. Knepley   PetscFunctionBegin;
441b16650c8SDave May   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
442b16650c8SDave May 
443fe39f135SDave May   ierr = DMSwarmGetCellDM(dm,&dmcell);CHKERRQ(ierr);
444fe39f135SDave May   if (!dmcell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid if cell DM provided");
445b16650c8SDave May   isdmda = PETSC_FALSE;
446b16650c8SDave May   PetscObjectTypeCompare((PetscObject)dmcell,DMDA,&isdmda);
447b16650c8SDave May   if (!isdmda) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only DMDA support for CollectBoundingBox");
448b16650c8SDave May 
4498dbd68bcSDave May   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
450b16650c8SDave May   sizeof_bbox_ctx = sizeof(CollectBBox);
451b16650c8SDave May   PetscMalloc1(1,&bbox);
452b16650c8SDave May   bbox->owner_rank = rank;
453b16650c8SDave May 
454b16650c8SDave May   /* compute the bounding box based on the overlapping / stenctil size */
455b16650c8SDave May   {
456b16650c8SDave May     Vec lcoor;
457b16650c8SDave May 
458b16650c8SDave May     ierr = DMGetCoordinatesLocal(dmcell,&lcoor);CHKERRQ(ierr);
459fe39f135SDave May     if (dim >= 1) {
460b16650c8SDave May       ierr = VecStrideMin(lcoor,0,NULL,&bbox->min[0]);CHKERRQ(ierr);
461b16650c8SDave May       ierr = VecStrideMax(lcoor,0,NULL,&bbox->max[0]);CHKERRQ(ierr);
462fe39f135SDave May     }
463fe39f135SDave May     if (dim >= 2) {
464fe39f135SDave May       ierr = VecStrideMin(lcoor,1,NULL,&bbox->min[1]);CHKERRQ(ierr);
465b16650c8SDave May       ierr = VecStrideMax(lcoor,1,NULL,&bbox->max[1]);CHKERRQ(ierr);
466b16650c8SDave May     }
467fe39f135SDave May     if (dim == 3) {
468fe39f135SDave May       ierr = VecStrideMin(lcoor,2,NULL,&bbox->min[2]);CHKERRQ(ierr);
469fe39f135SDave May       ierr = VecStrideMax(lcoor,2,NULL,&bbox->max[2]);CHKERRQ(ierr);
470fe39f135SDave May     }
471fe39f135SDave May   }
472b16650c8SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
473b16650c8SDave May   *globalsize = npoints;
47408056efcSDave May   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
475521f74f9SMatthew G. Knepley   ierr = DataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr);
476b16650c8SDave May   /* use DMDA neighbours */
477b16650c8SDave May   ierr = DMDAGetNeighbors(dmcell,&dmneighborranks);CHKERRQ(ierr);
4788dbd68bcSDave May   if (dim == 1) {
4798dbd68bcSDave May     neighbour_cells = 3;
4808dbd68bcSDave May   } else if (dim == 2) {
4818dbd68bcSDave May     neighbour_cells = 9;
4828dbd68bcSDave May   } else {
4838dbd68bcSDave May     neighbour_cells = 27;
4848dbd68bcSDave May   }
485b16650c8SDave May   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
486b16650c8SDave May   for (p=0; p<neighbour_cells; p++) {
487b16650c8SDave May     if ( (dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank) ) {
488b16650c8SDave May       ierr = DataExTopologyAddNeighbour(de,dmneighborranks[p]);CHKERRQ(ierr);
489b16650c8SDave May     }
490b16650c8SDave May   }
491b16650c8SDave May   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
492b16650c8SDave May   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
493b16650c8SDave May   for (p=0; p<neighbour_cells; p++) {
494b16650c8SDave May     if ( (dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank) ) {
495b16650c8SDave May       ierr = DataExAddToSendCount(de,dmneighborranks[p],1);CHKERRQ(ierr);
496b16650c8SDave May     }
497b16650c8SDave May   }
498b16650c8SDave May   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
499b16650c8SDave May   /* send bounding boxes */
500b16650c8SDave May   ierr = DataExPackInitialize(de,sizeof_bbox_ctx);CHKERRQ(ierr);
501b16650c8SDave May   for (p=0; p<neighbour_cells; p++) {
502b16650c8SDave May     nrank = dmneighborranks[p];
503b16650c8SDave May     if ( (nrank >= 0) && (nrank != rank) ) {
504b16650c8SDave May       /* insert bbox buffer into DataExchanger */
505b16650c8SDave May       ierr = DataExPackData(de,nrank,1,bbox);CHKERRQ(ierr);
506b16650c8SDave May     }
507b16650c8SDave May   }
508b16650c8SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
509b16650c8SDave May   /* recv bounding boxes */
510b16650c8SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
511b16650c8SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
512b16650c8SDave May   ierr = DataExGetRecvData(de,&n_bbox_recv,(void**)&recv_bbox);CHKERRQ(ierr);
513b16650c8SDave May   for (p=0; p<n_bbox_recv; p++) {
514b4b02483SDave May     PetscPrintf(PETSC_COMM_SELF,"[rank %d]: box from %d : range[%+1.4e,%+1.4e]x[%+1.4e,%+1.4e]\n",rank,recv_bbox[p].owner_rank,
515b4b02483SDave May            (double)recv_bbox[p].min[0],(double)recv_bbox[p].max[0],(double)recv_bbox[p].min[1],(double)recv_bbox[p].max[1]);
516b16650c8SDave May   }
517b16650c8SDave May   /* of course this is stupid as this "generic" function should have a better way to know what the coordinates are called */
518b16650c8SDave May   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
519b16650c8SDave May   for (pk=0; pk<n_bbox_recv; pk++) {
520b16650c8SDave May     PetscReal *array_x,*array_y;
521b16650c8SDave May 
522b16650c8SDave May     ierr = DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
523b16650c8SDave May     ierr = DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
524b16650c8SDave May     for (p=0; p<npoints; p++) {
525b16650c8SDave May       if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0]) ) {
526b16650c8SDave May         if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1]) ) {
527b16650c8SDave May           ierr = DataExAddToSendCount(de,recv_bbox[pk].owner_rank,1);CHKERRQ(ierr);
528b16650c8SDave May         }
529b16650c8SDave May       }
530b16650c8SDave May     }
531b16650c8SDave May     ierr = DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
532b16650c8SDave May     ierr = DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
533b16650c8SDave May   }
534b16650c8SDave May   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
535b16650c8SDave May   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
536b16650c8SDave May   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
537b16650c8SDave May   for (pk=0; pk<n_bbox_recv; pk++) {
538b16650c8SDave May     PetscReal *array_x,*array_y;
539b16650c8SDave May 
540b16650c8SDave May     ierr = DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
541b16650c8SDave May     ierr = DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
542b16650c8SDave May     for (p=0; p<npoints; p++) {
543b16650c8SDave May       if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0]) ) {
544b16650c8SDave May         if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1]) ) {
545521f74f9SMatthew G. Knepley           /* copy point into buffer */
546b16650c8SDave May           ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr);
547521f74f9SMatthew G. Knepley           /* insert point buffer into DataExchanger */
548b16650c8SDave May           ierr = DataExPackData(de,recv_bbox[pk].owner_rank,1,point_buffer);CHKERRQ(ierr);
549b16650c8SDave May         }
550b16650c8SDave May       }
551b16650c8SDave May     }
552b16650c8SDave May     ierr = DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
553b16650c8SDave May     ierr = DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
554b16650c8SDave May   }
555b16650c8SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
55608056efcSDave May   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
557b16650c8SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
558b16650c8SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
559b16650c8SDave May   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
560b16650c8SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
56165141ba8SDave May   ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr);
562b16650c8SDave May   for (p=0; p<n_points_recv; p++) {
563b16650c8SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
564b16650c8SDave May 
565b16650c8SDave May     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
566b16650c8SDave May   }
567b16650c8SDave May   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
568b16650c8SDave May   PetscFree(bbox);
569b16650c8SDave May   ierr = DataExView(de);CHKERRQ(ierr);
570b16650c8SDave May   ierr = DataExDestroy(de);CHKERRQ(ierr);
571b16650c8SDave May   PetscFunctionReturn(0);
572b16650c8SDave May }
573a9fd7477SDave 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 */
584a9fd7477SDave May PETSC_EXTERN PetscErrorCode DMSwarmCollect_General(DM dm,PetscErrorCode (*collect)(DM,void*,PetscInt*,PetscInt**),size_t ctx_size,void *ctx,PetscInt *globalsize)
585a9fd7477SDave May {
586a9fd7477SDave May   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
587a9fd7477SDave May   PetscErrorCode ierr;
588a9fd7477SDave May   DataEx         de;
589a9fd7477SDave May   PetscInt       p,r,npoints,n_points_recv;
590a9fd7477SDave May   PetscMPIInt    commsize,rank;
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;
597a9fd7477SDave May   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&commsize);CHKERRQ(ierr);
598a9fd7477SDave May   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
599a9fd7477SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
600a9fd7477SDave May   *globalsize = npoints;
601a9fd7477SDave May   /* Broadcast user context */
602a9fd7477SDave May   PetscMalloc(ctx_size*commsize,&ctxlist);
603a9fd7477SDave May   ierr = MPI_Allgather(ctx,ctx_size,MPI_CHAR,ctxlist,ctx_size,MPI_CHAR,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
604521f74f9SMatthew G. Knepley   ierr = PetscMalloc1(commsize,&n2collect);CHKERRQ(ierr);
605521f74f9SMatthew G. Knepley   ierr = PetscMalloc1(commsize,&collectlist);CHKERRQ(ierr);
606a9fd7477SDave May   for (r=0; r<commsize; 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 );
615a9fd7477SDave May       ierr = collect(dm,_ctx_r,&_n2collect,&_collectlist);CHKERRQ(ierr);
616a9fd7477SDave May     }
617a9fd7477SDave May     n2collect[r]   = _n2collect;
618a9fd7477SDave May     collectlist[r] = _collectlist;
619a9fd7477SDave May   }
620521f74f9SMatthew G. Knepley   ierr = DataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr);
621a9fd7477SDave May   /* Define topology */
622a9fd7477SDave May   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
623a9fd7477SDave May   for (r=0; r<commsize; r++) {
624a9fd7477SDave May     if (n2collect[r] > 0) {
625a9fd7477SDave May       ierr = DataExTopologyAddNeighbour(de,(PetscMPIInt)r);CHKERRQ(ierr);
626a9fd7477SDave May     }
627a9fd7477SDave May   }
628a9fd7477SDave May   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
629a9fd7477SDave May   /* Define send counts */
630a9fd7477SDave May   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
631a9fd7477SDave May   for (r=0; r<commsize; r++) {
632a9fd7477SDave May     if (n2collect[r] > 0) {
633a9fd7477SDave May       ierr = DataExAddToSendCount(de,r,n2collect[r]);CHKERRQ(ierr);
634a9fd7477SDave May     }
635a9fd7477SDave May   }
636a9fd7477SDave May   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
637a9fd7477SDave May   /* Pack data */
638a9fd7477SDave May   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
639a9fd7477SDave May   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
640a9fd7477SDave May   for (r=0; r<commsize; r++) {
641a9fd7477SDave May     for (p=0; p<n2collect[r]; p++) {
642a9fd7477SDave May       ierr = DataBucketFillPackedArray(swarm->db,collectlist[r][p],point_buffer);CHKERRQ(ierr);
643a9fd7477SDave May       /* insert point buffer into the data exchanger */
644a9fd7477SDave May       ierr = DataExPackData(de,r,1,point_buffer);CHKERRQ(ierr);
645a9fd7477SDave May     }
646a9fd7477SDave May   }
647a9fd7477SDave May   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
648a9fd7477SDave May   /* Scatter */
649a9fd7477SDave May   ierr = DataExBegin(de);CHKERRQ(ierr);
650a9fd7477SDave May   ierr = DataExEnd(de);CHKERRQ(ierr);
651a9fd7477SDave May   /* Collect data in DMSwarm container */
652a9fd7477SDave May   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
653a9fd7477SDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
65465141ba8SDave May   ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr);
655a9fd7477SDave May   for (p=0; p<n_points_recv; p++) {
656a9fd7477SDave May     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
657a9fd7477SDave May 
658a9fd7477SDave May     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
659a9fd7477SDave May   }
660a9fd7477SDave May   /* Release memory */
661a9fd7477SDave May   for (r=0; r<commsize; r++) {
662a9fd7477SDave May     if (collectlist[r]) PetscFree(collectlist[r]);
663a9fd7477SDave May   }
664521f74f9SMatthew G. Knepley   ierr = PetscFree(collectlist);CHKERRQ(ierr);
665521f74f9SMatthew G. Knepley   ierr = PetscFree(n2collect);CHKERRQ(ierr);
666521f74f9SMatthew G. Knepley   ierr = PetscFree(ctxlist);CHKERRQ(ierr);
667a9fd7477SDave May   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
668a9fd7477SDave May   ierr = DataExView(de);CHKERRQ(ierr);
669a9fd7477SDave May   ierr = DataExDestroy(de);CHKERRQ(ierr);
670a9fd7477SDave May   PetscFunctionReturn(0);
671a9fd7477SDave May }
672a9fd7477SDave May 
673