xref: /petsc/src/dm/impls/swarm/swarm_migrate.c (revision a4f09dd6775dcf4381066c075a763b56863c40b8)
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,DATA_BUCKET_BUFFER_DEFAULT);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,DATA_BUCKET_BUFFER_DEFAULT);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     DataField PField;
225     PetscInt npoints_curr;
226 
227     /* remove points which the domain */
228     ierr = DataBucketGetDataFieldByName(swarm->db,DMSwarmField_rank,&PField);CHKERRQ(ierr);
229     ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr);
230 
231     ierr = DataBucketGetSizes(swarm->db,&npoints_curr,NULL,NULL);CHKERRQ(ierr);
232     for (p=0; p<npoints_curr; p++) {
233       if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
234         /* kill point */
235         ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr);
236         ierr = DataBucketGetSizes(swarm->db,&npoints_curr,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */
237         ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); /* update date point increase realloc performed */
238         p--; /* check replacement point */
239       }
240     }
241     ierr = DMSwarmGetSize(dm,&npoints_prior_migration);CHKERRQ(ierr);
242 
243   }
244 
245   /* locate points newly recevied */
246   ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr);
247 
248 #if 0
249   { // safe alternative - however this performs two point locations on: (i) the intial points set and; (ii) the (intial + recieved) point set
250     PetscScalar *LA_coor;
251     PetscInt bs;
252     DataField PField;
253 
254     ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr);
255     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,bs*npoints2,(const PetscScalar*)LA_coor,&pos);CHKERRQ(ierr);
256     ierr = DMLocatePoints(dmcell,pos,DM_POINTLOCATION_NONE,&sfcell);CHKERRQ(ierr);
257 
258     ierr = VecDestroy(&pos);CHKERRQ(ierr);
259     ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr);
260 
261     ierr = PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell);CHKERRQ(ierr);
262     ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
263     for (p=0; p<npoints2; p++) {
264       rankval[p] = LA_sfcell[p].index;
265     }
266     ierr = PetscSFDestroy(&sfcell);CHKERRQ(ierr);
267     ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
268 
269     /* remove points which left processor */
270     ierr = DataBucketGetDataFieldByName(swarm->db,DMSwarmField_rank,&PField);CHKERRQ(ierr);
271     ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr);
272 
273     ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr);
274     for (p=0; p<npoints2; p++) {
275       if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
276         /* kill point */
277         ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr);
278         ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */
279         ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); /* update date point increase realloc performed */
280         p--; /* check replacement point */
281       }
282     }
283   }
284 #endif
285 
286   { // this performs two point locations: (i) on the intial points set prior to communication; and (ii) on the new (recieved) points
287     PetscScalar *LA_coor;
288     PetscInt npoints_from_neighbours,bs;
289     DataField PField;
290 
291     npoints_from_neighbours = npoints2 - npoints_prior_migration;
292 
293     ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr);
294     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,bs*npoints_from_neighbours,(const PetscScalar*)&LA_coor[bs*npoints_prior_migration],&pos);CHKERRQ(ierr);
295 
296     ierr = DMLocatePoints(dmcell,pos,DM_POINTLOCATION_NONE,&sfcell);CHKERRQ(ierr);
297 
298     ierr = VecDestroy(&pos);CHKERRQ(ierr);
299     ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,&bs,NULL,(void**)&LA_coor);CHKERRQ(ierr);
300 
301     ierr = PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell);CHKERRQ(ierr);
302     ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
303     for (p=0; p<npoints_from_neighbours; p++) {
304       rankval[npoints_prior_migration + p] = LA_sfcell[p].index;
305     }
306     ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
307     ierr = PetscSFDestroy(&sfcell);CHKERRQ(ierr);
308 
309     /* remove points which left processor */
310     ierr = DataBucketGetDataFieldByName(swarm->db,DMSwarmField_rank,&PField);CHKERRQ(ierr);
311     ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr);
312 
313     ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr);
314     for (p=npoints_prior_migration; p<npoints2; p++) {
315       if (rankval[p] == DMLOCATEPOINT_POINT_NOT_FOUND) {
316         /* kill point */
317         ierr = DataBucketRemovePointAtIndex(swarm->db,p);CHKERRQ(ierr);
318         ierr = DataBucketGetSizes(swarm->db,&npoints2,NULL,NULL);CHKERRQ(ierr); /* you need to update npoints as the list size decreases! */
319         ierr = DataFieldGetEntries(PField,(void**)&rankval);CHKERRQ(ierr); /* update date point increase realloc performed */
320         p--; /* check replacement point */
321       }
322     }
323   }
324 
325   /* check for error on removed points */
326   if (error_check) {
327     ierr = DMSwarmGetSize(dm,&npoints2g);CHKERRQ(ierr);
328     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);
329   }
330   PetscFunctionReturn(0);
331 }
332 
333 PetscErrorCode DMSwarmMigrate_CellDMExact(DM dm,PetscBool remove_sent_points)
334 {
335   PetscFunctionBegin;
336   PetscFunctionReturn(0);
337 }
338 
339 /*
340  Redundant as this assumes points can only be sent to a single rank
341 */
342 PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm,PetscInt *globalsize)
343 {
344   DM_Swarm *swarm = (DM_Swarm*)dm->data;
345   PetscErrorCode ierr;
346   DataEx de;
347   PetscInt p,npoints,*rankval,n_points_recv;
348   PetscMPIInt rank,nrank,negrank;
349   void *point_buffer,*recv_points;
350   size_t sizeof_dmswarm_point;
351 
352   PetscFunctionBegin;
353   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
354   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
355   *globalsize = npoints;
356   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
357   ierr = DataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr);
358   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
359   for (p=0; p<npoints; p++) {
360     negrank = rankval[p];
361     if (negrank < 0) {
362       nrank = -negrank - 1;
363       ierr = DataExTopologyAddNeighbour(de,nrank);CHKERRQ(ierr);
364     }
365   }
366   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
367   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
368   for (p=0; p<npoints; p++) {
369     negrank = rankval[p];
370     if (negrank < 0) {
371       nrank = -negrank - 1;
372       ierr = DataExAddToSendCount(de,nrank,1);CHKERRQ(ierr);
373     }
374   }
375   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
376   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
377   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
378   for (p=0; p<npoints; p++) {
379     negrank = rankval[p];
380     if (negrank < 0) {
381       nrank = -negrank - 1;
382       rankval[p] = nrank;
383       /* copy point into buffer */
384       ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr);
385       /* insert point buffer into DataExchanger */
386       ierr = DataExPackData(de,nrank,1,point_buffer);CHKERRQ(ierr);
387       rankval[p] = negrank;
388     }
389   }
390   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
391   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
392   ierr = DataExBegin(de);CHKERRQ(ierr);
393   ierr = DataExEnd(de);CHKERRQ(ierr);
394   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
395   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
396   ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr);
397   for (p=0; p<n_points_recv; p++) {
398     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
399 
400     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
401   }
402   ierr = DataExView(de);CHKERRQ(ierr);
403   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
404   ierr = DataExDestroy(de);CHKERRQ(ierr);
405   PetscFunctionReturn(0);
406 }
407 
408 typedef struct {
409   PetscMPIInt owner_rank;
410   PetscReal min[3],max[3];
411 } CollectBBox;
412 
413 PETSC_EXTERN PetscErrorCode DMSwarmCollect_DMDABoundingBox(DM dm,PetscInt *globalsize)
414 {
415   DM_Swarm *swarm = (DM_Swarm*)dm->data;
416   PetscErrorCode ierr;
417   DataEx de;
418   PetscInt p,pk,npoints,*rankval,n_points_recv,n_bbox_recv,dim,neighbour_cells;
419   PetscMPIInt rank,nrank;
420   void *point_buffer,*recv_points;
421   size_t sizeof_dmswarm_point,sizeof_bbox_ctx;
422   PetscBool isdmda;
423   CollectBBox *bbox,*recv_bbox;
424   const PetscMPIInt *dmneighborranks;
425   DM dmcell;
426 
427   PetscFunctionBegin;
428   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
429 
430   ierr = DMSwarmGetCellDM(dm,&dmcell);CHKERRQ(ierr);
431   if (!dmcell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid if cell DM provided");
432   isdmda = PETSC_FALSE;
433   PetscObjectTypeCompare((PetscObject)dmcell,DMDA,&isdmda);
434   if (!isdmda) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only DMDA support for CollectBoundingBox");
435 
436   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
437   sizeof_bbox_ctx = sizeof(CollectBBox);
438   PetscMalloc1(1,&bbox);
439   bbox->owner_rank = rank;
440 
441   /* compute the bounding box based on the overlapping / stenctil size */
442   {
443     Vec lcoor;
444 
445     ierr = DMGetCoordinatesLocal(dmcell,&lcoor);CHKERRQ(ierr);
446     if (dim >= 1) {
447       ierr = VecStrideMin(lcoor,0,NULL,&bbox->min[0]);CHKERRQ(ierr);
448       ierr = VecStrideMax(lcoor,0,NULL,&bbox->max[0]);CHKERRQ(ierr);
449     }
450     if (dim >= 2) {
451       ierr = VecStrideMin(lcoor,1,NULL,&bbox->min[1]);CHKERRQ(ierr);
452       ierr = VecStrideMax(lcoor,1,NULL,&bbox->max[1]);CHKERRQ(ierr);
453     }
454     if (dim == 3) {
455       ierr = VecStrideMin(lcoor,2,NULL,&bbox->min[2]);CHKERRQ(ierr);
456       ierr = VecStrideMax(lcoor,2,NULL,&bbox->max[2]);CHKERRQ(ierr);
457     }
458   }
459   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
460   *globalsize = npoints;
461   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
462   ierr = DataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr);
463   /* use DMDA neighbours */
464   ierr = DMDAGetNeighbors(dmcell,&dmneighborranks);CHKERRQ(ierr);
465   if (dim == 1) {
466     neighbour_cells = 3;
467   } else if (dim == 2) {
468     neighbour_cells = 9;
469   } else {
470     neighbour_cells = 27;
471   }
472   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
473   for (p=0; p<neighbour_cells; p++) {
474     if ( (dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank) ) {
475       ierr = DataExTopologyAddNeighbour(de,dmneighborranks[p]);CHKERRQ(ierr);
476     }
477   }
478   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
479   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
480   for (p=0; p<neighbour_cells; p++) {
481     if ( (dmneighborranks[p] >= 0) && (dmneighborranks[p] != rank) ) {
482       ierr = DataExAddToSendCount(de,dmneighborranks[p],1);CHKERRQ(ierr);
483     }
484   }
485   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
486   /* send bounding boxes */
487   ierr = DataExPackInitialize(de,sizeof_bbox_ctx);CHKERRQ(ierr);
488   for (p=0; p<neighbour_cells; p++) {
489     nrank = dmneighborranks[p];
490     if ( (nrank >= 0) && (nrank != rank) ) {
491       /* insert bbox buffer into DataExchanger */
492       ierr = DataExPackData(de,nrank,1,bbox);CHKERRQ(ierr);
493     }
494   }
495   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
496   /* recv bounding boxes */
497   ierr = DataExBegin(de);CHKERRQ(ierr);
498   ierr = DataExEnd(de);CHKERRQ(ierr);
499   ierr = DataExGetRecvData(de,&n_bbox_recv,(void**)&recv_bbox);CHKERRQ(ierr);
500   for (p=0; p<n_bbox_recv; p++) {
501     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,
502            (double)recv_bbox[p].min[0],(double)recv_bbox[p].max[0],(double)recv_bbox[p].min[1],(double)recv_bbox[p].max[1]);
503   }
504   /* of course this is stupid as this "generic" function should have a better way to know what the coordinates are called */
505   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
506   for (pk=0; pk<n_bbox_recv; pk++) {
507     PetscReal *array_x,*array_y;
508 
509     ierr = DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
510     ierr = DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
511     for (p=0; p<npoints; p++) {
512       if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0]) ) {
513         if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1]) ) {
514           ierr = DataExAddToSendCount(de,recv_bbox[pk].owner_rank,1);CHKERRQ(ierr);
515         }
516       }
517     }
518     ierr = DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
519     ierr = DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
520   }
521   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
522   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
523   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
524   for (pk=0; pk<n_bbox_recv; pk++) {
525     PetscReal *array_x,*array_y;
526 
527     ierr = DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
528     ierr = DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
529     for (p=0; p<npoints; p++) {
530       if ((array_x[p] >= recv_bbox[pk].min[0]) && (array_x[p] <= recv_bbox[pk].max[0]) ) {
531         if ((array_y[p] >= recv_bbox[pk].min[1]) && (array_y[p] <= recv_bbox[pk].max[1]) ) {
532           /* copy point into buffer */
533           ierr = DataBucketFillPackedArray(swarm->db,p,point_buffer);CHKERRQ(ierr);
534           /* insert point buffer into DataExchanger */
535           ierr = DataExPackData(de,recv_bbox[pk].owner_rank,1,point_buffer);CHKERRQ(ierr);
536         }
537       }
538     }
539     ierr = DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
540     ierr = DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
541   }
542   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
543   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
544   ierr = DataExBegin(de);CHKERRQ(ierr);
545   ierr = DataExEnd(de);CHKERRQ(ierr);
546   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
547   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
548   ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr);
549   for (p=0; p<n_points_recv; p++) {
550     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
551 
552     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
553   }
554   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
555   PetscFree(bbox);
556   ierr = DataExView(de);CHKERRQ(ierr);
557   ierr = DataExDestroy(de);CHKERRQ(ierr);
558   PetscFunctionReturn(0);
559 }
560 
561 
562 /* General collection when no order, or neighbour information is provided */
563 /*
564  User provides context and collect() method
565  Broadcast user context
566 
567  for each context / rank {
568    collect(swarm,context,n,list)
569  }
570 */
571 PETSC_EXTERN PetscErrorCode DMSwarmCollect_General(DM dm,PetscErrorCode (*collect)(DM,void*,PetscInt*,PetscInt**),size_t ctx_size,void *ctx,PetscInt *globalsize)
572 {
573   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
574   PetscErrorCode ierr;
575   DataEx         de;
576   PetscInt       p,r,npoints,n_points_recv;
577   PetscMPIInt    commsize,rank;
578   void           *point_buffer,*recv_points;
579   void           *ctxlist;
580   PetscInt       *n2collect,**collectlist;
581   size_t         sizeof_dmswarm_point;
582 
583   PetscFunctionBegin;
584   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&commsize);CHKERRQ(ierr);
585   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
586   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
587   *globalsize = npoints;
588   /* Broadcast user context */
589   PetscMalloc(ctx_size*commsize,&ctxlist);
590   ierr = MPI_Allgather(ctx,ctx_size,MPI_CHAR,ctxlist,ctx_size,MPI_CHAR,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
591   ierr = PetscMalloc1(commsize,&n2collect);CHKERRQ(ierr);
592   ierr = PetscMalloc1(commsize,&collectlist);CHKERRQ(ierr);
593   for (r=0; r<commsize; r++) {
594     PetscInt _n2collect;
595     PetscInt *_collectlist;
596     void     *_ctx_r;
597 
598     _n2collect   = 0;
599     _collectlist = NULL;
600     if (r != rank) { /* don't collect data from yourself */
601       _ctx_r = (void*)( (char*)ctxlist + r * ctx_size );
602       ierr = collect(dm,_ctx_r,&_n2collect,&_collectlist);CHKERRQ(ierr);
603     }
604     n2collect[r]   = _n2collect;
605     collectlist[r] = _collectlist;
606   }
607   ierr = DataExCreate(PetscObjectComm((PetscObject)dm),0,&de);CHKERRQ(ierr);
608   /* Define topology */
609   ierr = DataExTopologyInitialize(de);CHKERRQ(ierr);
610   for (r=0; r<commsize; r++) {
611     if (n2collect[r] > 0) {
612       ierr = DataExTopologyAddNeighbour(de,(PetscMPIInt)r);CHKERRQ(ierr);
613     }
614   }
615   ierr = DataExTopologyFinalize(de);CHKERRQ(ierr);
616   /* Define send counts */
617   ierr = DataExInitializeSendCount(de);CHKERRQ(ierr);
618   for (r=0; r<commsize; r++) {
619     if (n2collect[r] > 0) {
620       ierr = DataExAddToSendCount(de,r,n2collect[r]);CHKERRQ(ierr);
621     }
622   }
623   ierr = DataExFinalizeSendCount(de);CHKERRQ(ierr);
624   /* Pack data */
625   ierr = DataBucketCreatePackedArray(swarm->db,&sizeof_dmswarm_point,&point_buffer);CHKERRQ(ierr);
626   ierr = DataExPackInitialize(de,sizeof_dmswarm_point);CHKERRQ(ierr);
627   for (r=0; r<commsize; r++) {
628     for (p=0; p<n2collect[r]; p++) {
629       ierr = DataBucketFillPackedArray(swarm->db,collectlist[r][p],point_buffer);CHKERRQ(ierr);
630       /* insert point buffer into the data exchanger */
631       ierr = DataExPackData(de,r,1,point_buffer);CHKERRQ(ierr);
632     }
633   }
634   ierr = DataExPackFinalize(de);CHKERRQ(ierr);
635   /* Scatter */
636   ierr = DataExBegin(de);CHKERRQ(ierr);
637   ierr = DataExEnd(de);CHKERRQ(ierr);
638   /* Collect data in DMSwarm container */
639   ierr = DataExGetRecvData(de,&n_points_recv,(void**)&recv_points);CHKERRQ(ierr);
640   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
641   ierr = DataBucketSetSizes(swarm->db,npoints + n_points_recv,DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr);
642   for (p=0; p<n_points_recv; p++) {
643     void *data_p = (void*)( (char*)recv_points + p*sizeof_dmswarm_point );
644 
645     ierr = DataBucketInsertPackedArray(swarm->db,npoints+p,data_p);CHKERRQ(ierr);
646   }
647   /* Release memory */
648   for (r=0; r<commsize; r++) {
649     if (collectlist[r]) PetscFree(collectlist[r]);
650   }
651   ierr = PetscFree(collectlist);CHKERRQ(ierr);
652   ierr = PetscFree(n2collect);CHKERRQ(ierr);
653   ierr = PetscFree(ctxlist);CHKERRQ(ierr);
654   ierr = DataBucketDestroyPackedArray(swarm->db,&point_buffer);CHKERRQ(ierr);
655   ierr = DataExView(de);CHKERRQ(ierr);
656   ierr = DataExDestroy(de);CHKERRQ(ierr);
657   PetscFunctionReturn(0);
658 }
659 
660