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