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