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