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