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