xref: /petsc/src/dm/impls/swarm/swarm_migrate.c (revision 609bdbee21ea3be08735c64dbe00a9ab27759925)
1 
2 #include <petscsf.h>
3 #include <petscdmswarm.h>
4 #include <petsc/private/dmswarmimpl.h>    /*I   "petscdmswarm.h"   I*/
5 #include "data_bucket.h"
6 #include "data_ex.h"
7 
8 
9 /*
10  User loads desired location (MPI rank) into field DMSwarm_rank
11 */
12 PetscErrorCode DMSwarmMigrate_Push_Basic(DM dm,PetscBool remove_sent_points)
13 {
14   DM_Swarm *swarm = (DM_Swarm*)dm->data;
15   PetscErrorCode ierr;
16   DataEx de;
17   PetscInt p,npoints,*rankval,n_points_recv;
18   PetscMPIInt rank,nrank;
19   void *point_buffer,*recv_points;
20   size_t sizeof_dmswarm_point;
21 
22   PetscFunctionBegin;
23   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
24 
25   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
26   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
27   ierr = DataExCreate(PetscObjectComm((PetscObject)dm),0, &de);CHKERRQ(ierr);
28   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
29   for (p = 0; p < npoints; ++p) {
30     nrank = rankval[p];
31     if (nrank != rank) {
32       ierr = DataExTopologyAddNeighbour(de,nrank);CHKERRQ(ierr);
33     }
34   }
35   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
36   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
37   for (p=0; p<npoints; p++) {
38     nrank = rankval[p];
39     if (nrank != rank) {
40       ierr = DataExAddToSendCount(de,nrank,1);CHKERRQ(ierr);
41     }
42   }
43   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
44   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
45   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
46   for (p=0; p<npoints; p++) {
47     nrank = rankval[p];
48     if (nrank != rank) {
49       /* copy point into buffer */
50       ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr);
51       /* insert point buffer into DataExchanger */
52       ierr = DataExPackData(de,nrank,1,point_buffer);CHKERRQ(ierr);
53     }
54   }
55   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
56   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
57 
58   if (remove_sent_points) {
59     DataField gfield;
60 
61     ierr = DataBucketGetDataFieldByName(swarm->db,DMSwarmField_rank,&gfield);CHKERRQ(ierr);
62     ierr = DataFieldGetAccess(gfield);CHKERRQ(ierr);
63     ierr = DataFieldGetEntries(gfield,(void**)&rankval);CHKERRQ(ierr);
64 
65     /* remove points which left processor */
66     ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
67     for (p=0; p<npoints; p++) {
68       nrank = rankval[p];
69       if (nrank != rank) {
70         /* kill point */
71         ierr = DataFieldRestoreAccess(gfield);CHKERRQ(ierr);
72 
73         ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr);
74 
75         ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */
76         ierr = DataFieldGetAccess(gfield);CHKERRQ(ierr);
77         ierr = DataFieldGetEntries(gfield,(void**)&rankval);CHKERRQ(ierr);
78         p--; /* check replacement point */
79       }
80     }
81     ierr = DataFieldRestoreEntries(gfield,(void**)&rankval);CHKERRQ(ierr);
82     ierr = DataFieldRestoreAccess(gfield);CHKERRQ(ierr);
83   }
84   ierr = DataExBegin(de);CHKERRQ(ierr);
85   ierr = DataExEnd(de);CHKERRQ(ierr);
86   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
87   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
88   ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr);
89   for (p=0; p<n_points_recv; p++) {
90     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
91 
92     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
93   }
94   ierr = DataExView(de);CHKERRQ(ierr);
95   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
96   ierr = DataExDestroy(de);CHKERRQ(ierr);
97   PetscFunctionReturn(0);
98 }
99 
100 PetscErrorCode DMSwarmMigrate_DMNeighborScatter(DM dm,DM dmcell,PetscBool remove_sent_points,PetscInt *npoints_prior_migration)
101 {
102   DM_Swarm *swarm = (DM_Swarm*)dm->data;
103   PetscErrorCode ierr;
104   DataEx de;
105   PetscInt r,p,npoints,*rankval,n_points_recv;
106   PetscMPIInt rank,_rank;
107   const PetscMPIInt *neighbourranks;
108   void *point_buffer,*recv_points;
109   size_t sizeof_dmswarm_point;
110   PetscInt nneighbors;
111   PetscMPIInt mynneigh,*myneigh;
112 
113   PetscFunctionBegin;
114   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
115   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
116   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
117   ierr = DataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr);
118   ierr = DMGetNeighbors(dmcell,&nneighbors,&neighbourranks);CHKERRQ(ierr);
119   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
120   for (r=0; r<nneighbors; r++) {
121     _rank = neighbourranks[r];
122     if ((_rank != rank) && (_rank > 0)) {
123       ierr = DataExTopologyAddNeighbour(de,_rank);CHKERRQ(ierr);
124     }
125   }
126   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
127   ierr = DataExTopologyGetNeighbours(de,&mynneigh,&myneigh);CHKERRQ(ierr);
128   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
129   for (p=0; p<npoints; p++) {
130     if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
131       for (r=0; r<mynneigh; r++) {
132         _rank = myneigh[r];
133         ierr = DataExAddToSendCount(de,_rank,1);CHKERRQ(ierr);
134       }
135     }
136   }
137   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
138   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
139   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
140   for (p=0; p<npoints; p++) {
141     if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
142       for (r=0; r<mynneigh; r++) {
143         _rank = myneigh[r];
144         /* copy point into buffer */
145         ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr);
146         /* insert point buffer into DataExchanger */
147         ierr = DataExPackData(de,_rank,1,point_buffer);CHKERRQ(ierr);
148       }
149     }
150   }
151   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
152   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
153   if (remove_sent_points) {
154     DataField PField;
155 
156     ierr = DataBucketGetDataFieldByName(swarm->db,DMSwarmField_rank,&PField);CHKERRQ(ierr);
157     ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr);
158     /* remove points which left processor */
159     ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
160     for (p=0; p<npoints; p++) {
161       if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
162         /* kill point */
163         ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr);
164         ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */
165         ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); /* update date point increase realloc performed */
166         p--; /* check replacement point */
167       }
168     }
169   }
170   ierr = DataBucketGetSizes(swarm->db,npoints_prior_migration,NULL,NULL);CHKERRQ(ierr);
171   ierr = DataExBegin(de);CHKERRQ(ierr);
172   ierr = DataExEnd(de);CHKERRQ(ierr);
173   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
174   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
175   ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr);
176   for (p=0; p<n_points_recv; p++) {
177     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
178 
179     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
180   }
181   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
182   ierr = DataExDestroy(de);CHKERRQ(ierr);
183   PetscFunctionReturn(0);
184 }
185 
186 PetscErrorCode DMSwarmMigrate_CellDMScatter(DM dm,PetscBool remove_sent_points)
187 {
188   DM_Swarm *swarm = (DM_Swarm*)dm->data;
189   PetscErrorCode ierr;
190   PetscInt p,npoints,npointsg=0,npoints2,npoints2g,*rankval,npoints_prior_migration;
191   PetscSF sfcell = NULL;
192   const PetscSFNode *LA_sfcell;
193   DM dmcell;
194   Vec pos;
195   PetscBool error_check = swarm->migrate_error_on_missing_point;
196   PetscMPIInt commsize,rank;
197 
198   PetscFunctionBegin;
199   ierr = DMSwarmGetCellDM(dm,&dmcell);CHKERRQ(ierr);
200   if (!dmcell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid if cell DM provided");
201 
202   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&commsize);CHKERRQ(ierr);
203   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
204 
205   ierr = DMSwarmCreateLocalVectorFromField(dm, DMSwarmPICField_coor, &pos);CHKERRQ(ierr);
206   ierr = DMLocatePoints(dmcell, pos, DM_POINTLOCATION_NONE, &sfcell);CHKERRQ(ierr);
207   ierr = DMSwarmDestroyLocalVectorFromField(dm, DMSwarmPICField_coor, &pos);CHKERRQ(ierr);
208 
209   if (error_check) {
210     ierr = DMSwarmGetSize(dm,&npointsg);CHKERRQ(ierr);
211   }
212   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
213   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
214   ierr = PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell);CHKERRQ(ierr);
215   for (p=0; p<npoints; p++) {
216     rankval[p] = LA_sfcell[p].index;
217   }
218   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
219   ierr = PetscSFDestroy(&sfcell);CHKERRQ(ierr);
220 
221   if (commsize > 1) {
222     ierr = DMSwarmMigrate_DMNeighborScatter(dm,dmcell,remove_sent_points,&npoints_prior_migration);CHKERRQ(ierr);
223   } else {
224     ierr = DMSwarmGetSize(dm,&npoints_prior_migration);CHKERRQ(ierr);
225   }
226 
227   /* locate points newly recevied */
228   ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr);
229 #if 0
230   len = npoints2 - npoints_prior_migration;
231   if (len > 0) {
232     PetscScalar *LA_coor;
233     PetscInt bs;
234 
235     ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr);
236     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,bs*len,(const PetscScalar*)&LA_coor[bs*npoints_prior_migration],&pos);CHKERRQ(ierr);
237 
238     ierr = DMLocatePoints(dmcell,pos,&iscell);CHKERRQ(ierr);
239 
240     ierr = VecDestroy(&pos);CHKERRQ(ierr);
241     ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr);
242 
243     ierr = ISGetIndices(iscell,&LA_iscell);CHKERRQ(ierr);
244     ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
245     for (p=0; p<len; p++) {
246       rankval[npoints_prior_migration+p] = LA_iscell[p];
247     }
248     ierr = ISRestoreIndices(iscell,&LA_iscell);CHKERRQ(ierr);
249     ierr = ISDestroy(&iscell);CHKERRQ(ierr);
250     ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
251 
252     /* remove points which left processor */
253     ierr = DataBucketGetDataFieldByName(swarm->db,DMSwarmField_rank,&PField);CHKERRQ(ierr);
254     ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr);
255 
256     ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
257     for (p=npoints_prior_migration; p<npoints2; p++) {
258       if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
259         /* kill point */
260         ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr);
261         ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */
262         ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); /* update date point increase realloc performed */
263         p--; /* check replacement point */
264       }
265     }
266 
267   }
268 #endif
269 
270   {
271     PetscScalar *LA_coor;
272     PetscInt bs;
273     DataField PField;
274 
275     ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr);
276     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,bs*npoints2,(const PetscScalar*)LA_coor,&pos);CHKERRQ(ierr);
277     ierr = DMLocatePoints(dmcell,pos,DM_POINTLOCATION_NONE,&sfcell);CHKERRQ(ierr);
278 
279     ierr = VecDestroy(&pos);CHKERRQ(ierr);
280     ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr);
281 
282     ierr = PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell);CHKERRQ(ierr);
283     ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
284     for (p=0; p<npoints2; p++) {
285       rankval[p] = LA_sfcell[p].index;
286     }
287     ierr = PetscSFDestroy(&sfcell);CHKERRQ(ierr);
288     ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
289 
290     /* remove points which left processor */
291     ierr = DataBucketGetDataFieldByName(swarm->db,DMSwarmField_rank,&PField);CHKERRQ(ierr);
292     ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr);
293 
294     ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr);
295     for (p=0; p<npoints2; p++) {
296       if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
297         /* kill point */
298         ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr);
299         ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */
300         ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); /* update date point increase realloc performed */
301         p--; /* check replacement point */
302       }
303     }
304   }
305   /* check for error on removed points */
306   if (error_check) {
307     ierr = DMSwarmGetSize(dm,&npoints2g);CHKERRQ(ierr);
308     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);
309   }
310   PetscFunctionReturn(0);
311 }
312 
313 PetscErrorCode DMSwarmMigrate_CellDMExact(DM dm,PetscBool remove_sent_points)
314 {
315   PetscFunctionBegin;
316   PetscFunctionReturn(0);
317 }
318 
319 /*
320  Redundant as this assumes points can only be sent to a single rank
321 */
322 PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm,PetscInt *globalsize)
323 {
324   DM_Swarm *swarm = (DM_Swarm*)dm->data;
325   PetscErrorCode ierr;
326   DataEx de;
327   PetscInt p,npoints,*rankval,n_points_recv;
328   PetscMPIInt rank,nrank,negrank;
329   void *point_buffer,*recv_points;
330   size_t sizeof_dmswarm_point;
331 
332   PetscFunctionBegin;
333   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
334   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
335   *globalsize = npoints;
336   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
337   ierr = DataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr);
338   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
339   for (p=0; p<npoints; p++) {
340     negrank = rankval[p];
341     if (negrank < 0) {
342       nrank = -negrank - 1;
343       ierr = DataExTopologyAddNeighbour(de,nrank);CHKERRQ(ierr);
344     }
345   }
346   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
347   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
348   for (p=0; p<npoints; p++) {
349     negrank = rankval[p];
350     if (negrank < 0) {
351       nrank = -negrank - 1;
352       ierr = DataExAddToSendCount(de,nrank,1);CHKERRQ(ierr);
353     }
354   }
355   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
356   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
357   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
358   for (p=0; p<npoints; p++) {
359     negrank = rankval[p];
360     if (negrank < 0) {
361       nrank = -negrank - 1;
362       rankval[p] = nrank;
363       /* copy point into buffer */
364       ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr);
365       /* insert point buffer into DataExchanger */
366       ierr = DataExPackData(de,nrank,1,point_buffer);CHKERRQ(ierr);
367       rankval[p] = negrank;
368     }
369   }
370   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
371   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
372   ierr = DataExBegin(de);CHKERRQ(ierr);
373   ierr = DataExEnd(de);CHKERRQ(ierr);
374   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
375   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
376   ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr);
377   for (p=0; p<n_points_recv; p++) {
378     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
379 
380     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
381   }
382   ierr = DataExView(de);CHKERRQ(ierr);
383   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
384   ierr = DataExDestroy(de);CHKERRQ(ierr);
385   PetscFunctionReturn(0);
386 }
387 
388 typedef struct {
389   PetscMPIInt owner_rank;
390   PetscReal min[3],max[3];
391 } CollectBBox;
392 
393 PETSC_EXTERN PetscErrorCode DMSwarmCollect_DMDABoundingBox(DM dm,PetscInt *globalsize)
394 {
395   DM_Swarm *swarm = (DM_Swarm*)dm->data;
396   PetscErrorCode ierr;
397   DataEx de;
398   PetscInt p,pk,npoints,*rankval,n_points_recv,n_bbox_recv,dim,neighbour_cells;
399   PetscMPIInt rank,nrank;
400   void *point_buffer,*recv_points;
401   size_t sizeof_dmswarm_point,sizeof_bbox_ctx;
402   PetscBool isdmda;
403   CollectBBox *bbox,*recv_bbox;
404   const PetscMPIInt *dmneighborranks;
405   DM dmcell;
406 
407   PetscFunctionBegin;
408   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
409 
410   ierr = DMSwarmGetCellDM(dm,&dmcell);CHKERRQ(ierr);
411   if (!dmcell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid if cell DM provided");
412   isdmda = PETSC_FALSE;
413   PetscObjectTypeCompare((PetscObject)dmcell,DMDA,&isdmda);
414   if (!isdmda) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only DMDA support for CollectBoundingBox");
415 
416   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
417   sizeof_bbox_ctx = sizeof(CollectBBox);
418   PetscMalloc1(1,&bbox);
419   bbox->owner_rank = rank;
420 
421   /* compute the bounding box based on the overlapping / stenctil size */
422   {
423     Vec lcoor;
424 
425     ierr = DMGetCoordinatesLocal(dmcell,&lcoor);CHKERRQ(ierr);
426     if (dim >= 1) {
427       ierr = VecStrideMin(lcoor,0,NULL,&bbox->min[0]);CHKERRQ(ierr);
428       ierr = VecStrideMax(lcoor,0,NULL,&bbox->max[0]);CHKERRQ(ierr);
429     }
430     if (dim >= 2) {
431       ierr = VecStrideMin(lcoor,1,NULL,&bbox->min[1]);CHKERRQ(ierr);
432       ierr = VecStrideMax(lcoor,1,NULL,&bbox->max[1]);CHKERRQ(ierr);
433     }
434     if (dim == 3) {
435       ierr = VecStrideMin(lcoor,2,NULL,&bbox->min[2]);CHKERRQ(ierr);
436       ierr = VecStrideMax(lcoor,2,NULL,&bbox->max[2]);CHKERRQ(ierr);
437     }
438   }
439   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
440   *globalsize = npoints;
441   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
442   ierr = DataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr);
443   /* use DMDA neighbours */
444   ierr = DMDAGetNeighbors(dmcell,&dmneighborranks);CHKERRQ(ierr);
445   if (dim == 1) {
446     neighbour_cells = 3;
447   } else if (dim == 2) {
448     neighbour_cells = 9;
449   } else {
450     neighbour_cells = 27;
451   }
452   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
453   for (p=0; p<neighbour_cells; p++) {
454     if ( (dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank) ) {
455       ierr = DataExTopologyAddNeighbour(de,dmneighborranks[p]);CHKERRQ(ierr);
456     }
457   }
458   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
459   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
460   for (p=0; p<neighbour_cells; p++) {
461     if ( (dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank) ) {
462       ierr = DataExAddToSendCount(de,dmneighborranks[p],1);CHKERRQ(ierr);
463     }
464   }
465   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
466   /* send bounding boxes */
467   ierr = DataExPackInitialize(de,sizeof_bbox_ctx);CHKERRQ(ierr);
468   for (p=0; p<neighbour_cells; p++) {
469     nrank = dmneighborranks[p];
470     if ( (nrank >= 0) && (nrank != rank) ) {
471       /* insert bbox buffer into DataExchanger */
472       ierr = DataExPackData(de,nrank,1,bbox);CHKERRQ(ierr);
473     }
474   }
475   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
476   /* recv bounding boxes */
477   ierr = DataExBegin(de);CHKERRQ(ierr);
478   ierr = DataExEnd(de);CHKERRQ(ierr);
479   ierr = DataExGetRecvData(de,&n_bbox_recv,(void**)&recv_bbox);CHKERRQ(ierr);
480   for (p=0; p<n_bbox_recv; p++) {
481     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,
482            (double)recv_bbox[p].min[0],(double)recv_bbox[p].max[0],(double)recv_bbox[p].min[1],(double)recv_bbox[p].max[1]);
483   }
484   /* of course this is stupid as this "generic" function should have a better way to know what the coordinates are called */
485   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
486   for (pk=0; pk<n_bbox_recv; pk++) {
487     PetscReal *array_x,*array_y;
488 
489     ierr = DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
490     ierr = DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
491     for (p=0; p<npoints; p++) {
492       if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0]) ) {
493         if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1]) ) {
494           ierr = DataExAddToSendCount(de,recv_bbox[pk].owner_rank,1);CHKERRQ(ierr);
495         }
496       }
497     }
498     ierr = DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
499     ierr = DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
500   }
501   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
502   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
503   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
504   for (pk=0; pk<n_bbox_recv; pk++) {
505     PetscReal *array_x,*array_y;
506 
507     ierr = DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
508     ierr = DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
509     for (p=0; p<npoints; p++) {
510       if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0]) ) {
511         if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1]) ) {
512           /* copy point into buffer */
513           ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr);
514           /* insert point buffer into DataExchanger */
515           ierr = DataExPackData(de,recv_bbox[pk].owner_rank,1,point_buffer);CHKERRQ(ierr);
516         }
517       }
518     }
519     ierr = DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
520     ierr = DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
521   }
522   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
523   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
524   ierr = DataExBegin(de);CHKERRQ(ierr);
525   ierr = DataExEnd(de);CHKERRQ(ierr);
526   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
527   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
528   ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr);
529   for (p=0; p<n_points_recv; p++) {
530     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
531 
532     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
533   }
534   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
535   PetscFree(bbox);
536   ierr = DataExView(de);CHKERRQ(ierr);
537   ierr = DataExDestroy(de);CHKERRQ(ierr);
538   PetscFunctionReturn(0);
539 }
540 
541 
542 /* General collection when no order, or neighbour information is provided */
543 /*
544  User provides context and collect() method
545  Broadcast user context
546 
547  for each context / rank {
548    collect(swarm,context,n,list)
549  }
550 */
551 PETSC_EXTERN PetscErrorCode DMSwarmCollect_General(DM dm,PetscErrorCode (*collect)(DM,void*,PetscInt*,PetscInt**),size_t ctx_size,void *ctx,PetscInt *globalsize)
552 {
553   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
554   PetscErrorCode ierr;
555   DataEx         de;
556   PetscInt       p,r,npoints,n_points_recv;
557   PetscMPIInt    commsize,rank;
558   void           *point_buffer,*recv_points;
559   void           *ctxlist;
560   PetscInt       *n2collect,**collectlist;
561   size_t         sizeof_dmswarm_point;
562 
563   PetscFunctionBegin;
564   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&commsize);CHKERRQ(ierr);
565   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
566   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
567   *globalsize = npoints;
568   /* Broadcast user context */
569   PetscMalloc(ctx_size*commsize,&ctxlist);
570   ierr = MPI_Allgather(ctx,ctx_size,MPI_CHAR,ctxlist,ctx_size,MPI_CHAR,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
571   ierr = PetscMalloc1(commsize,&n2collect);CHKERRQ(ierr);
572   ierr = PetscMalloc1(commsize,&collectlist);CHKERRQ(ierr);
573   for (r=0; r<commsize; r++) {
574     PetscInt _n2collect;
575     PetscInt *_collectlist;
576     void     *_ctx_r;
577 
578     _n2collect   = 0;
579     _collectlist = NULL;
580     if (r != rank) { /* don't collect data from yourself */
581       _ctx_r = (void*)( (char*)ctxlist + r * ctx_size );
582       ierr = collect(dm,_ctx_r,&_n2collect,&_collectlist);CHKERRQ(ierr);
583     }
584     n2collect[r]   = _n2collect;
585     collectlist[r] = _collectlist;
586   }
587   ierr = DataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr);
588   /* Define topology */
589   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
590   for (r=0; r<commsize; r++) {
591     if (n2collect[r] > 0) {
592       ierr = DataExTopologyAddNeighbour(de,(PetscMPIInt)r);CHKERRQ(ierr);
593     }
594   }
595   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
596   /* Define send counts */
597   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
598   for (r=0; r<commsize; r++) {
599     if (n2collect[r] > 0) {
600       ierr = DataExAddToSendCount(de,r,n2collect[r]);CHKERRQ(ierr);
601     }
602   }
603   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
604   /* Pack data */
605   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
606   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
607   for (r=0; r<commsize; r++) {
608     for (p=0; p<n2collect[r]; p++) {
609       ierr = DataBucketFillPackedArray(swarm->db,collectlist[r][p],point_buffer);CHKERRQ(ierr);
610       /* insert point buffer into the data exchanger */
611       ierr = DataExPackData(de,r,1,point_buffer);CHKERRQ(ierr);
612     }
613   }
614   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
615   /* Scatter */
616   ierr = DataExBegin(de);CHKERRQ(ierr);
617   ierr = DataExEnd(de);CHKERRQ(ierr);
618   /* Collect data in DMSwarm container */
619   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
620   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
621   ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,-1);CHKERRQ(ierr);
622   for (p=0; p<n_points_recv; p++) {
623     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
624 
625     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
626   }
627   /* Release memory */
628   for (r=0; r<commsize; r++) {
629     if (collectlist[r]) PetscFree(collectlist[r]);
630   }
631   ierr = PetscFree(collectlist);CHKERRQ(ierr);
632   ierr = PetscFree(n2collect);CHKERRQ(ierr);
633   ierr = PetscFree(ctxlist);CHKERRQ(ierr);
634   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
635   ierr = DataExView(de);CHKERRQ(ierr);
636   ierr = DataExDestroy(de);CHKERRQ(ierr);
637   PetscFunctionReturn(0);
638 }
639 
640