xref: /petsc/src/dm/impls/swarm/data_bucket.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
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   PetscCheck(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 DMSwarmDataBucketGetDMSwarmDataFieldIdByName(DMSwarmDataBucket db,const char name[],PetscInt *idx)
160 {
161   PetscBool found;
162 
163   PetscFunctionBegin;
164   *idx = -1;
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   PetscFunctionReturn(0);
169 }
170 
171 PetscErrorCode DMSwarmDataBucketGetDMSwarmDataFieldByName(DMSwarmDataBucket db,const char name[],DMSwarmDataField *gfield)
172 {
173   PetscInt       idx;
174   PetscBool      found;
175 
176   PetscFunctionBegin;
177   PetscCall(DMSwarmDataFieldStringInList(name,db->nfields,(const DMSwarmDataField*)db->field,&found));
178   PetscCheck(found,PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot find DMSwarmDataField with name %s",name);
179   PetscCall(DMSwarmDataFieldStringFindInList(name,db->nfields,(const DMSwarmDataField*)db->field,&idx));
180   *gfield = db->field[idx];
181   PetscFunctionReturn(0);
182 }
183 
184 PetscErrorCode DMSwarmDataBucketQueryDMSwarmDataFieldByName(DMSwarmDataBucket db,const char name[],PetscBool *found)
185 {
186   PetscFunctionBegin;
187   *found = PETSC_FALSE;
188   PetscCall(DMSwarmDataFieldStringInList(name,db->nfields,(const DMSwarmDataField*)db->field,found));
189   PetscFunctionReturn(0);
190 }
191 
192 PetscErrorCode DMSwarmDataBucketFinalize(DMSwarmDataBucket db)
193 {
194   PetscFunctionBegin;
195   db->finalised = PETSC_TRUE;
196   PetscFunctionReturn(0);
197 }
198 
199 PetscErrorCode DMSwarmDataFieldGetNumEntries(DMSwarmDataField df,PetscInt *sum)
200 {
201   PetscFunctionBegin;
202   *sum = df->L;
203   PetscFunctionReturn(0);
204 }
205 
206 PetscErrorCode DMSwarmDataFieldSetBlockSize(DMSwarmDataField df,PetscInt blocksize)
207 {
208   PetscFunctionBegin;
209   df->bs = blocksize;
210   PetscFunctionReturn(0);
211 }
212 
213 PetscErrorCode DMSwarmDataFieldSetSize(DMSwarmDataField df,const PetscInt new_L)
214 {
215   PetscFunctionBegin;
216   PetscCheck(new_L >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot set size of DMSwarmDataField to be < 0");
217   if (new_L == df->L) PetscFunctionReturn(0);
218   if (new_L > df->L) {
219     PetscCall(PetscRealloc(df->atomic_size * (new_L), &df->data));
220     /* init new contents */
221     PetscCall(PetscMemzero(( ((char*)df->data)+df->L*df->atomic_size), (new_L-df->L)*df->atomic_size));
222   } else {
223     /* reallocate pointer list, add +1 in case new_L = 0 */
224     PetscCall(PetscRealloc(df->atomic_size * (new_L+1), &df->data));
225   }
226   df->L = new_L;
227   PetscFunctionReturn(0);
228 }
229 
230 PetscErrorCode DMSwarmDataFieldZeroBlock(DMSwarmDataField df,const PetscInt start,const PetscInt end)
231 {
232   PetscFunctionBegin;
233   PetscCheck(start <= end,PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if start(%" PetscInt_FMT ") > end(%" PetscInt_FMT ")",start,end);
234   PetscCheck(start >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if start(%" PetscInt_FMT ") < 0",start);
235   PetscCheck(end <= df->L,PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if end(%" PetscInt_FMT ") >= array size(%" PetscInt_FMT ")",end,df->L);
236   PetscCall(PetscMemzero((((char*)df->data)+start*df->atomic_size), (end-start)*df->atomic_size));
237   PetscFunctionReturn(0);
238 }
239 
240 /*
241  A negative buffer value will simply be ignored and the old buffer value will be used.
242  */
243 PetscErrorCode DMSwarmDataBucketSetSizes(DMSwarmDataBucket db,const PetscInt L,const PetscInt buffer)
244 {
245   PetscInt       current_allocated,new_used,new_unused,new_buffer,new_allocated,f;
246   PetscBool      any_active_fields;
247 
248   PetscFunctionBegin;
249   PetscCheck(db->finalised != PETSC_FALSE,PETSC_COMM_SELF,PETSC_ERR_USER,"You must call DMSwarmDataBucketFinalize() before DMSwarmDataBucketSetSizes()");
250   PetscCall(DMSwarmDataBucketQueryForActiveFields(db,&any_active_fields));
251   PetscCheck(!any_active_fields,PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot safely re-size as at least one DMSwarmDataField is currently being accessed");
252 
253   current_allocated = db->allocated;
254   new_used   = L;
255   new_unused = current_allocated - new_used;
256   new_buffer = db->buffer;
257   if (buffer >= 0) { /* update the buffer value */
258     new_buffer = buffer;
259   }
260   new_allocated = new_used + new_buffer;
261   /* action */
262   if (new_allocated > current_allocated) {
263     /* increase size to new_used + new_buffer */
264     for (f=0; f<db->nfields; f++) {
265       PetscCall(DMSwarmDataFieldSetSize(db->field[f], new_allocated));
266     }
267     db->L         = new_used;
268     db->buffer    = new_buffer;
269     db->allocated = new_used + new_buffer;
270   } else {
271     if (new_unused > 2 * new_buffer) {
272       /* shrink array to new_used + new_buffer */
273       for (f = 0; f < db->nfields; ++f) {
274         PetscCall(DMSwarmDataFieldSetSize(db->field[f], new_allocated));
275       }
276       db->L         = new_used;
277       db->buffer    = new_buffer;
278       db->allocated = new_used + new_buffer;
279     } else {
280       db->L      = new_used;
281       db->buffer = new_buffer;
282     }
283   }
284   /* zero all entries from db->L to db->allocated */
285   for (f = 0; f < db->nfields; ++f) {
286     DMSwarmDataField field = db->field[f];
287     PetscCall(DMSwarmDataFieldZeroBlock(field, db->L,db->allocated));
288   }
289   PetscFunctionReturn(0);
290 }
291 
292 PetscErrorCode DMSwarmDataBucketSetInitialSizes(DMSwarmDataBucket db,const PetscInt L,const PetscInt buffer)
293 {
294   PetscInt       f;
295 
296   PetscFunctionBegin;
297   PetscCall(DMSwarmDataBucketSetSizes(db,L,buffer));
298   for (f = 0; f < db->nfields; ++f) {
299     DMSwarmDataField field = db->field[f];
300     PetscCall(DMSwarmDataFieldZeroBlock(field,0,db->allocated));
301   }
302   PetscFunctionReturn(0);
303 }
304 
305 PetscErrorCode DMSwarmDataBucketGetSizes(DMSwarmDataBucket db,PetscInt *L,PetscInt *buffer,PetscInt *allocated)
306 {
307   PetscFunctionBegin;
308   if (L) {*L = db->L;}
309   if (buffer) {*buffer = db->buffer;}
310   if (allocated) {*allocated = db->allocated;}
311   PetscFunctionReturn(0);
312 }
313 
314 PetscErrorCode DMSwarmDataBucketGetGlobalSizes(MPI_Comm comm,DMSwarmDataBucket db,PetscInt *L,PetscInt *buffer,PetscInt *allocated)
315 {
316   PetscFunctionBegin;
317   if (L) PetscCallMPI(MPI_Allreduce(&db->L,L,1,MPIU_INT,MPI_SUM,comm));
318   if (buffer) PetscCallMPI(MPI_Allreduce(&db->buffer,buffer,1,MPIU_INT,MPI_SUM,comm));
319   if (allocated) PetscCallMPI(MPI_Allreduce(&db->allocated,allocated,1,MPIU_INT,MPI_SUM,comm));
320   PetscFunctionReturn(0);
321 }
322 
323 PetscErrorCode DMSwarmDataBucketGetDMSwarmDataFields(DMSwarmDataBucket db,PetscInt *L,DMSwarmDataField *fields[])
324 {
325   PetscFunctionBegin;
326   if (L)      {*L      = db->nfields;}
327   if (fields) {*fields = db->field;}
328   PetscFunctionReturn(0);
329 }
330 
331 PetscErrorCode DMSwarmDataFieldGetAccess(const DMSwarmDataField gfield)
332 {
333   PetscFunctionBegin;
334   PetscCheck(!gfield->active,PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is already active. You must call DMSwarmDataFieldRestoreAccess()",gfield->name);
335   gfield->active = PETSC_TRUE;
336   PetscFunctionReturn(0);
337 }
338 
339 PetscErrorCode DMSwarmDataFieldAccessPoint(const DMSwarmDataField gfield,const PetscInt pid,void **ctx_p)
340 {
341   PetscFunctionBegin;
342   *ctx_p = NULL;
343 #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
344   /* debug mode */
345   /* check point is valid */
346   PetscCheck(pid >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
347   PetscCheck(pid < gfield->L,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %" PetscInt_FMT,gfield->L);
348   PetscCheck(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);
349 #endif
350   *ctx_p = DMSWARM_DATAFIELD_point_access(gfield->data,pid,gfield->atomic_size);
351   PetscFunctionReturn(0);
352 }
353 
354 PetscErrorCode DMSwarmDataFieldAccessPointOffset(const DMSwarmDataField gfield,const size_t offset,const PetscInt pid,void **ctx_p)
355 {
356   PetscFunctionBegin;
357 #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
358   /* debug mode */
359   /* check point is valid */
360   /* PetscCheck(offset >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be >= 0");*/
361   /* Note compiler realizes this can never happen with an unsigned PetscInt */
362   PetscCheck(offset < gfield->atomic_size,PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be < %zu",gfield->atomic_size);
363   /* check point is valid */
364   PetscCheck(pid >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
365   PetscCheck(pid < gfield->L,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %" PetscInt_FMT,gfield->L);
366   PetscCheck(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);
367 #endif
368   *ctx_p = DMSWARM_DATAFIELD_point_access_offset(gfield->data,pid,gfield->atomic_size,offset);
369   PetscFunctionReturn(0);
370 }
371 
372 PetscErrorCode DMSwarmDataFieldRestoreAccess(DMSwarmDataField gfield)
373 {
374   PetscFunctionBegin;
375   PetscCheck(gfield->active != PETSC_FALSE,PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is not active. You must call DMSwarmDataFieldGetAccess()", gfield->name);
376   gfield->active = PETSC_FALSE;
377   PetscFunctionReturn(0);
378 }
379 
380 PetscErrorCode DMSwarmDataFieldVerifyAccess(const DMSwarmDataField gfield,const size_t size)
381 {
382   PetscFunctionBegin;
383 #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
384   PetscCheck(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);
385 #endif
386   PetscFunctionReturn(0);
387 }
388 
389 PetscErrorCode DMSwarmDataFieldGetAtomicSize(const DMSwarmDataField gfield,size_t *size)
390 {
391   PetscFunctionBegin;
392   if (size) {*size = gfield->atomic_size;}
393   PetscFunctionReturn(0);
394 }
395 
396 PetscErrorCode DMSwarmDataFieldGetEntries(const DMSwarmDataField gfield,void **data)
397 {
398   PetscFunctionBegin;
399   if (data) {*data = gfield->data;}
400   PetscFunctionReturn(0);
401 }
402 
403 PetscErrorCode DMSwarmDataFieldRestoreEntries(const DMSwarmDataField gfield,void **data)
404 {
405   PetscFunctionBegin;
406   if (data) {*data = NULL;}
407   PetscFunctionReturn(0);
408 }
409 
410 /* y = x */
411 PetscErrorCode DMSwarmDataBucketCopyPoint(const DMSwarmDataBucket xb,const PetscInt pid_x,const DMSwarmDataBucket yb,const PetscInt pid_y)
412 {
413   PetscInt       f;
414 
415   PetscFunctionBegin;
416   for (f = 0; f < xb->nfields; ++f) {
417     void *dest;
418     void *src;
419 
420     PetscCall(DMSwarmDataFieldGetAccess(xb->field[f]));
421     if (xb != yb) PetscCall(DMSwarmDataFieldGetAccess( yb->field[f]));
422     PetscCall(DMSwarmDataFieldAccessPoint(xb->field[f],pid_x, &src));
423     PetscCall(DMSwarmDataFieldAccessPoint(yb->field[f],pid_y, &dest));
424     PetscCall(PetscMemcpy(dest, src, xb->field[f]->atomic_size));
425     PetscCall(DMSwarmDataFieldRestoreAccess(xb->field[f]));
426     if (xb != yb) PetscCall(DMSwarmDataFieldRestoreAccess(yb->field[f]));
427   }
428   PetscFunctionReturn(0);
429 }
430 
431 PetscErrorCode DMSwarmDataBucketCreateFromSubset(DMSwarmDataBucket DBIn,const PetscInt N,const PetscInt list[],DMSwarmDataBucket *DB)
432 {
433   PetscInt         nfields;
434   DMSwarmDataField *fields;
435   PetscInt         f,L,buffer,allocated,p;
436 
437   PetscFunctionBegin;
438   PetscCall(DMSwarmDataBucketCreate(DB));
439   /* copy contents of DBIn */
440   PetscCall(DMSwarmDataBucketGetDMSwarmDataFields(DBIn,&nfields,&fields));
441   PetscCall(DMSwarmDataBucketGetSizes(DBIn,&L,&buffer,&allocated));
442   for (f = 0; f < nfields; ++f) {
443     PetscCall(DMSwarmDataBucketRegisterField(*DB,"DMSwarmDataBucketCreateFromSubset",fields[f]->name,fields[f]->atomic_size,NULL));
444   }
445   PetscCall(DMSwarmDataBucketFinalize(*DB));
446   PetscCall(DMSwarmDataBucketSetSizes(*DB,L,buffer));
447   /* now copy the desired guys from DBIn => DB */
448   for (p = 0; p < N; ++p) {
449     PetscCall(DMSwarmDataBucketCopyPoint(DBIn,list[p], *DB,list[p]));
450   }
451   PetscFunctionReturn(0);
452 }
453 
454 /* insert into an exisitng location */
455 PetscErrorCode DMSwarmDataFieldInsertPoint(const DMSwarmDataField field,const PetscInt index,const void *ctx)
456 {
457   PetscFunctionBegin;
458 #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
459   /* check point is valid */
460   PetscCheck(index >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
461   PetscCheck(index < field->L,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %" PetscInt_FMT,field->L);
462 #endif
463   PetscCall(PetscMemcpy(DMSWARM_DATAFIELD_point_access(field->data,index,field->atomic_size), ctx, field->atomic_size));
464   PetscFunctionReturn(0);
465 }
466 
467 /* remove data at index - replace with last point */
468 PetscErrorCode DMSwarmDataBucketRemovePointAtIndex(const DMSwarmDataBucket db,const PetscInt index)
469 {
470   PetscInt       f;
471   PetscBool      any_active_fields;
472 
473   PetscFunctionBegin;
474 #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
475   /* check point is valid */
476   PetscCheck(index >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
477   PetscCheck(index < db->allocated,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %" PetscInt_FMT,db->L+db->buffer);
478 #endif
479   PetscCall(DMSwarmDataBucketQueryForActiveFields(db,&any_active_fields));
480   PetscCheck(!any_active_fields,PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot safely remove point as at least one DMSwarmDataField is currently being accessed");
481   if (index >= db->L) { /* this point is not in the list - no need to error, but I will anyway */
482     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"You should not be trying to remove point at index=%" PetscInt_FMT " since it's < db->L = %" PetscInt_FMT, index, db->L);
483   }
484   if (index != db->L-1) { /* not last point in list */
485     for (f = 0; f < db->nfields; ++f) {
486       DMSwarmDataField field = db->field[f];
487 
488       /* copy then remove */
489       PetscCall(DMSwarmDataFieldCopyPoint(db->L-1, field, index, field));
490       /* DMSwarmDataFieldZeroPoint(field,index); */
491     }
492   }
493   /* decrement size */
494   /* this will zero out an crap at the end of the list */
495   PetscCall(DMSwarmDataBucketRemovePoint(db));
496   PetscFunctionReturn(0);
497 }
498 
499 /* copy x into y */
500 PetscErrorCode DMSwarmDataFieldCopyPoint(const PetscInt pid_x,const DMSwarmDataField field_x,const PetscInt pid_y,const DMSwarmDataField field_y)
501 {
502   PetscFunctionBegin;
503 #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
504   /* check point is valid */
505   PetscCheck(pid_x >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"(IN) index must be >= 0");
506   PetscCheck(pid_x < field_x->L,PETSC_COMM_SELF,PETSC_ERR_USER,"(IN) index must be < %" PetscInt_FMT,field_x->L);
507   PetscCheck(pid_y >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"(OUT) index must be >= 0");
508   PetscCheck(pid_y < field_y->L,PETSC_COMM_SELF,PETSC_ERR_USER,"(OUT) index must be < %" PetscInt_FMT,field_y->L);
509   PetscCheck(field_y->atomic_size == field_x->atomic_size,PETSC_COMM_SELF,PETSC_ERR_USER,"atomic size must match");
510 #endif
511   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));
512   PetscFunctionReturn(0);
513 }
514 
515 /* zero only the datafield at this point */
516 PetscErrorCode DMSwarmDataFieldZeroPoint(const DMSwarmDataField field,const PetscInt index)
517 {
518   PetscFunctionBegin;
519 #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
520   /* check point is valid */
521   PetscCheck(index >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
522   PetscCheck(index < field->L,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %" PetscInt_FMT,field->L);
523 #endif
524   PetscCall(PetscMemzero(DMSWARM_DATAFIELD_point_access(field->data,index,field->atomic_size), field->atomic_size));
525   PetscFunctionReturn(0);
526 }
527 
528 /* zero ALL data for this point */
529 PetscErrorCode DMSwarmDataBucketZeroPoint(const DMSwarmDataBucket db,const PetscInt index)
530 {
531   PetscInt       f;
532 
533   PetscFunctionBegin;
534   /* check point is valid */
535   PetscCheck(index >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
536   PetscCheck(index < db->allocated,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %" PetscInt_FMT,db->allocated);
537   for (f = 0; f < db->nfields; ++f) {
538     DMSwarmDataField field = db->field[f];
539     PetscCall(DMSwarmDataFieldZeroPoint(field,index));
540   }
541   PetscFunctionReturn(0);
542 }
543 
544 /* increment */
545 PetscErrorCode DMSwarmDataBucketAddPoint(DMSwarmDataBucket db)
546 {
547   PetscFunctionBegin;
548   PetscCall(DMSwarmDataBucketSetSizes(db,db->L+1,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
549   PetscFunctionReturn(0);
550 }
551 
552 /* decrement */
553 PetscErrorCode DMSwarmDataBucketRemovePoint(DMSwarmDataBucket db)
554 {
555   PetscFunctionBegin;
556   PetscCall(DMSwarmDataBucketSetSizes(db,db->L-1,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
557   PetscFunctionReturn(0);
558 }
559 
560 /*  Should be redone to user PetscViewer */
561 PetscErrorCode DMSwarmDataBucketView_stdout(MPI_Comm comm,DMSwarmDataBucket db)
562 {
563   PetscInt       f;
564   double         memory_usage_total,memory_usage_total_local = 0.0;
565 
566   PetscFunctionBegin;
567   PetscCall(PetscPrintf(comm,"DMSwarmDataBucketView: \n"));
568   PetscCall(PetscPrintf(comm,"  L                  = %" PetscInt_FMT " \n", db->L));
569   PetscCall(PetscPrintf(comm,"  buffer             = %" PetscInt_FMT " \n", db->buffer));
570   PetscCall(PetscPrintf(comm,"  allocated          = %" PetscInt_FMT " \n", db->allocated));
571   PetscCall(PetscPrintf(comm,"  nfields registered = %" PetscInt_FMT " \n", db->nfields));
572 
573   for (f = 0; f < db->nfields; ++f) {
574     double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
575     memory_usage_total_local += memory_usage_f;
576   }
577   PetscCallMPI(MPI_Allreduce(&memory_usage_total_local,&memory_usage_total,1,MPI_DOUBLE,MPI_SUM,comm));
578 
579   for (f = 0; f < db->nfields; ++f) {
580     double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
581     PetscCall(PetscPrintf(comm,"    [%3" PetscInt_FMT "] %15s : Mem. usage       = %1.2e (MB) [rank0]\n", f, db->field[f]->name, memory_usage_f));
582     PetscCall(PetscPrintf(comm,"                            blocksize        = %" PetscInt_FMT " \n", db->field[f]->bs));
583     if (db->field[f]->bs != 1) {
584       PetscCall(PetscPrintf(comm,"                            atomic size      = %zu [full block, bs=%" PetscInt_FMT "]\n", db->field[f]->atomic_size,db->field[f]->bs));
585       PetscCall(PetscPrintf(comm,"                            atomic size/item = %zu \n", (size_t)(db->field[f]->atomic_size/db->field[f]->bs)));
586     } else {
587       PetscCall(PetscPrintf(comm,"                            atomic size      = %zu \n", db->field[f]->atomic_size));
588     }
589   }
590   PetscCall(PetscPrintf(comm,"  Total mem. usage                           = %1.2e (MB) (collective)\n", memory_usage_total));
591   PetscFunctionReturn(0);
592 }
593 
594 PetscErrorCode DMSwarmDataBucketView_Seq(MPI_Comm comm,DMSwarmDataBucket db,const char filename[],DMSwarmDataBucketViewType type)
595 {
596   PetscFunctionBegin;
597   switch (type) {
598   case DATABUCKET_VIEW_STDOUT:
599     PetscCall(DMSwarmDataBucketView_stdout(PETSC_COMM_SELF,db));
600     break;
601   case DATABUCKET_VIEW_ASCII:
602     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for ascii output");
603   case DATABUCKET_VIEW_BINARY:
604     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for binary output");
605   case DATABUCKET_VIEW_HDF5:
606     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for HDF5 output");
607   default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer method requested");
608   }
609   PetscFunctionReturn(0);
610 }
611 
612 PetscErrorCode DMSwarmDataBucketView_MPI(MPI_Comm comm,DMSwarmDataBucket db,const char filename[],DMSwarmDataBucketViewType type)
613 {
614   PetscFunctionBegin;
615   switch (type) {
616   case DATABUCKET_VIEW_STDOUT:
617     PetscCall(DMSwarmDataBucketView_stdout(comm,db));
618     break;
619   case DATABUCKET_VIEW_ASCII:
620     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for ascii output");
621   case DATABUCKET_VIEW_BINARY:
622     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for binary output");
623   case DATABUCKET_VIEW_HDF5:
624     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for HDF5 output");
625   default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer method requested");
626   }
627   PetscFunctionReturn(0);
628 }
629 
630 PetscErrorCode DMSwarmDataBucketView(MPI_Comm comm,DMSwarmDataBucket db,const char filename[],DMSwarmDataBucketViewType type)
631 {
632   PetscMPIInt    size;
633 
634   PetscFunctionBegin;
635   PetscCallMPI(MPI_Comm_size(comm,&size));
636   if (size == 1) {
637     PetscCall(DMSwarmDataBucketView_Seq(comm,db,filename,type));
638   } else {
639     PetscCall(DMSwarmDataBucketView_MPI(comm,db,filename,type));
640   }
641   PetscFunctionReturn(0);
642 }
643 
644 PetscErrorCode DMSwarmDataBucketDuplicateFields(DMSwarmDataBucket dbA,DMSwarmDataBucket *dbB)
645 {
646   DMSwarmDataBucket db2;
647   PetscInt          f;
648 
649   PetscFunctionBegin;
650   PetscCall(DMSwarmDataBucketCreate(&db2));
651   /* copy contents from dbA into db2 */
652   for (f = 0; f < dbA->nfields; ++f) {
653     DMSwarmDataField field;
654     size_t           atomic_size;
655     char             *name;
656 
657     field = dbA->field[f];
658     atomic_size = field->atomic_size;
659     name        = field->name;
660     PetscCall(DMSwarmDataBucketRegisterField(db2,"DMSwarmDataBucketDuplicateFields",name,atomic_size,NULL));
661   }
662   PetscCall(DMSwarmDataBucketFinalize(db2));
663   PetscCall(DMSwarmDataBucketSetInitialSizes(db2,0,1000));
664   *dbB = db2;
665   PetscFunctionReturn(0);
666 }
667 
668 /*
669  Insert points from db2 into db1
670  db1 <<== db2
671  */
672 PetscErrorCode DMSwarmDataBucketInsertValues(DMSwarmDataBucket db1,DMSwarmDataBucket db2)
673 {
674   PetscInt       n_mp_points1,n_mp_points2;
675   PetscInt       n_mp_points1_new,p;
676 
677   PetscFunctionBegin;
678   PetscCall(DMSwarmDataBucketGetSizes(db1,&n_mp_points1,NULL,NULL));
679   PetscCall(DMSwarmDataBucketGetSizes(db2,&n_mp_points2,NULL,NULL));
680   n_mp_points1_new = n_mp_points1 + n_mp_points2;
681   PetscCall(DMSwarmDataBucketSetSizes(db1,n_mp_points1_new,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
682   for (p = 0; p < n_mp_points2; ++p) {
683     /* db1 <<== db2 */
684     PetscCall(DMSwarmDataBucketCopyPoint(db2,p, db1,(n_mp_points1 + p)));
685   }
686   PetscFunctionReturn(0);
687 }
688 
689 /* helpers for parallel send/recv */
690 PetscErrorCode DMSwarmDataBucketCreatePackedArray(DMSwarmDataBucket db,size_t *bytes,void **buf)
691 {
692   PetscInt       f;
693   size_t         sizeof_marker_contents;
694   void          *buffer;
695 
696   PetscFunctionBegin;
697   sizeof_marker_contents = 0;
698   for (f = 0; f < db->nfields; ++f) {
699     DMSwarmDataField df = db->field[f];
700     sizeof_marker_contents += df->atomic_size;
701   }
702   PetscCall(PetscMalloc(sizeof_marker_contents, &buffer));
703   PetscCall(PetscMemzero(buffer, sizeof_marker_contents));
704   if (bytes) {*bytes = sizeof_marker_contents;}
705   if (buf)   {*buf   = buffer;}
706   PetscFunctionReturn(0);
707 }
708 
709 PetscErrorCode DMSwarmDataBucketDestroyPackedArray(DMSwarmDataBucket db,void **buf)
710 {
711   PetscFunctionBegin;
712   if (buf) {
713     PetscCall(PetscFree(*buf));
714     *buf = NULL;
715   }
716   PetscFunctionReturn(0);
717 }
718 
719 PetscErrorCode DMSwarmDataBucketFillPackedArray(DMSwarmDataBucket db,const PetscInt index,void *buf)
720 {
721   PetscInt       f;
722   void          *data, *data_p;
723   size_t         asize, offset;
724 
725   PetscFunctionBegin;
726   offset = 0;
727   for (f = 0; f < db->nfields; ++f) {
728     DMSwarmDataField df = db->field[f];
729 
730     asize = df->atomic_size;
731     data = (void*)( df->data);
732     data_p = (void*)( (char*)data + index*asize);
733     PetscCall(PetscMemcpy((void*)((char*)buf + offset), data_p, asize));
734     offset = offset + asize;
735   }
736   PetscFunctionReturn(0);
737 }
738 
739 PetscErrorCode DMSwarmDataBucketInsertPackedArray(DMSwarmDataBucket db,const PetscInt idx,void *data)
740 {
741   PetscInt       f;
742   void           *data_p;
743   size_t         offset;
744 
745   PetscFunctionBegin;
746   offset = 0;
747   for (f = 0; f < db->nfields; ++f) {
748     DMSwarmDataField df = db->field[f];
749 
750     data_p = (void*)( (char*)data + offset);
751     PetscCall(DMSwarmDataFieldInsertPoint(df, idx, (void*)data_p));
752     offset = offset + df->atomic_size;
753   }
754   PetscFunctionReturn(0);
755 }
756