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