xref: /petsc/src/dm/impls/swarm/data_bucket.c (revision fc2c8c6ab10fa0cecf029e6ef74b247b4e841624)
1 
2 #include "data_bucket.h"
3 
4 /* string helpers */
5 #undef __FUNCT__
6 #define __FUNCT__ "StringInList"
7 PetscErrorCode StringInList(const char name[],const PetscInt N,const DataField gfield[],PetscBool *val)
8 {
9   PetscInt       i;
10   PetscErrorCode ierr;
11 
12   PetscFunctionBegin;
13   *val = PETSC_FALSE;
14   for (i = 0; i < N; ++i) {
15     PetscBool flg;
16     ierr = PetscStrcmp(name, gfield[i]->name, &flg);CHKERRQ(ierr);
17     if (flg) {
18       *val = PETSC_TRUE;
19       PetscFunctionReturn(0);
20     }
21   }
22   PetscFunctionReturn(0);
23 }
24 
25 #undef __FUNCT__
26 #define __FUNCT__ "StringFindInList"
27 PetscErrorCode StringFindInList(const char name[],const PetscInt N,const DataField gfield[],PetscInt *index)
28 {
29   PetscInt       i;
30   PetscErrorCode ierr;
31 
32   PetscFunctionBegin;
33   *index = -1;
34   for (i = 0; i < N; ++i) {
35     PetscBool flg;
36     ierr = PetscStrcmp(name, gfield[i]->name, &flg);CHKERRQ(ierr);
37     if (flg) {
38       *index = i;
39       PetscFunctionReturn(0);
40     }
41   }
42   PetscFunctionReturn(0);
43 }
44 
45 #undef __FUNCT__
46 #define __FUNCT__ "DataFieldCreate"
47 PetscErrorCode DataFieldCreate(const char registeration_function[],const char name[],const size_t size,const PetscInt L,DataField *DF)
48 {
49   DataField      df;
50   PetscErrorCode ierr;
51 
52   PetscFunctionBegin;
53   ierr = PetscMalloc(sizeof(struct _p_DataField), &df);CHKERRQ(ierr);
54   ierr = PetscMemzero(df, sizeof(struct _p_DataField));CHKERRQ(ierr);
55   ierr = PetscStrallocpy(registeration_function, &df->registeration_function);CHKERRQ(ierr);
56   ierr = PetscStrallocpy(name, &df->name);CHKERRQ(ierr);
57   df->atomic_size = size;
58   df->L  = L;
59   df->bs = 1;
60   /* allocate something so we don't have to reallocate */
61   ierr = PetscMalloc(size * L, &df->data);CHKERRQ(ierr);
62   ierr = PetscMemzero(df->data, size * L);CHKERRQ(ierr);
63   *DF = df;
64   PetscFunctionReturn(0);
65 }
66 
67 #undef __FUNCT__
68 #define __FUNCT__ "DataFieldDestroy"
69 PetscErrorCode DataFieldDestroy(DataField *DF)
70 {
71   DataField      df = *DF;
72   PetscErrorCode ierr;
73 
74   PetscFunctionBegin;
75   ierr = PetscFree(df->registeration_function);CHKERRQ(ierr);
76   ierr = PetscFree(df->name);CHKERRQ(ierr);
77   ierr = PetscFree(df->data);CHKERRQ(ierr);
78   ierr = PetscFree(df);CHKERRQ(ierr);
79   *DF  = NULL;
80   PetscFunctionReturn(0);
81 }
82 
83 /* data bucket */
84 #undef __FUNCT__
85 #define __FUNCT__ "DataBucketCreate"
86 PetscErrorCode DataBucketCreate(DataBucket *DB)
87 {
88   DataBucket     db;
89   PetscErrorCode ierr;
90 
91   PetscFunctionBegin;
92   ierr = PetscMalloc(sizeof(struct _p_DataBucket), &db);CHKERRQ(ierr);
93   ierr = PetscMemzero(db, sizeof(struct _p_DataBucket));CHKERRQ(ierr);
94 
95   db->finalised = PETSC_FALSE;
96   /* create empty spaces for fields */
97   db->L         = -1;
98   db->buffer    = 1;
99   db->allocated = 1;
100   db->nfields   = 0;
101   ierr = PetscMalloc1(1, &db->field);CHKERRQ(ierr);
102   *DB  = db;
103   PetscFunctionReturn(0);
104 }
105 
106 #undef __FUNCT__
107 #define __FUNCT__ "DataBucketDestroy"
108 PetscErrorCode DataBucketDestroy(DataBucket *DB)
109 {
110   DataBucket     db = *DB;
111   PetscInt       f;
112   PetscErrorCode ierr;
113 
114   PetscFunctionBegin;
115   /* release fields */
116   for (f = 0; f < db->nfields; ++f) {
117     ierr = DataFieldDestroy(&db->field[f]);CHKERRQ(ierr);
118   }
119   /* this will catch the initially allocated objects in the event that no fields are registered */
120   if (db->field != NULL) {
121     ierr = PetscFree(db->field);CHKERRQ(ierr);
122   }
123   ierr = PetscFree(db);CHKERRQ(ierr);
124   *DB = NULL;
125   PetscFunctionReturn(0);
126 }
127 
128 #undef __FUNCT__
129 #define __FUNCT__ "DataBucketQueryForActiveFields"
130 PetscErrorCode DataBucketQueryForActiveFields(DataBucket db,PetscBool *any_active_fields)
131 {
132   PetscInt f;
133 
134   PetscFunctionBegin;
135   *any_active_fields = PETSC_FALSE;
136   for (f = 0; f < db->nfields; ++f) {
137     if (db->field[f]->active) {
138       *any_active_fields = PETSC_TRUE;
139       PetscFunctionReturn(0);
140     }
141   }
142   PetscFunctionReturn(0);
143 }
144 
145 #undef __FUNCT__
146 #define __FUNCT__ "DataBucketRegisterField"
147 PetscErrorCode DataBucketRegisterField(
148                               DataBucket db,
149                               const char registeration_function[],
150                               const char field_name[],
151                               size_t atomic_size, DataField *_gfield)
152 {
153   PetscBool val;
154   DataField *field,fp;
155   PetscErrorCode ierr;
156 
157   PetscFunctionBegin;
158 	/* check we haven't finalised the registration of fields */
159 	/*
160    if(db->finalised==PETSC_TRUE) {
161    printf("ERROR: DataBucketFinalize() has been called. Cannot register more fields\n");
162    ERROR();
163    }
164    */
165   /* check for repeated name */
166   ierr = StringInList(field_name, db->nfields, (const DataField*) db->field, &val);CHKERRQ(ierr);
167   if (val == PETSC_TRUE) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field %s already exists. Cannot add same field twice",field_name);
168   /* create new space for data */
169   field = (DataField *) realloc(db->field, sizeof(DataField)*(db->nfields+1));
170   db->field = field;
171   /* add field */
172   ierr = DataFieldCreate(registeration_function, field_name, atomic_size, db->allocated, &fp);CHKERRQ(ierr);
173   db->field[db->nfields] = fp;
174   db->nfields++;
175   if (_gfield != NULL) {
176     *_gfield = fp;
177   }
178   PetscFunctionReturn(0);
179 }
180 
181 /*
182  #define DataBucketRegisterField(db,name,size,k) {\
183  char *location;\
184  asprintf(&location,"Registered by %s() at line %d within file %s", __FUNCTION__, __LINE__, __FILE__);\
185  _DataBucketRegisterField( (db), location, (name), (size), (k) );\
186  ierr = PetscFree(location);\
187  }
188  */
189 
190 #undef __FUNCT__
191 #define __FUNCT__ "DataBucketGetDataFieldByName"
192 PetscErrorCode DataBucketGetDataFieldByName(DataBucket db,const char name[],DataField *gfield)
193 {
194   PetscInt       idx;
195   PetscBool      found;
196   PetscErrorCode ierr;
197 
198   PetscFunctionBegin;
199   ierr = StringInList(name,db->nfields,(const DataField*)db->field,&found);CHKERRQ(ierr);
200   if (!found) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot find DataField with name %s",name);
201   ierr = StringFindInList(name,db->nfields,(const DataField*)db->field,&idx);CHKERRQ(ierr);
202   *gfield = db->field[idx];
203   PetscFunctionReturn(0);
204 }
205 
206 #undef __FUNCT__
207 #define __FUNCT__ "DataBucketQueryDataFieldByName"
208 PetscErrorCode DataBucketQueryDataFieldByName(DataBucket db,const char name[],PetscBool *found)
209 {
210   PetscErrorCode ierr;
211 
212   PetscFunctionBegin;
213   *found = PETSC_FALSE;
214   ierr = StringInList(name,db->nfields,(const DataField*)db->field,found);CHKERRQ(ierr);
215   PetscFunctionReturn(0);
216 }
217 
218 #undef __FUNCT__
219 #define __FUNCT__ "DataBucketFinalize"
220 PetscErrorCode DataBucketFinalize(DataBucket db)
221 {
222   PetscFunctionBegin;
223   db->finalised = PETSC_TRUE;
224   PetscFunctionReturn(0);
225 }
226 
227 #undef __FUNCT__
228 #define __FUNCT__ "DataFieldGetNumEntries"
229 PetscErrorCode DataFieldGetNumEntries(DataField df,PetscInt *sum)
230 {
231   PetscFunctionBegin;
232   *sum = df->L;
233   PetscFunctionReturn(0);
234 }
235 
236 #undef __FUNCT__
237 #define __FUNCT__ "DataFieldSetBlockSize"
238 PetscErrorCode DataFieldSetBlockSize(DataField df,PetscInt blocksize)
239 {
240   PetscFunctionBegin;
241   df->bs = blocksize;
242   PetscFunctionReturn(0);
243 }
244 
245 #undef __FUNCT__
246 #define __FUNCT__ "DataFieldSetSize"
247 PetscErrorCode DataFieldSetSize(DataField df,const PetscInt new_L)
248 {
249   void          *tmp_data;
250   PetscErrorCode ierr;
251 
252   PetscFunctionBegin;
253   if (new_L <= 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot set size of DataField to be <= 0");
254   if (new_L == df->L) PetscFunctionReturn(0);
255   if (new_L > df->L) {
256     tmp_data = realloc( df->data, df->atomic_size * (new_L) );
257     df->data = tmp_data;
258     /* init new contents */
259     ierr = PetscMemzero(( ((char*)df->data)+df->L*df->atomic_size), (new_L-df->L)*df->atomic_size);CHKERRQ(ierr);
260   } else {
261     /* reallocate pointer list, add +1 in case new_L = 0 */
262     tmp_data = realloc( df->data, df->atomic_size * (new_L+1) );
263     df->data = tmp_data;
264   }
265   df->L = new_L;
266   PetscFunctionReturn(0);
267 }
268 
269 #undef __FUNCT__
270 #define __FUNCT__ "DataFieldZeroBlock"
271 PetscErrorCode DataFieldZeroBlock(DataField df,const PetscInt start,const PetscInt end)
272 {
273   PetscErrorCode ierr;
274 
275   PetscFunctionBegin;
276   if (start > end) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if start(%D) > end(%D)",start,end);
277   if (start < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if start(%D) < 0",start);
278   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);
279   ierr = PetscMemzero((((char*)df->data)+start*df->atomic_size), (end-start)*df->atomic_size);CHKERRQ(ierr);
280   PetscFunctionReturn(0);
281 }
282 
283 /*
284  A negative buffer value will simply be ignored and the old buffer value will be used.
285  */
286 #undef __FUNCT__
287 #define __FUNCT__ "DataBucketSetSizes"
288 PetscErrorCode DataBucketSetSizes(DataBucket db,const PetscInt L,const PetscInt buffer)
289 {
290   PetscInt       current_allocated,new_used,new_unused,new_buffer,new_allocated,f;
291   PetscBool      any_active_fields;
292   PetscErrorCode ierr;
293 
294   PetscFunctionBegin;
295   if (db->finalised == PETSC_FALSE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"You must call DataBucketFinalize() before DataBucketSetSizes()");
296   ierr = DataBucketQueryForActiveFields(db,&any_active_fields);CHKERRQ(ierr);
297   if (any_active_fields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot safely re-size as at least one DataField is currently being accessed");
298 
299   current_allocated = db->allocated;
300   new_used   = L;
301   new_unused = current_allocated - new_used;
302   new_buffer = db->buffer;
303   if (buffer >= 0) { /* update the buffer value */
304     new_buffer = buffer;
305   }
306   new_allocated = new_used + new_buffer;
307   /* action */
308   if (new_allocated > current_allocated) {
309     /* increase size to new_used + new_buffer */
310     for (f=0; f<db->nfields; f++) {
311       ierr = DataFieldSetSize(db->field[f], new_allocated);CHKERRQ(ierr);
312     }
313     db->L         = new_used;
314     db->buffer    = new_buffer;
315     db->allocated = new_used + new_buffer;
316   } else {
317     if (new_unused > 2 * new_buffer) {
318       /* shrink array to new_used + new_buffer */
319       for (f = 0; f < db->nfields; ++f) {
320         ierr = DataFieldSetSize(db->field[f], new_allocated);CHKERRQ(ierr);
321       }
322       db->L         = new_used;
323       db->buffer    = new_buffer;
324       db->allocated = new_used + new_buffer;
325     } else {
326       db->L      = new_used;
327       db->buffer = new_buffer;
328     }
329   }
330   /* zero all entries from db->L to db->allocated */
331   for (f = 0; f < db->nfields; ++f) {
332     DataField field = db->field[f];
333     ierr = DataFieldZeroBlock(field, db->L,db->allocated);CHKERRQ(ierr);
334   }
335   PetscFunctionReturn(0);
336 }
337 
338 #undef __FUNCT__
339 #define __FUNCT__ "DataBucketSetInitialSizes"
340 PetscErrorCode DataBucketSetInitialSizes(DataBucket db,const PetscInt L,const PetscInt buffer)
341 {
342   PetscInt       f;
343   PetscErrorCode ierr;
344 
345   PetscFunctionBegin;
346   ierr = DataBucketSetSizes(db,L,buffer);CHKERRQ(ierr);
347   for (f = 0; f < db->nfields; ++f) {
348     DataField field = db->field[f];
349     ierr = DataFieldZeroBlock(field,0,db->allocated);CHKERRQ(ierr);
350   }
351   PetscFunctionReturn(0);
352 }
353 
354 #undef __FUNCT__
355 #define __FUNCT__ "DataBucketGetSizes"
356 PetscErrorCode DataBucketGetSizes(DataBucket db,PetscInt *L,PetscInt *buffer,PetscInt *allocated)
357 {
358   PetscFunctionBegin;
359   if (L) {*L = db->L;}
360   if (buffer) {*buffer = db->buffer;}
361   if (allocated) {*allocated = db->allocated;}
362   PetscFunctionReturn(0);
363 }
364 
365 #undef __FUNCT__
366 #define __FUNCT__ "DataBucketGetGlobalSizes"
367 PetscErrorCode DataBucketGetGlobalSizes(MPI_Comm comm,DataBucket db,PetscInt *L,PetscInt *buffer,PetscInt *allocated)
368 {
369   PetscInt _L,_buffer,_allocated;
370   PetscInt ierr;
371 
372   PetscFunctionBegin;
373   _L = db->L;
374   _buffer = db->buffer;
375   _allocated = db->allocated;
376 
377   if (L) {         ierr = MPI_Allreduce(&_L,L,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); }
378   if (buffer) {    ierr = MPI_Allreduce(&_buffer,buffer,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); }
379   if (allocated) { ierr = MPI_Allreduce(&_allocated,allocated,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); }
380   PetscFunctionReturn(0);
381 }
382 
383 #undef __FUNCT__
384 #define __FUNCT__ "DataBucketGetDataFields"
385 PetscErrorCode DataBucketGetDataFields(DataBucket db,PetscInt *L,DataField *fields[])
386 {
387   PetscFunctionBegin;
388   if (L)      {*L      = db->nfields;}
389   if (fields) {*fields = db->field;}
390   PetscFunctionReturn(0);
391 }
392 
393 #undef __FUNCT__
394 #define __FUNCT__ "DataFieldGetAccess"
395 PetscErrorCode DataFieldGetAccess(const DataField gfield)
396 {
397   PetscFunctionBegin;
398   if (gfield->active) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is already active. You must call DataFieldRestoreAccess()",gfield->name);
399   gfield->active = PETSC_TRUE;
400   PetscFunctionReturn(0);
401 }
402 
403 #undef __FUNCT__
404 #define __FUNCT__ "DataFieldAccessPoint"
405 PetscErrorCode DataFieldAccessPoint(const DataField gfield,const PetscInt pid,void **ctx_p)
406 {
407   PetscFunctionBegin;
408 #ifdef DATAFIELD_POINT_ACCESS_GUARD
409   /* debug mode */
410   /* check point is valid */
411   if (pid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
412   if (pid >= gfield->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",gfield->L);
413   if (gfield->active == PETSC_FALSE) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is not active. You must call DataFieldGetAccess() before point data can be retrivied",gfield->name);
414 #endif
415   *ctx_p = __DATATFIELD_point_access(gfield->data,pid,gfield->atomic_size);
416   PetscFunctionReturn(0);
417 }
418 
419 #undef __FUNCT__
420 #define __FUNCT__ "DataFieldAccessPointOffset"
421 PetscErrorCode DataFieldAccessPointOffset(const DataField gfield,const size_t offset,const PetscInt pid,void **ctx_p)
422 {
423   PetscFunctionBegin;
424 #ifdef DATAFIELD_POINT_ACCESS_GUARD
425   /* debug mode */
426   /* check point is valid */
427   /* if( offset < 0 ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be >= 0");*/
428   /* Note compiler realizes this can never happen with an unsigned PetscInt */
429   if (offset >= gfield->atomic_size) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be < %zu",gfield->atomic_size);
430   /* check point is valid */
431   if (pid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
432   if (pid >= gfield->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",gfield->L);
433   if (gfield->active == PETSC_FALSE) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is not active. You must call DataFieldGetAccess() before point data can be retrivied",gfield->name);
434 #endif
435   *ctx_p = __DATATFIELD_point_access_offset(gfield->data,pid,gfield->atomic_size,offset);
436   PetscFunctionReturn(0);
437 }
438 
439 #undef __FUNCT__
440 #define __FUNCT__ "DataFieldRestoreAccess"
441 PetscErrorCode DataFieldRestoreAccess(DataField gfield)
442 {
443   PetscFunctionBegin;
444   if (gfield->active == PETSC_FALSE) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is not active. You must call DataFieldGetAccess()", gfield->name);
445   gfield->active = PETSC_FALSE;
446   PetscFunctionReturn(0);
447 }
448 
449 #undef __FUNCT__
450 #define __FUNCT__ "DataFieldVerifyAccess"
451 PetscErrorCode DataFieldVerifyAccess(const DataField gfield,const size_t size)
452 {
453   PetscFunctionBegin;
454 #ifdef DATAFIELD_POINT_ACCESS_GUARD
455   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 );
456 #endif
457   PetscFunctionReturn(0);
458 }
459 
460 #undef __FUNCT__
461 #define __FUNCT__ "DataFieldGetAtomicSize"
462 PetscErrorCode DataFieldGetAtomicSize(const DataField gfield,size_t *size)
463 {
464   PetscFunctionBegin;
465   if (size) {*size = gfield->atomic_size;}
466   PetscFunctionReturn(0);
467 }
468 
469 #undef __FUNCT__
470 #define __FUNCT__ "DataFieldGetEntries"
471 PetscErrorCode DataFieldGetEntries(const DataField gfield,void **data)
472 {
473   PetscFunctionBegin;
474   if (data) {*data = gfield->data;}
475   PetscFunctionReturn(0);
476 }
477 
478 #undef __FUNCT__
479 #define __FUNCT__ "DataFieldRestoreEntries"
480 PetscErrorCode DataFieldRestoreEntries(const DataField gfield,void **data)
481 {
482   PetscFunctionBegin;
483   if (data) {*data = NULL;}
484   PetscFunctionReturn(0);
485 }
486 
487 /* y = x */
488 #undef __FUNCT__
489 #define __FUNCT__ "DataBucketCopyPoint"
490 PetscErrorCode DataBucketCopyPoint(const DataBucket xb,const PetscInt pid_x,
491                          const DataBucket yb,const PetscInt pid_y)
492 {
493   PetscInt f;
494   PetscErrorCode ierr;
495 
496   PetscFunctionBegin;
497   for (f = 0; f < xb->nfields; ++f) {
498     void *dest;
499     void *src;
500 
501     ierr = DataFieldGetAccess(xb->field[f]);CHKERRQ(ierr);
502     if (xb != yb) { ierr = DataFieldGetAccess( yb->field[f]);CHKERRQ(ierr); }
503     ierr = DataFieldAccessPoint(xb->field[f],pid_x, &src);CHKERRQ(ierr);
504     ierr = DataFieldAccessPoint(yb->field[f],pid_y, &dest);CHKERRQ(ierr);
505     ierr = PetscMemcpy(dest, src, xb->field[f]->atomic_size);CHKERRQ(ierr);
506     ierr = DataFieldRestoreAccess(xb->field[f]);CHKERRQ(ierr);
507     if (xb != yb) {ierr = DataFieldRestoreAccess(yb->field[f]);CHKERRQ(ierr);}
508   }
509   PetscFunctionReturn(0);
510 }
511 
512 #undef __FUNCT__
513 #define __FUNCT__ "DataBucketCreateFromSubset"
514 PetscErrorCode DataBucketCreateFromSubset(DataBucket DBIn,const PetscInt N,const PetscInt list[],DataBucket *DB)
515 {
516   PetscInt nfields;
517   DataField *fields;
518   PetscInt f,L,buffer,allocated,p;
519   PetscErrorCode ierr;
520 
521   PetscFunctionBegin;
522   ierr = DataBucketCreate(DB);CHKERRQ(ierr);
523   /* copy contents of DBIn */
524   ierr = DataBucketGetDataFields(DBIn,&nfields,&fields);CHKERRQ(ierr);
525   ierr = DataBucketGetSizes(DBIn,&L,&buffer,&allocated);CHKERRQ(ierr);
526   for (f = 0; f < nfields; ++f) {
527     ierr = DataBucketRegisterField(*DB,"DataBucketCreateFromSubset",fields[f]->name,fields[f]->atomic_size,NULL);CHKERRQ(ierr);
528   }
529   ierr = DataBucketFinalize(*DB);CHKERRQ(ierr);
530   ierr = DataBucketSetSizes(*DB,L,buffer);CHKERRQ(ierr);
531   /* now copy the desired guys from DBIn => DB */
532   for (p = 0; p < N; ++p) {
533     ierr = DataBucketCopyPoint(DBIn,list[p], *DB,p);CHKERRQ(ierr);
534   }
535   PetscFunctionReturn(0);
536 }
537 
538 /* insert into an exisitng location */
539 #undef __FUNCT__
540 #define __FUNCT__ "DataFieldInsertPoint"
541 PetscErrorCode DataFieldInsertPoint(const DataField field,const PetscInt index,const void *ctx)
542 {
543   PetscErrorCode ierr;
544 
545   PetscFunctionBegin;
546 #ifdef DATAFIELD_POINT_ACCESS_GUARD
547   /* check point is valid */
548   if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
549   if (index >= field->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",field->L);
550 #endif
551   ierr = PetscMemcpy(__DATATFIELD_point_access(field->data,index,field->atomic_size), ctx, field->atomic_size);CHKERRQ(ierr);
552   PetscFunctionReturn(0);
553 }
554 
555 /* remove data at index - replace with last point */
556 #undef __FUNCT__
557 #define __FUNCT__ "DataBucketRemovePointAtIndex"
558 PetscErrorCode DataBucketRemovePointAtIndex(const DataBucket db,const PetscInt index)
559 {
560   PetscInt       f;
561   PetscBool      any_active_fields;
562   PetscErrorCode ierr;
563 
564   PetscFunctionBegin;
565 #ifdef DATAFIELD_POINT_ACCESS_GUARD
566   /* check point is valid */
567   if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
568   if (index >= db->allocated) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",db->L+db->buffer);
569 #endif
570   ierr = DataBucketQueryForActiveFields(db,&any_active_fields);CHKERRQ(ierr);
571   if (any_active_fields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot safely remove point as at least one DataField is currently being accessed");
572   if (index >= db->L) { /* this point is not in the list - no need to error, but I will anyway */
573     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);
574   }
575   if (index != db->L-1) { /* not last point in list */
576     for (f = 0; f < db->nfields; ++f) {
577       DataField field = db->field[f];
578 
579       /* copy then remove */
580       ierr = DataFieldCopyPoint(db->L-1, field, index, field);CHKERRQ(ierr);
581       /* DataFieldZeroPoint(field,index); */
582     }
583   }
584   /* decrement size */
585   /* this will zero out an crap at the end of the list */
586   ierr = DataBucketRemovePoint(db);CHKERRQ(ierr);
587   PetscFunctionReturn(0);
588 }
589 
590 /* copy x into y */
591 #undef __FUNCT__
592 #define __FUNCT__ "DataFieldCopyPoint"
593 PetscErrorCode DataFieldCopyPoint(const PetscInt pid_x,const DataField field_x,
594                         const PetscInt pid_y,const DataField field_y )
595 {
596   PetscErrorCode ierr;
597 
598   PetscFunctionBegin;
599 #ifdef DATAFIELD_POINT_ACCESS_GUARD
600   /* check point is valid */
601   if (pid_x < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"(IN) index must be >= 0");
602   if (pid_x >= field_x->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"(IN) index must be < %D",field_x->L);
603   if (pid_y < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"(OUT) index must be >= 0");
604   if (pid_y >= field_y->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"(OUT) index must be < %D",field_y->L);
605   if( field_y->atomic_size != field_x->atomic_size ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"atomic size must match");
606 #endif
607   ierr = PetscMemcpy(__DATATFIELD_point_access(field_y->data,pid_y,field_y->atomic_size),
608                      __DATATFIELD_point_access(field_x->data,pid_x,field_x->atomic_size),
609                      field_y->atomic_size);CHKERRQ(ierr);
610   PetscFunctionReturn(0);
611 }
612 
613 
614 /* zero only the datafield at this point */
615 #undef __FUNCT__
616 #define __FUNCT__ "DataFieldZeroPoint"
617 PetscErrorCode DataFieldZeroPoint(const DataField field,const PetscInt index)
618 {
619   PetscErrorCode ierr;
620 
621   PetscFunctionBegin;
622 #ifdef DATAFIELD_POINT_ACCESS_GUARD
623   /* check point is valid */
624   if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
625   if (index >= field->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",field->L);
626 #endif
627   ierr = PetscMemzero(__DATATFIELD_point_access(field->data,index,field->atomic_size), field->atomic_size);CHKERRQ(ierr);
628   PetscFunctionReturn(0);
629 }
630 
631 /* zero ALL data for this point */
632 #undef __FUNCT__
633 #define __FUNCT__ "DataBucketZeroPoint"
634 PetscErrorCode DataBucketZeroPoint(const DataBucket db,const PetscInt index)
635 {
636   PetscInt f;
637   PetscErrorCode ierr;
638 
639   PetscFunctionBegin;
640   /* check point is valid */
641   if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
642   if (index >= db->allocated) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",db->allocated);
643   for (f = 0; f < db->nfields; ++f) {
644     DataField field = db->field[f];
645     ierr = DataFieldZeroPoint(field,index);CHKERRQ(ierr);
646   }
647   PetscFunctionReturn(0);
648 }
649 
650 /* increment */
651 #undef __FUNCT__
652 #define __FUNCT__ "DataBucketAddPoint"
653 PetscErrorCode DataBucketAddPoint(DataBucket db)
654 {
655   PetscErrorCode ierr;
656 
657   PetscFunctionBegin;
658   ierr = DataBucketSetSizes(db,db->L+1,-1);CHKERRQ(ierr);
659   PetscFunctionReturn(0);
660 }
661 
662 /* decrement */
663 #undef __FUNCT__
664 #define __FUNCT__ "DataBucketRemovePoint"
665 PetscErrorCode DataBucketRemovePoint(DataBucket db)
666 {
667   PetscErrorCode ierr;
668 
669   PetscFunctionBegin;
670   ierr = DataBucketSetSizes(db,db->L-1,-1);CHKERRQ(ierr);
671   PetscFunctionReturn(0);
672 }
673 
674 #undef __FUNCT__
675 #define __FUNCT__ "_DataFieldViewBinary"
676 PetscErrorCode _DataFieldViewBinary(DataField field,FILE *fp)
677 {
678   PetscErrorCode ierr;
679 
680   PetscFunctionBegin;
681   ierr = PetscFPrintf(PETSC_COMM_SELF, fp,"<DataField>\n");CHKERRQ(ierr);
682   ierr = PetscFPrintf(PETSC_COMM_SELF, fp,"%D\n", field->L);CHKERRQ(ierr);
683   ierr = PetscFPrintf(PETSC_COMM_SELF, fp,"%zu\n",field->atomic_size);CHKERRQ(ierr);
684   ierr = PetscFPrintf(PETSC_COMM_SELF, fp,"%s\n", field->registeration_function);CHKERRQ(ierr);
685   ierr = PetscFPrintf(PETSC_COMM_SELF, fp,"%s\n", field->name);CHKERRQ(ierr);
686   fwrite(field->data, field->atomic_size, field->L, fp);
687   ierr = PetscFPrintf(PETSC_COMM_SELF, fp,"\n</DataField>\n");CHKERRQ(ierr);
688   PetscFunctionReturn(0);
689 }
690 
691 #undef __FUNCT__
692 #define __FUNCT__ "_DataBucketRegisterFieldFromFile"
693 PetscErrorCode _DataBucketRegisterFieldFromFile(FILE *fp,DataBucket db)
694 {
695   PetscBool      val;
696   DataField     *field;
697   DataField      gfield;
698   char           dummy[100];
699   char           registeration_function[5000];
700   char           field_name[5000];
701   PetscInt       L;
702   size_t         atomic_size,strL;
703   PetscErrorCode ierr;
704 
705   PetscFunctionBegin;
706   /* check we haven't finalised the registration of fields */
707   /*
708    if(db->finalised==PETSC_TRUE) {
709    printf("ERROR: DataBucketFinalize() has been called. Cannot register more fields\n");
710    ERROR();
711    }
712    */
713   /* read file contents */
714   fgets(dummy,99,fp);
715   fscanf(fp, "%" PetscInt_FMT "\n",&L);
716   fscanf(fp, "%zu\n",&atomic_size);
717   fgets(registeration_function,4999,fp);
718   strL = strlen(registeration_function);
719   if (strL > 1) {
720     registeration_function[strL-1] = 0;
721   }
722   fgets(field_name,4999,fp);
723   strL = strlen(field_name);
724   if (strL > 1) {
725     field_name[strL-1] = 0;
726   }
727 
728 #ifdef DATA_BUCKET_LOG
729   ierr = PetscPrintf(PETSC_COMM_SELF,"  ** read L=%D; atomic_size=%zu; reg_func=\"%s\"; name=\"%s\" \n", L,atomic_size,registeration_function,field_name);CHKERRQ(ierr);
730 #endif
731   /* check for repeated name */
732   ierr = StringInList( field_name, db->nfields, (const DataField*)db->field, &val );CHKERRQ(ierr);
733   if (val == PETSC_TRUE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot add same field twice");
734   /* create new space for data */
735   field = (DataField *) realloc( db->field,     sizeof(DataField)*(db->nfields+1));
736   db->field     = field;
737   /* add field */
738   ierr = DataFieldCreate( registeration_function, field_name, atomic_size, L, &gfield );CHKERRQ(ierr);
739   /* copy contents of file */
740   fread(gfield->data, gfield->atomic_size, gfield->L, fp);
741 #ifdef DATA_BUCKET_LOG
742   ierr = PetscPrintf(PETSC_COMM_SELF,"  ** read %zu bytes for DataField \"%s\" \n", gfield->atomic_size * gfield->L, field_name);CHKERRQ(ierr);
743 #endif
744   /* finish reading meta data */
745   fgets(dummy,99,fp);
746   fgets(dummy,99,fp);
747   db->field[db->nfields] = gfield;
748   db->nfields++;
749   PetscFunctionReturn(0);
750 }
751 
752 #undef __FUNCT__
753 #define __FUNCT__ "_DataBucketViewAscii_HeaderWrite_v00"
754 PetscErrorCode _DataBucketViewAscii_HeaderWrite_v00(FILE *fp)
755 {
756   PetscErrorCode ierr;
757 
758   PetscFunctionBegin;
759   ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"<DataBucketHeader>\n");CHKERRQ(ierr);
760   ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"type=DataBucket\n");CHKERRQ(ierr);
761   ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"format=ascii\n");CHKERRQ(ierr);
762   ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"version=0.0\n");CHKERRQ(ierr);
763   ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"options=\n");CHKERRQ(ierr);
764   ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"</DataBucketHeader>\n");CHKERRQ(ierr);
765   PetscFunctionReturn(0);
766 }
767 
768 #undef __FUNCT__
769 #define __FUNCT__ "_DataBucketViewAscii_HeaderRead_v00"
770 PetscErrorCode _DataBucketViewAscii_HeaderRead_v00(FILE *fp)
771 {
772   char           dummy[100];
773   size_t         strL;
774   PetscBool      flg;
775   PetscErrorCode ierr;
776 
777   PetscFunctionBegin;
778   /* header open */
779   fgets(dummy,99,fp);
780 
781   /* type */
782   fgets(dummy,99,fp);
783   strL = strlen(dummy);
784   if (strL > 1) {dummy[strL-1] = 0;}
785   ierr = PetscStrcmp(dummy, "type=DataBucket", &flg);CHKERRQ(ierr);
786   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Data file doesn't contain a DataBucket type");
787   /* format */
788   fgets(dummy,99,fp);
789   /* version */
790   fgets(dummy,99,fp);
791   strL = strlen(dummy);
792   if (strL > 1) { dummy[strL-1] = 0; }
793   ierr = PetscStrcmp(dummy, "version=0.0", &flg);CHKERRQ(ierr);
794   if (!flg) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"DataBucket file must be parsed with version=0.0 : You tried %s", dummy);
795   /* options */
796   fgets(dummy,99,fp);
797   /* header close */
798   fgets(dummy,99,fp);
799   PetscFunctionReturn(0);
800 }
801 
802 #undef __FUNCT__
803 #define __FUNCT__ "_DataBucketLoadFromFileBinary_SEQ"
804 PetscErrorCode _DataBucketLoadFromFileBinary_SEQ(const char filename[],DataBucket *_db)
805 {
806   DataBucket db;
807   FILE *fp;
808   PetscInt L,buffer,f,nfields;
809   PetscErrorCode ierr;
810 
811   PetscFunctionBegin;
812 #ifdef DATA_BUCKET_LOG
813   ierr = PetscPrintf(PETSC_COMM_SELF,"** DataBucketLoadFromFile **\n");CHKERRQ(ierr);
814 #endif
815   /* open file */
816   fp = fopen(filename,"rb");
817   if (fp == NULL) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file with name %s", filename);
818   /* read header */
819   ierr = _DataBucketViewAscii_HeaderRead_v00(fp);CHKERRQ(ierr);
820   fscanf(fp,"%" PetscInt_FMT "\n%" PetscInt_FMT "\n%" PetscInt_FMT "\n",&L,&buffer,&nfields);
821   ierr = DataBucketCreate(&db);CHKERRQ(ierr);
822   for (f = 0; f < nfields; ++f) {
823     ierr = _DataBucketRegisterFieldFromFile(fp,db);CHKERRQ(ierr);
824   }
825   fclose(fp);
826   ierr = DataBucketFinalize(db);CHKERRQ(ierr);
827   /*
828    DataBucketSetSizes(db,L,buffer);
829    */
830   db->L = L;
831   db->buffer = buffer;
832   db->allocated = L + buffer;
833   *_db = db;
834   PetscFunctionReturn(0);
835 }
836 
837 #undef __FUNCT__
838 #define __FUNCT__ "DataBucketLoadFromFile"
839 PetscErrorCode DataBucketLoadFromFile(MPI_Comm comm,const char filename[],DataBucketViewType type,DataBucket *db)
840 {
841   PetscMPIInt    nproc,rank;
842   PetscErrorCode ierr;
843 
844   PetscFunctionBegin;
845   ierr = MPI_Comm_size(comm,&nproc);CHKERRQ(ierr);
846   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
847 #ifdef DATA_BUCKET_LOG
848   ierr = PetscPrintf(PETSC_COMM_SELF,"** DataBucketLoadFromFile **\n");CHKERRQ(ierr);
849 #endif
850   if (type == DATABUCKET_VIEW_STDOUT) {
851   } else if (type == DATABUCKET_VIEW_ASCII) {
852     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot be implemented as we don't know the underlying particle data structure");
853   } else if (type == DATABUCKET_VIEW_BINARY) {
854     if (nproc == 1) {
855       ierr = _DataBucketLoadFromFileBinary_SEQ(filename,db);CHKERRQ(ierr);
856     } else {
857       char name[PETSC_MAX_PATH_LEN];
858 
859       ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "%s_p%1.5d", filename, rank);CHKERRQ(ierr);
860       ierr = _DataBucketLoadFromFileBinary_SEQ(name, db);CHKERRQ(ierr);
861     }
862   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown viewer requested");
863   PetscFunctionReturn(0);
864 }
865 
866 #undef __FUNCT__
867 #define __FUNCT__ "_DataBucketViewBinary"
868 PetscErrorCode _DataBucketViewBinary(DataBucket db,const char filename[])
869 {
870   FILE          *fp = NULL;
871   PetscInt       f;
872   PetscErrorCode ierr;
873 
874   PetscFunctionBegin;
875   fp = fopen(filename,"wb");
876   if (fp == NULL) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot open file with name %s", filename);
877   /* db header */
878   ierr =_DataBucketViewAscii_HeaderWrite_v00(fp);CHKERRQ(ierr);
879   /* meta-data */
880   ierr = PetscFPrintf(PETSC_COMM_SELF, fp, "%D\n%D\n%D\n", db->L,db->buffer,db->nfields);CHKERRQ(ierr);
881   /* load datafields */
882   for (f = 0; f < db->nfields; ++f) {
883     ierr = _DataFieldViewBinary(db->field[f],fp);CHKERRQ(ierr);
884   }
885   fclose(fp);
886   PetscFunctionReturn(0);
887 }
888 
889 #undef __FUNCT__
890 #define __FUNCT__ "DataBucketView_SEQ"
891 PetscErrorCode DataBucketView_SEQ(DataBucket db,const char filename[],DataBucketViewType type)
892 {
893   PetscErrorCode ierr;
894 
895   PetscFunctionBegin;
896   switch (type) {
897   case DATABUCKET_VIEW_STDOUT:
898   {
899     PetscInt f;
900     double memory_usage_total = 0.0;
901 
902     ierr = PetscPrintf(PETSC_COMM_SELF,"DataBucketView(SEQ): (\"%s\")\n",filename);CHKERRQ(ierr);
903     ierr = PetscPrintf(PETSC_COMM_SELF,"  L                  = %D \n", db->L);CHKERRQ(ierr);
904     ierr = PetscPrintf(PETSC_COMM_SELF,"  buffer             = %D \n", db->buffer);CHKERRQ(ierr);
905     ierr = PetscPrintf(PETSC_COMM_SELF,"  allocated          = %D \n", db->allocated);CHKERRQ(ierr);
906     ierr = PetscPrintf(PETSC_COMM_SELF,"  nfields registered = %D \n", db->nfields);CHKERRQ(ierr);
907     for (f = 0; f < db->nfields; ++f) {
908       double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
909       ierr = PetscPrintf(PETSC_COMM_SELF,"    [%3D]: field name  ==>> %30s : Mem. usage = %1.2e (MB) \n", f, db->field[f]->name, memory_usage_f );CHKERRQ(ierr);
910       ierr = PetscPrintf(PETSC_COMM_SELF,"           blocksize          = %D \n", db->field[f]->bs);CHKERRQ(ierr);
911       if (db->field[f]->bs != 1) {
912         ierr = PetscPrintf(PETSC_COMM_SELF,"           atomic size        = %zu [full block, bs=%D]\n", db->field[f]->atomic_size,db->field[f]->bs);CHKERRQ(ierr);
913         ierr = PetscPrintf(PETSC_COMM_SELF,"           atomic size/item   = %zu \n", db->field[f]->atomic_size/db->field[f]->bs);CHKERRQ(ierr);
914       } else {
915         ierr = PetscPrintf(PETSC_COMM_SELF,"           atomic size        = %zu \n", db->field[f]->atomic_size);CHKERRQ(ierr);
916       }
917       memory_usage_total += memory_usage_f;
918     }
919     ierr = PetscPrintf(PETSC_COMM_SELF,"  Total mem. usage                                                      = %1.2e (MB) \n", memory_usage_total);CHKERRQ(ierr);
920   }
921   break;
922   case DATABUCKET_VIEW_ASCII:
923   {
924     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot be implemented as we don't know the underlying particle data structure");
925   }
926   break;
927   case DATABUCKET_VIEW_BINARY:
928   {
929     ierr = _DataBucketViewBinary(db,filename);CHKERRQ(ierr);
930   }
931   break;
932   case DATABUCKET_VIEW_HDF5:
933   {
934     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No HDF5 support");
935   }
936   break;
937   default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer method requested");
938   }
939   PetscFunctionReturn(0);
940 }
941 
942 #undef __FUNCT__
943 #define __FUNCT__ "DataBucketView_MPI"
944 PetscErrorCode DataBucketView_MPI(MPI_Comm comm,DataBucket db,const char filename[],DataBucketViewType type)
945 {
946   PetscErrorCode ierr;
947 
948   PetscFunctionBegin;
949   switch (type) {
950   case DATABUCKET_VIEW_STDOUT:
951   {
952     PetscInt f;
953     PetscInt L,buffer,allocated;
954     double memory_usage_total,memory_usage_total_local = 0.0;
955     PetscMPIInt rank;
956 
957     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
958     ierr = DataBucketGetGlobalSizes(comm,db,&L,&buffer,&allocated);CHKERRQ(ierr);
959     for (f = 0; f < db->nfields; ++f) {
960       double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
961       memory_usage_total_local += memory_usage_f;
962     }
963     ierr = MPI_Allreduce(&memory_usage_total_local,&memory_usage_total,1,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
964     ierr = PetscPrintf(comm,"DataBucketView(MPI): (\"%s\")\n",filename);CHKERRQ(ierr);
965     ierr = PetscPrintf(comm,"  L                  = %D \n", L);CHKERRQ(ierr);
966     ierr = PetscPrintf(comm,"  buffer (max)       = %D \n", buffer);CHKERRQ(ierr);
967     ierr = PetscPrintf(comm,"  allocated          = %D \n", allocated);CHKERRQ(ierr);
968     ierr = PetscPrintf(comm,"  nfields registered = %D \n", db->nfields);CHKERRQ(ierr);
969     for (f=0; f<db->nfields; f++) {
970       double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
971       ierr = PetscPrintf(PETSC_COMM_SELF,"    [%3D]: field name  ==>> %30s : Mem. usage = %1.2e (MB) : rank0\n", f, db->field[f]->name, memory_usage_f);CHKERRQ(ierr);
972     }
973     ierr = PetscPrintf(PETSC_COMM_SELF,"  Total mem. usage                                                      = %1.2e (MB) : collective\n", memory_usage_total);CHKERRQ(ierr);
974   }
975   break;
976   case DATABUCKET_VIEW_ASCII:
977   {
978     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot be implemented as we don't know the underlying data structure");
979   }
980   break;
981   case DATABUCKET_VIEW_BINARY:
982   {
983     char        name[PETSC_MAX_PATH_LEN];
984     PetscMPIInt rank;
985 
986     /* create correct extension */
987     ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
988     ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "%s_p%1.5d", filename, rank);CHKERRQ(ierr);
989     ierr = _DataBucketViewBinary(db, name);CHKERRQ(ierr);
990   }
991   break;
992   case DATABUCKET_VIEW_HDF5:
993   {
994     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for HDF5");
995   }
996   break;
997   default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer method requested");
998   }
999   PetscFunctionReturn(0);
1000 }
1001 
1002 #undef __FUNCT__
1003 #define __FUNCT__ "DataBucketView"
1004 PetscErrorCode DataBucketView(MPI_Comm comm,DataBucket db,const char filename[],DataBucketViewType type)
1005 {
1006   PetscMPIInt nproc;
1007   PetscErrorCode ierr;
1008 
1009   PetscFunctionBegin;
1010   ierr = MPI_Comm_size(comm,&nproc);CHKERRQ(ierr);
1011   if (nproc == 1) {
1012     ierr = DataBucketView_SEQ(db,filename,type);CHKERRQ(ierr);
1013   } else {
1014     ierr = DataBucketView_MPI(comm,db,filename,type);CHKERRQ(ierr);
1015   }
1016   PetscFunctionReturn(0);
1017 }
1018 
1019 #undef __FUNCT__
1020 #define __FUNCT__ "DataBucketDuplicateFields"
1021 PetscErrorCode DataBucketDuplicateFields(DataBucket dbA,DataBucket *dbB)
1022 {
1023   DataBucket db2;
1024   PetscInt f;
1025   PetscErrorCode ierr;
1026 
1027   PetscFunctionBegin;
1028   ierr = DataBucketCreate(&db2);CHKERRQ(ierr);
1029   /* copy contents from dbA into db2 */
1030   for (f = 0; f < dbA->nfields; ++f) {
1031     DataField field;
1032     size_t    atomic_size;
1033     char      *name;
1034 
1035     field = dbA->field[f];
1036     atomic_size = field->atomic_size;
1037     name        = field->name;
1038     ierr = DataBucketRegisterField(db2,"DataBucketDuplicateFields",name,atomic_size,NULL);CHKERRQ(ierr);
1039   }
1040   ierr = DataBucketFinalize(db2);CHKERRQ(ierr);
1041   ierr = DataBucketSetInitialSizes(db2,0,1000);CHKERRQ(ierr);
1042   *dbB = db2;
1043   PetscFunctionReturn(0);
1044 }
1045 
1046 /*
1047  Insert points from db2 into db1
1048  db1 <<== db2
1049  */
1050 #undef __FUNCT__
1051 #define __FUNCT__ "DataBucketInsertValues"
1052 PetscErrorCode DataBucketInsertValues(DataBucket db1,DataBucket db2)
1053 {
1054   PetscInt n_mp_points1,n_mp_points2;
1055   PetscInt n_mp_points1_new,p;
1056   PetscErrorCode ierr;
1057 
1058   PetscFunctionBegin;
1059   ierr = DataBucketGetSizes(db1,&n_mp_points1,0,0);CHKERRQ(ierr);
1060   ierr = DataBucketGetSizes(db2,&n_mp_points2,0,0);CHKERRQ(ierr);
1061   n_mp_points1_new = n_mp_points1 + n_mp_points2;
1062   ierr = DataBucketSetSizes(db1,n_mp_points1_new,-1);CHKERRQ(ierr);
1063   for (p = 0; p < n_mp_points2; ++p) {
1064     /* db1 <<== db2 */
1065     ierr = DataBucketCopyPoint(db2,p, db1,(n_mp_points1 + p));CHKERRQ(ierr);
1066   }
1067   PetscFunctionReturn(0);
1068 }
1069 
1070 /* helpers for parallel send/recv */
1071 #undef __FUNCT__
1072 #define __FUNCT__ "DataBucketCreatePackedArray"
1073 PetscErrorCode DataBucketCreatePackedArray(DataBucket db,size_t *bytes,void **buf)
1074 {
1075   PetscInt       f;
1076   size_t         sizeof_marker_contents;
1077   void          *buffer;
1078   PetscErrorCode ierr;
1079 
1080   PetscFunctionBegin;
1081   sizeof_marker_contents = 0;
1082   for (f = 0; f < db->nfields; ++f) {
1083     DataField df = db->field[f];
1084     sizeof_marker_contents += df->atomic_size;
1085   }
1086   ierr = PetscMalloc(sizeof_marker_contents, &buffer);CHKERRQ(ierr);
1087   ierr = PetscMemzero(buffer, sizeof_marker_contents);CHKERRQ(ierr);
1088   if (bytes) {*bytes = sizeof_marker_contents;}
1089   if (buf)   {*buf   = buffer;}
1090   PetscFunctionReturn(0);
1091 }
1092 
1093 #undef __FUNCT__
1094 #define __FUNCT__ "DataBucketDestroyPackedArray"
1095 PetscErrorCode DataBucketDestroyPackedArray(DataBucket db,void **buf)
1096 {
1097   PetscErrorCode ierr;
1098 
1099   PetscFunctionBegin;
1100   if (buf) {
1101     ierr = PetscFree(*buf);CHKERRQ(ierr);
1102     *buf = NULL;
1103   }
1104   PetscFunctionReturn(0);
1105 }
1106 
1107 #undef __FUNCT__
1108 #define __FUNCT__ "DataBucketFillPackedArray"
1109 PetscErrorCode DataBucketFillPackedArray(DataBucket db,const PetscInt index,void *buf)
1110 {
1111   PetscInt       f;
1112   void          *data, *data_p;
1113   size_t         asize, offset;
1114   PetscErrorCode ierr;
1115 
1116   PetscFunctionBegin;
1117   offset = 0;
1118   for (f = 0; f < db->nfields; ++f) {
1119     DataField df = db->field[f];
1120 
1121     asize = df->atomic_size;
1122     data = (void*)( df->data );
1123     data_p = (void*)( (char*)data + index*asize );
1124     ierr = PetscMemcpy((void*)((char*)buf + offset), data_p, asize);CHKERRQ(ierr);
1125     offset = offset + asize;
1126   }
1127   PetscFunctionReturn(0);
1128 }
1129 
1130 #undef __FUNCT__
1131 #define __FUNCT__ "DataBucketInsertPackedArray"
1132 PetscErrorCode DataBucketInsertPackedArray(DataBucket db,const PetscInt idx,void *data)
1133 {
1134   PetscInt f;
1135   void *data_p;
1136   size_t offset;
1137   PetscErrorCode ierr;
1138 
1139   PetscFunctionBegin;
1140   offset = 0;
1141   for (f = 0; f < db->nfields; ++f) {
1142     DataField df = db->field[f];
1143 
1144     data_p = (void*)( (char*)data + offset );
1145     ierr = DataFieldInsertPoint(df, idx, (void*)data_p);CHKERRQ(ierr);
1146     offset = offset + df->atomic_size;
1147   }
1148   PetscFunctionReturn(0);
1149 }
1150