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