1 #include "../src/dm/impls/swarm/data_bucket.h"
2 
3 /* string helpers */
DMSwarmDataFieldStringInList(const char name[],const PetscInt N,const DMSwarmDataField gfield[],PetscBool * val)4 PetscErrorCode DMSwarmDataFieldStringInList(const char name[],const PetscInt N,const DMSwarmDataField gfield[],PetscBool *val)
5 {
6   PetscInt       i;
7   PetscErrorCode ierr;
8 
9   PetscFunctionBegin;
10   *val = PETSC_FALSE;
11   for (i = 0; i < N; ++i) {
12     PetscBool flg;
13     ierr = PetscStrcmp(name, gfield[i]->name, &flg);CHKERRQ(ierr);
14     if (flg) {
15       *val = PETSC_TRUE;
16       PetscFunctionReturn(0);
17     }
18   }
19   PetscFunctionReturn(0);
20 }
21 
DMSwarmDataFieldStringFindInList(const char name[],const PetscInt N,const DMSwarmDataField gfield[],PetscInt * index)22 PetscErrorCode DMSwarmDataFieldStringFindInList(const char name[],const PetscInt N,const DMSwarmDataField gfield[],PetscInt *index)
23 {
24   PetscInt       i;
25   PetscErrorCode ierr;
26 
27   PetscFunctionBegin;
28   *index = -1;
29   for (i = 0; i < N; ++i) {
30     PetscBool flg;
31     ierr = PetscStrcmp(name, gfield[i]->name, &flg);CHKERRQ(ierr);
32     if (flg) {
33       *index = i;
34       PetscFunctionReturn(0);
35     }
36   }
37   PetscFunctionReturn(0);
38 }
39 
DMSwarmDataFieldCreate(const char registration_function[],const char name[],const size_t size,const PetscInt L,DMSwarmDataField * DF)40 PetscErrorCode DMSwarmDataFieldCreate(const char registration_function[],const char name[],const size_t size,const PetscInt L,DMSwarmDataField *DF)
41 {
42   DMSwarmDataField df;
43   PetscErrorCode   ierr;
44 
45   PetscFunctionBegin;
46   ierr = PetscMalloc(sizeof(struct _p_DMSwarmDataField), &df);CHKERRQ(ierr);
47   ierr = PetscMemzero(df, sizeof(struct _p_DMSwarmDataField));CHKERRQ(ierr);
48   ierr = PetscStrallocpy(registration_function, &df->registration_function);CHKERRQ(ierr);
49   ierr = PetscStrallocpy(name, &df->name);CHKERRQ(ierr);
50   df->atomic_size = size;
51   df->L  = L;
52   df->bs = 1;
53   /* allocate something so we don't have to reallocate */
54   ierr = PetscMalloc(size * L, &df->data);CHKERRQ(ierr);
55   ierr = PetscMemzero(df->data, size * L);CHKERRQ(ierr);
56   *DF = df;
57   PetscFunctionReturn(0);
58 }
59 
DMSwarmDataFieldDestroy(DMSwarmDataField * DF)60 PetscErrorCode DMSwarmDataFieldDestroy(DMSwarmDataField *DF)
61 {
62   DMSwarmDataField df = *DF;
63   PetscErrorCode   ierr;
64 
65   PetscFunctionBegin;
66   ierr = PetscFree(df->registration_function);CHKERRQ(ierr);
67   ierr = PetscFree(df->name);CHKERRQ(ierr);
68   ierr = PetscFree(df->data);CHKERRQ(ierr);
69   ierr = PetscFree(df);CHKERRQ(ierr);
70   *DF  = NULL;
71   PetscFunctionReturn(0);
72 }
73 
74 /* data bucket */
DMSwarmDataBucketCreate(DMSwarmDataBucket * DB)75 PetscErrorCode DMSwarmDataBucketCreate(DMSwarmDataBucket *DB)
76 {
77   DMSwarmDataBucket db;
78   PetscErrorCode    ierr;
79 
80   PetscFunctionBegin;
81   ierr = PetscMalloc(sizeof(struct _p_DMSwarmDataBucket), &db);CHKERRQ(ierr);
82   ierr = PetscMemzero(db, sizeof(struct _p_DMSwarmDataBucket));CHKERRQ(ierr);
83 
84   db->finalised = PETSC_FALSE;
85   /* create empty spaces for fields */
86   db->L         = -1;
87   db->buffer    = 1;
88   db->allocated = 1;
89   db->nfields   = 0;
90   ierr = PetscMalloc1(1, &db->field);CHKERRQ(ierr);
91   *DB  = db;
92   PetscFunctionReturn(0);
93 }
94 
DMSwarmDataBucketDestroy(DMSwarmDataBucket * DB)95 PetscErrorCode DMSwarmDataBucketDestroy(DMSwarmDataBucket *DB)
96 {
97   DMSwarmDataBucket db = *DB;
98   PetscInt          f;
99   PetscErrorCode    ierr;
100 
101   PetscFunctionBegin;
102   /* release fields */
103   for (f = 0; f < db->nfields; ++f) {
104     ierr = DMSwarmDataFieldDestroy(&db->field[f]);CHKERRQ(ierr);
105   }
106   /* this will catch the initially allocated objects in the event that no fields are registered */
107   if (db->field != NULL) {
108     ierr = PetscFree(db->field);CHKERRQ(ierr);
109   }
110   ierr = PetscFree(db);CHKERRQ(ierr);
111   *DB = NULL;
112   PetscFunctionReturn(0);
113 }
114 
DMSwarmDataBucketQueryForActiveFields(DMSwarmDataBucket db,PetscBool * any_active_fields)115 PetscErrorCode DMSwarmDataBucketQueryForActiveFields(DMSwarmDataBucket db,PetscBool *any_active_fields)
116 {
117   PetscInt f;
118 
119   PetscFunctionBegin;
120   *any_active_fields = PETSC_FALSE;
121   for (f = 0; f < db->nfields; ++f) {
122     if (db->field[f]->active) {
123       *any_active_fields = PETSC_TRUE;
124       PetscFunctionReturn(0);
125     }
126   }
127   PetscFunctionReturn(0);
128 }
129 
DMSwarmDataBucketRegisterField(DMSwarmDataBucket db,const char registration_function[],const char field_name[],size_t atomic_size,DMSwarmDataField * _gfield)130 PetscErrorCode DMSwarmDataBucketRegisterField(
131                               DMSwarmDataBucket db,
132                               const char registration_function[],
133                               const char field_name[],
134                               size_t atomic_size, DMSwarmDataField *_gfield)
135 {
136   PetscBool val;
137   DMSwarmDataField fp;
138   PetscErrorCode ierr;
139 
140   PetscFunctionBegin;
141         /* check we haven't finalised the registration of fields */
142         /*
143    if (db->finalised==PETSC_TRUE) {
144    printf("ERROR: DMSwarmDataBucketFinalize() has been called. Cannot register more fields\n");
145    ERROR();
146    }
147    */
148   /* check for repeated name */
149   ierr = DMSwarmDataFieldStringInList(field_name, db->nfields, (const DMSwarmDataField*) db->field, &val);CHKERRQ(ierr);
150   if (val == PETSC_TRUE) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field %s already exists. Cannot add same field twice",field_name);
151   /* create new space for data */
152   ierr = PetscRealloc(sizeof(DMSwarmDataField)*(db->nfields+1), &db->field);CHKERRQ(ierr);
153   /* add field */
154   ierr = DMSwarmDataFieldCreate(registration_function, field_name, atomic_size, db->allocated, &fp);CHKERRQ(ierr);
155   db->field[db->nfields] = fp;
156   db->nfields++;
157   if (_gfield != NULL) {
158     *_gfield = fp;
159   }
160   PetscFunctionReturn(0);
161 }
162 
163 /*
164  #define DMSwarmDataBucketRegisterField(db,name,size,k) {\
165  char *location;\
166  asprintf(&location,"Registered by %s() at line %d within file %s", __FUNCTION__, __LINE__, __FILE__);\
167  _DMSwarmDataBucketRegisterField( (db), location, (name), (size), (k));\
168  ierr = PetscFree(location);\
169  }
170  */
171 
DMSwarmDataBucketGetDMSwarmDataFieldByName(DMSwarmDataBucket db,const char name[],DMSwarmDataField * gfield)172 PetscErrorCode DMSwarmDataBucketGetDMSwarmDataFieldByName(DMSwarmDataBucket db,const char name[],DMSwarmDataField *gfield)
173 {
174   PetscInt       idx;
175   PetscBool      found;
176   PetscErrorCode ierr;
177 
178   PetscFunctionBegin;
179   ierr = DMSwarmDataFieldStringInList(name,db->nfields,(const DMSwarmDataField*)db->field,&found);CHKERRQ(ierr);
180   if (!found) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot find DMSwarmDataField with name %s",name);
181   ierr = DMSwarmDataFieldStringFindInList(name,db->nfields,(const DMSwarmDataField*)db->field,&idx);CHKERRQ(ierr);
182   *gfield = db->field[idx];
183   PetscFunctionReturn(0);
184 }
185 
DMSwarmDataBucketQueryDMSwarmDataFieldByName(DMSwarmDataBucket db,const char name[],PetscBool * found)186 PetscErrorCode DMSwarmDataBucketQueryDMSwarmDataFieldByName(DMSwarmDataBucket db,const char name[],PetscBool *found)
187 {
188   PetscErrorCode ierr;
189 
190   PetscFunctionBegin;
191   *found = PETSC_FALSE;
192   ierr = DMSwarmDataFieldStringInList(name,db->nfields,(const DMSwarmDataField*)db->field,found);CHKERRQ(ierr);
193   PetscFunctionReturn(0);
194 }
195 
DMSwarmDataBucketFinalize(DMSwarmDataBucket db)196 PetscErrorCode DMSwarmDataBucketFinalize(DMSwarmDataBucket db)
197 {
198   PetscFunctionBegin;
199   db->finalised = PETSC_TRUE;
200   PetscFunctionReturn(0);
201 }
202 
DMSwarmDataFieldGetNumEntries(DMSwarmDataField df,PetscInt * sum)203 PetscErrorCode DMSwarmDataFieldGetNumEntries(DMSwarmDataField df,PetscInt *sum)
204 {
205   PetscFunctionBegin;
206   *sum = df->L;
207   PetscFunctionReturn(0);
208 }
209 
DMSwarmDataFieldSetBlockSize(DMSwarmDataField df,PetscInt blocksize)210 PetscErrorCode DMSwarmDataFieldSetBlockSize(DMSwarmDataField df,PetscInt blocksize)
211 {
212   PetscFunctionBegin;
213   df->bs = blocksize;
214   PetscFunctionReturn(0);
215 }
216 
DMSwarmDataFieldSetSize(DMSwarmDataField df,const PetscInt new_L)217 PetscErrorCode DMSwarmDataFieldSetSize(DMSwarmDataField df,const PetscInt new_L)
218 {
219   PetscErrorCode ierr;
220 
221   PetscFunctionBegin;
222   if (new_L < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot set size of DMSwarmDataField to be < 0");
223   if (new_L == df->L) PetscFunctionReturn(0);
224   if (new_L > df->L) {
225     ierr = PetscRealloc(df->atomic_size * (new_L), &df->data);CHKERRQ(ierr);
226     /* init new contents */
227     ierr = PetscMemzero(( ((char*)df->data)+df->L*df->atomic_size), (new_L-df->L)*df->atomic_size);CHKERRQ(ierr);
228   } else {
229     /* reallocate pointer list, add +1 in case new_L = 0 */
230     ierr = PetscRealloc(df->atomic_size * (new_L+1), &df->data);CHKERRQ(ierr);
231   }
232   df->L = new_L;
233   PetscFunctionReturn(0);
234 }
235 
DMSwarmDataFieldZeroBlock(DMSwarmDataField df,const PetscInt start,const PetscInt end)236 PetscErrorCode DMSwarmDataFieldZeroBlock(DMSwarmDataField df,const PetscInt start,const PetscInt end)
237 {
238   PetscErrorCode ierr;
239 
240   PetscFunctionBegin;
241   if (start > end) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if start(%D) > end(%D)",start,end);
242   if (start < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if start(%D) < 0",start);
243   if (end > df->L) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if end(%D) >= array size(%D)",end,df->L);
244   ierr = PetscMemzero((((char*)df->data)+start*df->atomic_size), (end-start)*df->atomic_size);CHKERRQ(ierr);
245   PetscFunctionReturn(0);
246 }
247 
248 /*
249  A negative buffer value will simply be ignored and the old buffer value will be used.
250  */
DMSwarmDataBucketSetSizes(DMSwarmDataBucket db,const PetscInt L,const PetscInt buffer)251 PetscErrorCode DMSwarmDataBucketSetSizes(DMSwarmDataBucket db,const PetscInt L,const PetscInt buffer)
252 {
253   PetscInt       current_allocated,new_used,new_unused,new_buffer,new_allocated,f;
254   PetscBool      any_active_fields;
255   PetscErrorCode ierr;
256 
257   PetscFunctionBegin;
258   if (db->finalised == PETSC_FALSE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"You must call DMSwarmDataBucketFinalize() before DMSwarmDataBucketSetSizes()");
259   ierr = DMSwarmDataBucketQueryForActiveFields(db,&any_active_fields);CHKERRQ(ierr);
260   if (any_active_fields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot safely re-size as at least one DMSwarmDataField is currently being accessed");
261 
262   current_allocated = db->allocated;
263   new_used   = L;
264   new_unused = current_allocated - new_used;
265   new_buffer = db->buffer;
266   if (buffer >= 0) { /* update the buffer value */
267     new_buffer = buffer;
268   }
269   new_allocated = new_used + new_buffer;
270   /* action */
271   if (new_allocated > current_allocated) {
272     /* increase size to new_used + new_buffer */
273     for (f=0; f<db->nfields; f++) {
274       ierr = DMSwarmDataFieldSetSize(db->field[f], new_allocated);CHKERRQ(ierr);
275     }
276     db->L         = new_used;
277     db->buffer    = new_buffer;
278     db->allocated = new_used + new_buffer;
279   } else {
280     if (new_unused > 2 * new_buffer) {
281       /* shrink array to new_used + new_buffer */
282       for (f = 0; f < db->nfields; ++f) {
283         ierr = DMSwarmDataFieldSetSize(db->field[f], new_allocated);CHKERRQ(ierr);
284       }
285       db->L         = new_used;
286       db->buffer    = new_buffer;
287       db->allocated = new_used + new_buffer;
288     } else {
289       db->L      = new_used;
290       db->buffer = new_buffer;
291     }
292   }
293   /* zero all entries from db->L to db->allocated */
294   for (f = 0; f < db->nfields; ++f) {
295     DMSwarmDataField field = db->field[f];
296     ierr = DMSwarmDataFieldZeroBlock(field, db->L,db->allocated);CHKERRQ(ierr);
297   }
298   PetscFunctionReturn(0);
299 }
300 
DMSwarmDataBucketSetInitialSizes(DMSwarmDataBucket db,const PetscInt L,const PetscInt buffer)301 PetscErrorCode DMSwarmDataBucketSetInitialSizes(DMSwarmDataBucket db,const PetscInt L,const PetscInt buffer)
302 {
303   PetscInt       f;
304   PetscErrorCode ierr;
305 
306   PetscFunctionBegin;
307   ierr = DMSwarmDataBucketSetSizes(db,L,buffer);CHKERRQ(ierr);
308   for (f = 0; f < db->nfields; ++f) {
309     DMSwarmDataField field = db->field[f];
310     ierr = DMSwarmDataFieldZeroBlock(field,0,db->allocated);CHKERRQ(ierr);
311   }
312   PetscFunctionReturn(0);
313 }
314 
DMSwarmDataBucketGetSizes(DMSwarmDataBucket db,PetscInt * L,PetscInt * buffer,PetscInt * allocated)315 PetscErrorCode DMSwarmDataBucketGetSizes(DMSwarmDataBucket db,PetscInt *L,PetscInt *buffer,PetscInt *allocated)
316 {
317   PetscFunctionBegin;
318   if (L) {*L = db->L;}
319   if (buffer) {*buffer = db->buffer;}
320   if (allocated) {*allocated = db->allocated;}
321   PetscFunctionReturn(0);
322 }
323 
DMSwarmDataBucketGetGlobalSizes(MPI_Comm comm,DMSwarmDataBucket db,PetscInt * L,PetscInt * buffer,PetscInt * allocated)324 PetscErrorCode DMSwarmDataBucketGetGlobalSizes(MPI_Comm comm,DMSwarmDataBucket db,PetscInt *L,PetscInt *buffer,PetscInt *allocated)
325 {
326   PetscInt ierr;
327 
328   PetscFunctionBegin;
329   if (L) {         ierr = MPI_Allreduce(&db->L,L,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); }
330   if (buffer) {    ierr = MPI_Allreduce(&db->buffer,buffer,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); }
331   if (allocated) { ierr = MPI_Allreduce(&db->allocated,allocated,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); }
332   PetscFunctionReturn(0);
333 }
334 
DMSwarmDataBucketGetDMSwarmDataFields(DMSwarmDataBucket db,PetscInt * L,DMSwarmDataField * fields[])335 PetscErrorCode DMSwarmDataBucketGetDMSwarmDataFields(DMSwarmDataBucket db,PetscInt *L,DMSwarmDataField *fields[])
336 {
337   PetscFunctionBegin;
338   if (L)      {*L      = db->nfields;}
339   if (fields) {*fields = db->field;}
340   PetscFunctionReturn(0);
341 }
342 
DMSwarmDataFieldGetAccess(const DMSwarmDataField gfield)343 PetscErrorCode DMSwarmDataFieldGetAccess(const DMSwarmDataField gfield)
344 {
345   PetscFunctionBegin;
346   if (gfield->active) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is already active. You must call DMSwarmDataFieldRestoreAccess()",gfield->name);
347   gfield->active = PETSC_TRUE;
348   PetscFunctionReturn(0);
349 }
350 
DMSwarmDataFieldAccessPoint(const DMSwarmDataField gfield,const PetscInt pid,void ** ctx_p)351 PetscErrorCode DMSwarmDataFieldAccessPoint(const DMSwarmDataField gfield,const PetscInt pid,void **ctx_p)
352 {
353   PetscFunctionBegin;
354   *ctx_p = NULL;
355 #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
356   /* debug mode */
357   /* check point is valid */
358   if (pid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
359   if (pid >= gfield->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",gfield->L);
360   if (gfield->active == PETSC_FALSE) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is not active. You must call DMSwarmDataFieldGetAccess() before point data can be retrivied",gfield->name);
361 #endif
362   *ctx_p = DMSWARM_DATAFIELD_point_access(gfield->data,pid,gfield->atomic_size);
363   PetscFunctionReturn(0);
364 }
365 
DMSwarmDataFieldAccessPointOffset(const DMSwarmDataField gfield,const size_t offset,const PetscInt pid,void ** ctx_p)366 PetscErrorCode DMSwarmDataFieldAccessPointOffset(const DMSwarmDataField gfield,const size_t offset,const PetscInt pid,void **ctx_p)
367 {
368   PetscFunctionBegin;
369 #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
370   /* debug mode */
371   /* check point is valid */
372   /* if (offset < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be >= 0");*/
373   /* Note compiler realizes this can never happen with an unsigned PetscInt */
374   if (offset >= gfield->atomic_size) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be < %zu",gfield->atomic_size);
375   /* check point is valid */
376   if (pid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
377   if (pid >= gfield->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",gfield->L);
378   if (gfield->active == PETSC_FALSE) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is not active. You must call DMSwarmDataFieldGetAccess() before point data can be retrivied",gfield->name);
379 #endif
380   *ctx_p = DMSWARM_DATAFIELD_point_access_offset(gfield->data,pid,gfield->atomic_size,offset);
381   PetscFunctionReturn(0);
382 }
383 
DMSwarmDataFieldRestoreAccess(DMSwarmDataField gfield)384 PetscErrorCode DMSwarmDataFieldRestoreAccess(DMSwarmDataField gfield)
385 {
386   PetscFunctionBegin;
387   if (gfield->active == PETSC_FALSE) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is not active. You must call DMSwarmDataFieldGetAccess()", gfield->name);
388   gfield->active = PETSC_FALSE;
389   PetscFunctionReturn(0);
390 }
391 
DMSwarmDataFieldVerifyAccess(const DMSwarmDataField gfield,const size_t size)392 PetscErrorCode DMSwarmDataFieldVerifyAccess(const DMSwarmDataField gfield,const size_t size)
393 {
394   PetscFunctionBegin;
395 #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
396   if (gfield->atomic_size != size) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" must be mapped to %zu bytes, your intended structure is %zu bytes in length.",gfield->name, gfield->atomic_size, size);
397 #endif
398   PetscFunctionReturn(0);
399 }
400 
DMSwarmDataFieldGetAtomicSize(const DMSwarmDataField gfield,size_t * size)401 PetscErrorCode DMSwarmDataFieldGetAtomicSize(const DMSwarmDataField gfield,size_t *size)
402 {
403   PetscFunctionBegin;
404   if (size) {*size = gfield->atomic_size;}
405   PetscFunctionReturn(0);
406 }
407 
DMSwarmDataFieldGetEntries(const DMSwarmDataField gfield,void ** data)408 PetscErrorCode DMSwarmDataFieldGetEntries(const DMSwarmDataField gfield,void **data)
409 {
410   PetscFunctionBegin;
411   if (data) {*data = gfield->data;}
412   PetscFunctionReturn(0);
413 }
414 
DMSwarmDataFieldRestoreEntries(const DMSwarmDataField gfield,void ** data)415 PetscErrorCode DMSwarmDataFieldRestoreEntries(const DMSwarmDataField gfield,void **data)
416 {
417   PetscFunctionBegin;
418   if (data) {*data = NULL;}
419   PetscFunctionReturn(0);
420 }
421 
422 /* y = x */
DMSwarmDataBucketCopyPoint(const DMSwarmDataBucket xb,const PetscInt pid_x,const DMSwarmDataBucket yb,const PetscInt pid_y)423 PetscErrorCode DMSwarmDataBucketCopyPoint(const DMSwarmDataBucket xb,const PetscInt pid_x,
424                          const DMSwarmDataBucket yb,const PetscInt pid_y)
425 {
426   PetscInt f;
427   PetscErrorCode ierr;
428 
429   PetscFunctionBegin;
430   for (f = 0; f < xb->nfields; ++f) {
431     void *dest;
432     void *src;
433 
434     ierr = DMSwarmDataFieldGetAccess(xb->field[f]);CHKERRQ(ierr);
435     if (xb != yb) { ierr = DMSwarmDataFieldGetAccess( yb->field[f]);CHKERRQ(ierr); }
436     ierr = DMSwarmDataFieldAccessPoint(xb->field[f],pid_x, &src);CHKERRQ(ierr);
437     ierr = DMSwarmDataFieldAccessPoint(yb->field[f],pid_y, &dest);CHKERRQ(ierr);
438     ierr = PetscMemcpy(dest, src, xb->field[f]->atomic_size);CHKERRQ(ierr);
439     ierr = DMSwarmDataFieldRestoreAccess(xb->field[f]);CHKERRQ(ierr);
440     if (xb != yb) {ierr = DMSwarmDataFieldRestoreAccess(yb->field[f]);CHKERRQ(ierr);}
441   }
442   PetscFunctionReturn(0);
443 }
444 
DMSwarmDataBucketCreateFromSubset(DMSwarmDataBucket DBIn,const PetscInt N,const PetscInt list[],DMSwarmDataBucket * DB)445 PetscErrorCode DMSwarmDataBucketCreateFromSubset(DMSwarmDataBucket DBIn,const PetscInt N,const PetscInt list[],DMSwarmDataBucket *DB)
446 {
447   PetscInt nfields;
448   DMSwarmDataField *fields;
449   PetscInt f,L,buffer,allocated,p;
450   PetscErrorCode ierr;
451 
452   PetscFunctionBegin;
453   ierr = DMSwarmDataBucketCreate(DB);CHKERRQ(ierr);
454   /* copy contents of DBIn */
455   ierr = DMSwarmDataBucketGetDMSwarmDataFields(DBIn,&nfields,&fields);CHKERRQ(ierr);
456   ierr = DMSwarmDataBucketGetSizes(DBIn,&L,&buffer,&allocated);CHKERRQ(ierr);
457   for (f = 0; f < nfields; ++f) {
458     ierr = DMSwarmDataBucketRegisterField(*DB,"DMSwarmDataBucketCreateFromSubset",fields[f]->name,fields[f]->atomic_size,NULL);CHKERRQ(ierr);
459   }
460   ierr = DMSwarmDataBucketFinalize(*DB);CHKERRQ(ierr);
461   ierr = DMSwarmDataBucketSetSizes(*DB,L,buffer);CHKERRQ(ierr);
462   /* now copy the desired guys from DBIn => DB */
463   for (p = 0; p < N; ++p) {
464     ierr = DMSwarmDataBucketCopyPoint(DBIn,list[p], *DB,p);CHKERRQ(ierr);
465   }
466   PetscFunctionReturn(0);
467 }
468 
469 /* insert into an exisitng location */
DMSwarmDataFieldInsertPoint(const DMSwarmDataField field,const PetscInt index,const void * ctx)470 PetscErrorCode DMSwarmDataFieldInsertPoint(const DMSwarmDataField field,const PetscInt index,const void *ctx)
471 {
472   PetscErrorCode ierr;
473 
474   PetscFunctionBegin;
475 #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
476   /* check point is valid */
477   if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
478   if (index >= field->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",field->L);
479 #endif
480   ierr = PetscMemcpy(DMSWARM_DATAFIELD_point_access(field->data,index,field->atomic_size), ctx, field->atomic_size);CHKERRQ(ierr);
481   PetscFunctionReturn(0);
482 }
483 
484 /* remove data at index - replace with last point */
DMSwarmDataBucketRemovePointAtIndex(const DMSwarmDataBucket db,const PetscInt index)485 PetscErrorCode DMSwarmDataBucketRemovePointAtIndex(const DMSwarmDataBucket db,const PetscInt index)
486 {
487   PetscInt       f;
488   PetscBool      any_active_fields;
489   PetscErrorCode ierr;
490 
491   PetscFunctionBegin;
492 #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
493   /* check point is valid */
494   if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
495   if (index >= db->allocated) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",db->L+db->buffer);
496 #endif
497   ierr = DMSwarmDataBucketQueryForActiveFields(db,&any_active_fields);CHKERRQ(ierr);
498   if (any_active_fields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot safely remove point as at least one DMSwarmDataField is currently being accessed");
499   if (index >= db->L) { /* this point is not in the list - no need to error, but I will anyway */
500     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"You should not be trying to remove point at index=%D since it's < db->L = %D", index, db->L);
501   }
502   if (index != db->L-1) { /* not last point in list */
503     for (f = 0; f < db->nfields; ++f) {
504       DMSwarmDataField field = db->field[f];
505 
506       /* copy then remove */
507       ierr = DMSwarmDataFieldCopyPoint(db->L-1, field, index, field);CHKERRQ(ierr);
508       /* DMSwarmDataFieldZeroPoint(field,index); */
509     }
510   }
511   /* decrement size */
512   /* this will zero out an crap at the end of the list */
513   ierr = DMSwarmDataBucketRemovePoint(db);CHKERRQ(ierr);
514   PetscFunctionReturn(0);
515 }
516 
517 /* copy x into y */
DMSwarmDataFieldCopyPoint(const PetscInt pid_x,const DMSwarmDataField field_x,const PetscInt pid_y,const DMSwarmDataField field_y)518 PetscErrorCode DMSwarmDataFieldCopyPoint(const PetscInt pid_x,const DMSwarmDataField field_x,
519                         const PetscInt pid_y,const DMSwarmDataField field_y)
520 {
521   PetscErrorCode ierr;
522 
523   PetscFunctionBegin;
524 #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
525   /* check point is valid */
526   if (pid_x < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"(IN) index must be >= 0");
527   if (pid_x >= field_x->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"(IN) index must be < %D",field_x->L);
528   if (pid_y < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"(OUT) index must be >= 0");
529   if (pid_y >= field_y->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"(OUT) index must be < %D",field_y->L);
530   if (field_y->atomic_size != field_x->atomic_size) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"atomic size must match");
531 #endif
532   ierr = PetscMemcpy(DMSWARM_DATAFIELD_point_access(field_y->data,pid_y,field_y->atomic_size),DMSWARM_DATAFIELD_point_access(field_x->data,pid_x,field_x->atomic_size),field_y->atomic_size);CHKERRQ(ierr);
533   PetscFunctionReturn(0);
534 }
535 
536 
537 /* zero only the datafield at this point */
DMSwarmDataFieldZeroPoint(const DMSwarmDataField field,const PetscInt index)538 PetscErrorCode DMSwarmDataFieldZeroPoint(const DMSwarmDataField field,const PetscInt index)
539 {
540   PetscErrorCode ierr;
541 
542   PetscFunctionBegin;
543 #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
544   /* check point is valid */
545   if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
546   if (index >= field->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",field->L);
547 #endif
548   ierr = PetscMemzero(DMSWARM_DATAFIELD_point_access(field->data,index,field->atomic_size), field->atomic_size);CHKERRQ(ierr);
549   PetscFunctionReturn(0);
550 }
551 
552 /* zero ALL data for this point */
DMSwarmDataBucketZeroPoint(const DMSwarmDataBucket db,const PetscInt index)553 PetscErrorCode DMSwarmDataBucketZeroPoint(const DMSwarmDataBucket db,const PetscInt index)
554 {
555   PetscInt f;
556   PetscErrorCode ierr;
557 
558   PetscFunctionBegin;
559   /* check point is valid */
560   if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
561   if (index >= db->allocated) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",db->allocated);
562   for (f = 0; f < db->nfields; ++f) {
563     DMSwarmDataField field = db->field[f];
564     ierr = DMSwarmDataFieldZeroPoint(field,index);CHKERRQ(ierr);
565   }
566   PetscFunctionReturn(0);
567 }
568 
569 /* increment */
DMSwarmDataBucketAddPoint(DMSwarmDataBucket db)570 PetscErrorCode DMSwarmDataBucketAddPoint(DMSwarmDataBucket db)
571 {
572   PetscErrorCode ierr;
573 
574   PetscFunctionBegin;
575   ierr = DMSwarmDataBucketSetSizes(db,db->L+1,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr);
576   PetscFunctionReturn(0);
577 }
578 
579 /* decrement */
DMSwarmDataBucketRemovePoint(DMSwarmDataBucket db)580 PetscErrorCode DMSwarmDataBucketRemovePoint(DMSwarmDataBucket db)
581 {
582   PetscErrorCode ierr;
583 
584   PetscFunctionBegin;
585   ierr = DMSwarmDataBucketSetSizes(db,db->L-1,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr);
586   PetscFunctionReturn(0);
587 }
588 
DMSwarmDataBucketView_stdout(MPI_Comm comm,DMSwarmDataBucket db)589 PetscErrorCode DMSwarmDataBucketView_stdout(MPI_Comm comm,DMSwarmDataBucket db)
590 {
591   PetscInt       f;
592   double         memory_usage_total,memory_usage_total_local = 0.0;
593   PetscErrorCode ierr;
594 
595   PetscFunctionBegin;
596   ierr = PetscPrintf(comm,"DMSwarmDataBucketView: \n");CHKERRQ(ierr);
597   ierr = PetscPrintf(comm,"  L                  = %D \n", db->L);CHKERRQ(ierr);
598   ierr = PetscPrintf(comm,"  buffer             = %D \n", db->buffer);CHKERRQ(ierr);
599   ierr = PetscPrintf(comm,"  allocated          = %D \n", db->allocated);CHKERRQ(ierr);
600   ierr = PetscPrintf(comm,"  nfields registered = %D \n", db->nfields);CHKERRQ(ierr);
601 
602   for (f = 0; f < db->nfields; ++f) {
603     double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
604     memory_usage_total_local += memory_usage_f;
605   }
606   ierr = MPI_Allreduce(&memory_usage_total_local,&memory_usage_total,1,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
607 
608   for (f = 0; f < db->nfields; ++f) {
609     double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
610     ierr = PetscPrintf(comm,"    [%3D] %15s : Mem. usage       = %1.2e (MB) [rank0]\n", f, db->field[f]->name, memory_usage_f);CHKERRQ(ierr);
611     ierr = PetscPrintf(comm,"                            blocksize        = %D \n", db->field[f]->bs);CHKERRQ(ierr);
612     if (db->field[f]->bs != 1) {
613       ierr = PetscPrintf(comm,"                            atomic size      = %zu [full block, bs=%D]\n", db->field[f]->atomic_size,db->field[f]->bs);CHKERRQ(ierr);
614       ierr = PetscPrintf(comm,"                            atomic size/item = %zu \n", db->field[f]->atomic_size/db->field[f]->bs);CHKERRQ(ierr);
615     } else {
616       ierr = PetscPrintf(comm,"                            atomic size      = %zu \n", db->field[f]->atomic_size);CHKERRQ(ierr);
617     }
618   }
619   ierr = PetscPrintf(comm,"  Total mem. usage                           = %1.2e (MB) (collective)\n", memory_usage_total);CHKERRQ(ierr);
620   PetscFunctionReturn(0);
621 }
622 
DMSwarmDataBucketView_SEQ(MPI_Comm comm,DMSwarmDataBucket db,const char filename[],DMSwarmDataBucketViewType type)623 PetscErrorCode DMSwarmDataBucketView_SEQ(MPI_Comm comm,DMSwarmDataBucket db,const char filename[],DMSwarmDataBucketViewType type)
624 {
625   PetscErrorCode ierr;
626 
627   PetscFunctionBegin;
628   switch (type) {
629   case DATABUCKET_VIEW_STDOUT:
630     ierr = DMSwarmDataBucketView_stdout(PETSC_COMM_SELF,db);CHKERRQ(ierr);
631     break;
632   case DATABUCKET_VIEW_ASCII:
633     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for ascii output");
634     break;
635   case DATABUCKET_VIEW_BINARY:
636     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for binary output");
637     break;
638   case DATABUCKET_VIEW_HDF5:
639     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for HDF5 output");
640     break;
641   default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer method requested");
642   }
643   PetscFunctionReturn(0);
644 }
645 
DMSwarmDataBucketView_MPI(MPI_Comm comm,DMSwarmDataBucket db,const char filename[],DMSwarmDataBucketViewType type)646 PetscErrorCode DMSwarmDataBucketView_MPI(MPI_Comm comm,DMSwarmDataBucket db,const char filename[],DMSwarmDataBucketViewType type)
647 {
648   PetscErrorCode ierr;
649 
650   PetscFunctionBegin;
651   switch (type) {
652   case DATABUCKET_VIEW_STDOUT:
653     ierr = DMSwarmDataBucketView_stdout(comm,db);CHKERRQ(ierr);
654     break;
655   case DATABUCKET_VIEW_ASCII:
656     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for ascii output");
657     break;
658   case DATABUCKET_VIEW_BINARY:
659     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for binary output");
660     break;
661   case DATABUCKET_VIEW_HDF5:
662     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for HDF5 output");
663     break;
664   default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer method requested");
665   }
666   PetscFunctionReturn(0);
667 }
668 
DMSwarmDataBucketView(MPI_Comm comm,DMSwarmDataBucket db,const char filename[],DMSwarmDataBucketViewType type)669 PetscErrorCode DMSwarmDataBucketView(MPI_Comm comm,DMSwarmDataBucket db,const char filename[],DMSwarmDataBucketViewType type)
670 {
671   PetscMPIInt size;
672   PetscErrorCode ierr;
673 
674   PetscFunctionBegin;
675   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
676   if (size == 1) {
677     ierr = DMSwarmDataBucketView_SEQ(comm,db,filename,type);CHKERRQ(ierr);
678   } else {
679     ierr = DMSwarmDataBucketView_MPI(comm,db,filename,type);CHKERRQ(ierr);
680   }
681   PetscFunctionReturn(0);
682 }
683 
DMSwarmDataBucketDuplicateFields(DMSwarmDataBucket dbA,DMSwarmDataBucket * dbB)684 PetscErrorCode DMSwarmDataBucketDuplicateFields(DMSwarmDataBucket dbA,DMSwarmDataBucket *dbB)
685 {
686   DMSwarmDataBucket db2;
687   PetscInt          f;
688   PetscErrorCode    ierr;
689 
690   PetscFunctionBegin;
691   ierr = DMSwarmDataBucketCreate(&db2);CHKERRQ(ierr);
692   /* copy contents from dbA into db2 */
693   for (f = 0; f < dbA->nfields; ++f) {
694     DMSwarmDataField field;
695     size_t    atomic_size;
696     char      *name;
697 
698     field = dbA->field[f];
699     atomic_size = field->atomic_size;
700     name        = field->name;
701     ierr = DMSwarmDataBucketRegisterField(db2,"DMSwarmDataBucketDuplicateFields",name,atomic_size,NULL);CHKERRQ(ierr);
702   }
703   ierr = DMSwarmDataBucketFinalize(db2);CHKERRQ(ierr);
704   ierr = DMSwarmDataBucketSetInitialSizes(db2,0,1000);CHKERRQ(ierr);
705   *dbB = db2;
706   PetscFunctionReturn(0);
707 }
708 
709 /*
710  Insert points from db2 into db1
711  db1 <<== db2
712  */
DMSwarmDataBucketInsertValues(DMSwarmDataBucket db1,DMSwarmDataBucket db2)713 PetscErrorCode DMSwarmDataBucketInsertValues(DMSwarmDataBucket db1,DMSwarmDataBucket db2)
714 {
715   PetscInt n_mp_points1,n_mp_points2;
716   PetscInt n_mp_points1_new,p;
717   PetscErrorCode ierr;
718 
719   PetscFunctionBegin;
720   ierr = DMSwarmDataBucketGetSizes(db1,&n_mp_points1,NULL,NULL);CHKERRQ(ierr);
721   ierr = DMSwarmDataBucketGetSizes(db2,&n_mp_points2,NULL,NULL);CHKERRQ(ierr);
722   n_mp_points1_new = n_mp_points1 + n_mp_points2;
723   ierr = DMSwarmDataBucketSetSizes(db1,n_mp_points1_new,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr);
724   for (p = 0; p < n_mp_points2; ++p) {
725     /* db1 <<== db2 */
726     ierr = DMSwarmDataBucketCopyPoint(db2,p, db1,(n_mp_points1 + p));CHKERRQ(ierr);
727   }
728   PetscFunctionReturn(0);
729 }
730 
731 /* helpers for parallel send/recv */
DMSwarmDataBucketCreatePackedArray(DMSwarmDataBucket db,size_t * bytes,void ** buf)732 PetscErrorCode DMSwarmDataBucketCreatePackedArray(DMSwarmDataBucket db,size_t *bytes,void **buf)
733 {
734   PetscInt       f;
735   size_t         sizeof_marker_contents;
736   void          *buffer;
737   PetscErrorCode ierr;
738 
739   PetscFunctionBegin;
740   sizeof_marker_contents = 0;
741   for (f = 0; f < db->nfields; ++f) {
742     DMSwarmDataField df = db->field[f];
743     sizeof_marker_contents += df->atomic_size;
744   }
745   ierr = PetscMalloc(sizeof_marker_contents, &buffer);CHKERRQ(ierr);
746   ierr = PetscMemzero(buffer, sizeof_marker_contents);CHKERRQ(ierr);
747   if (bytes) {*bytes = sizeof_marker_contents;}
748   if (buf)   {*buf   = buffer;}
749   PetscFunctionReturn(0);
750 }
751 
DMSwarmDataBucketDestroyPackedArray(DMSwarmDataBucket db,void ** buf)752 PetscErrorCode DMSwarmDataBucketDestroyPackedArray(DMSwarmDataBucket db,void **buf)
753 {
754   PetscErrorCode ierr;
755 
756   PetscFunctionBegin;
757   if (buf) {
758     ierr = PetscFree(*buf);CHKERRQ(ierr);
759     *buf = NULL;
760   }
761   PetscFunctionReturn(0);
762 }
763 
DMSwarmDataBucketFillPackedArray(DMSwarmDataBucket db,const PetscInt index,void * buf)764 PetscErrorCode DMSwarmDataBucketFillPackedArray(DMSwarmDataBucket db,const PetscInt index,void *buf)
765 {
766   PetscInt       f;
767   void          *data, *data_p;
768   size_t         asize, offset;
769   PetscErrorCode ierr;
770 
771   PetscFunctionBegin;
772   offset = 0;
773   for (f = 0; f < db->nfields; ++f) {
774     DMSwarmDataField df = db->field[f];
775 
776     asize = df->atomic_size;
777     data = (void*)( df->data);
778     data_p = (void*)( (char*)data + index*asize);
779     ierr = PetscMemcpy((void*)((char*)buf + offset), data_p, asize);CHKERRQ(ierr);
780     offset = offset + asize;
781   }
782   PetscFunctionReturn(0);
783 }
784 
DMSwarmDataBucketInsertPackedArray(DMSwarmDataBucket db,const PetscInt idx,void * data)785 PetscErrorCode DMSwarmDataBucketInsertPackedArray(DMSwarmDataBucket db,const PetscInt idx,void *data)
786 {
787   PetscInt f;
788   void *data_p;
789   size_t offset;
790   PetscErrorCode ierr;
791 
792   PetscFunctionBegin;
793   offset = 0;
794   for (f = 0; f < db->nfields; ++f) {
795     DMSwarmDataField df = db->field[f];
796 
797     data_p = (void*)( (char*)data + offset);
798     ierr = DMSwarmDataFieldInsertPoint(df, idx, (void*)data_p);CHKERRQ(ierr);
799     offset = offset + df->atomic_size;
800   }
801   PetscFunctionReturn(0);
802 }
803