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