xref: /petsc/src/dm/impls/swarm/data_bucket.c (revision b41ce5d507ea9a58bfa83cf403107a702e77a67d)
1 #include "data_bucket.h"
2 
3 /* string helpers */
4 PetscErrorCode StringInList(const char name[],const PetscInt N,const DataField 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 StringFindInList(const char name[],const PetscInt N,const DataField 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 DataFieldCreate(const char registeration_function[],const char name[],const size_t size,const PetscInt L,DataField *DF)
41 {
42   DataField      df;
43   PetscErrorCode ierr;
44 
45   PetscFunctionBegin;
46   ierr = PetscMalloc(sizeof(struct _p_DataField), &df);CHKERRQ(ierr);
47   ierr = PetscMemzero(df, sizeof(struct _p_DataField));CHKERRQ(ierr);
48   ierr = PetscStrallocpy(registeration_function, &df->registeration_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 DataFieldDestroy(DataField *DF)
61 {
62   DataField      df = *DF;
63   PetscErrorCode ierr;
64 
65   PetscFunctionBegin;
66   ierr = PetscFree(df->registeration_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 DataBucketCreate(DataBucket *DB)
76 {
77   DataBucket     db;
78   PetscErrorCode ierr;
79 
80   PetscFunctionBegin;
81   ierr = PetscMalloc(sizeof(struct _p_DataBucket), &db);CHKERRQ(ierr);
82   ierr = PetscMemzero(db, sizeof(struct _p_DataBucket));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 DataBucketDestroy(DataBucket *DB)
96 {
97   DataBucket     db = *DB;
98   PetscInt       f;
99   PetscErrorCode ierr;
100 
101   PetscFunctionBegin;
102   /* release fields */
103   for (f = 0; f < db->nfields; ++f) {
104     ierr = DataFieldDestroy(&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 DataBucketQueryForActiveFields(DataBucket 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 DataBucketRegisterField(
131                               DataBucket db,
132                               const char registeration_function[],
133                               const char field_name[],
134                               size_t atomic_size, DataField *_gfield)
135 {
136   PetscBool val;
137   DataField 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: DataBucketFinalize() has been called. Cannot register more fields\n");
145    ERROR();
146    }
147    */
148   /* check for repeated name */
149   ierr = StringInList(field_name, db->nfields, (const DataField*) 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(DataField)*(db->nfields+1), &db->field);CHKERRQ(ierr);
153   /* add field */
154   ierr = DataFieldCreate(registeration_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 DataBucketRegisterField(db,name,size,k) {\
165  char *location;\
166  asprintf(&location,"Registered by %s() at line %d within file %s", __FUNCTION__, __LINE__, __FILE__);\
167  _DataBucketRegisterField( (db), location, (name), (size), (k) );\
168  ierr = PetscFree(location);\
169  }
170  */
171 
172 PetscErrorCode DataBucketGetDataFieldByName(DataBucket db,const char name[],DataField *gfield)
173 {
174   PetscInt       idx;
175   PetscBool      found;
176   PetscErrorCode ierr;
177 
178   PetscFunctionBegin;
179   ierr = StringInList(name,db->nfields,(const DataField*)db->field,&found);CHKERRQ(ierr);
180   if (!found) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot find DataField with name %s",name);
181   ierr = StringFindInList(name,db->nfields,(const DataField*)db->field,&idx);CHKERRQ(ierr);
182   *gfield = db->field[idx];
183   PetscFunctionReturn(0);
184 }
185 
186 PetscErrorCode DataBucketQueryDataFieldByName(DataBucket db,const char name[],PetscBool *found)
187 {
188   PetscErrorCode ierr;
189 
190   PetscFunctionBegin;
191   *found = PETSC_FALSE;
192   ierr = StringInList(name,db->nfields,(const DataField*)db->field,found);CHKERRQ(ierr);
193   PetscFunctionReturn(0);
194 }
195 
196 PetscErrorCode DataBucketFinalize(DataBucket db)
197 {
198   PetscFunctionBegin;
199   db->finalised = PETSC_TRUE;
200   PetscFunctionReturn(0);
201 }
202 
203 PetscErrorCode DataFieldGetNumEntries(DataField df,PetscInt *sum)
204 {
205   PetscFunctionBegin;
206   *sum = df->L;
207   PetscFunctionReturn(0);
208 }
209 
210 PetscErrorCode DataFieldSetBlockSize(DataField df,PetscInt blocksize)
211 {
212   PetscFunctionBegin;
213   df->bs = blocksize;
214   PetscFunctionReturn(0);
215 }
216 
217 PetscErrorCode DataFieldSetSize(DataField 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 DataField 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 DataFieldZeroBlock(DataField 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 DataBucketSetSizes(DataBucket 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 DataBucketFinalize() before DataBucketSetSizes()");
259   ierr = DataBucketQueryForActiveFields(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 DataField 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 = DataFieldSetSize(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 = DataFieldSetSize(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     DataField field = db->field[f];
296     ierr = DataFieldZeroBlock(field, db->L,db->allocated);CHKERRQ(ierr);
297   }
298   PetscFunctionReturn(0);
299 }
300 
301 PetscErrorCode DataBucketSetInitialSizes(DataBucket db,const PetscInt L,const PetscInt buffer)
302 {
303   PetscInt       f;
304   PetscErrorCode ierr;
305 
306   PetscFunctionBegin;
307   ierr = DataBucketSetSizes(db,L,buffer);CHKERRQ(ierr);
308   for (f = 0; f < db->nfields; ++f) {
309     DataField field = db->field[f];
310     ierr = DataFieldZeroBlock(field,0,db->allocated);CHKERRQ(ierr);
311   }
312   PetscFunctionReturn(0);
313 }
314 
315 PetscErrorCode DataBucketGetSizes(DataBucket 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 DataBucketGetGlobalSizes(MPI_Comm comm,DataBucket 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 DataBucketGetDataFields(DataBucket db,PetscInt *L,DataField *fields[])
341 {
342   PetscFunctionBegin;
343   if (L)      {*L      = db->nfields;}
344   if (fields) {*fields = db->field;}
345   PetscFunctionReturn(0);
346 }
347 
348 PetscErrorCode DataFieldGetAccess(const DataField gfield)
349 {
350   PetscFunctionBegin;
351   if (gfield->active) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is already active. You must call DataFieldRestoreAccess()",gfield->name);
352   gfield->active = PETSC_TRUE;
353   PetscFunctionReturn(0);
354 }
355 
356 PetscErrorCode DataFieldAccessPoint(const DataField gfield,const PetscInt pid,void **ctx_p)
357 {
358   PetscFunctionBegin;
359   *ctx_p = NULL;
360 #ifdef 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 DataFieldGetAccess() before point data can be retrivied",gfield->name);
366 #endif
367   *ctx_p = __DATATFIELD_point_access(gfield->data,pid,gfield->atomic_size);
368   PetscFunctionReturn(0);
369 }
370 
371 PetscErrorCode DataFieldAccessPointOffset(const DataField gfield,const size_t offset,const PetscInt pid,void **ctx_p)
372 {
373   PetscFunctionBegin;
374 #ifdef 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 DataFieldGetAccess() before point data can be retrivied",gfield->name);
384 #endif
385   *ctx_p = __DATATFIELD_point_access_offset(gfield->data,pid,gfield->atomic_size,offset);
386   PetscFunctionReturn(0);
387 }
388 
389 PetscErrorCode DataFieldRestoreAccess(DataField 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 DataFieldGetAccess()", gfield->name);
393   gfield->active = PETSC_FALSE;
394   PetscFunctionReturn(0);
395 }
396 
397 PetscErrorCode DataFieldVerifyAccess(const DataField gfield,const size_t size)
398 {
399   PetscFunctionBegin;
400 #ifdef 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 DataFieldGetAtomicSize(const DataField gfield,size_t *size)
407 {
408   PetscFunctionBegin;
409   if (size) {*size = gfield->atomic_size;}
410   PetscFunctionReturn(0);
411 }
412 
413 PetscErrorCode DataFieldGetEntries(const DataField gfield,void **data)
414 {
415   PetscFunctionBegin;
416   if (data) {*data = gfield->data;}
417   PetscFunctionReturn(0);
418 }
419 
420 PetscErrorCode DataFieldRestoreEntries(const DataField gfield,void **data)
421 {
422   PetscFunctionBegin;
423   if (data) {*data = NULL;}
424   PetscFunctionReturn(0);
425 }
426 
427 /* y = x */
428 PetscErrorCode DataBucketCopyPoint(const DataBucket xb,const PetscInt pid_x,
429                          const DataBucket 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 = DataFieldGetAccess(xb->field[f]);CHKERRQ(ierr);
440     if (xb != yb) { ierr = DataFieldGetAccess( yb->field[f]);CHKERRQ(ierr); }
441     ierr = DataFieldAccessPoint(xb->field[f],pid_x, &src);CHKERRQ(ierr);
442     ierr = DataFieldAccessPoint(yb->field[f],pid_y, &dest);CHKERRQ(ierr);
443     ierr = PetscMemcpy(dest, src, xb->field[f]->atomic_size);CHKERRQ(ierr);
444     ierr = DataFieldRestoreAccess(xb->field[f]);CHKERRQ(ierr);
445     if (xb != yb) {ierr = DataFieldRestoreAccess(yb->field[f]);CHKERRQ(ierr);}
446   }
447   PetscFunctionReturn(0);
448 }
449 
450 PetscErrorCode DataBucketCreateFromSubset(DataBucket DBIn,const PetscInt N,const PetscInt list[],DataBucket *DB)
451 {
452   PetscInt nfields;
453   DataField *fields;
454   PetscInt f,L,buffer,allocated,p;
455   PetscErrorCode ierr;
456 
457   PetscFunctionBegin;
458   ierr = DataBucketCreate(DB);CHKERRQ(ierr);
459   /* copy contents of DBIn */
460   ierr = DataBucketGetDataFields(DBIn,&nfields,&fields);CHKERRQ(ierr);
461   ierr = DataBucketGetSizes(DBIn,&L,&buffer,&allocated);CHKERRQ(ierr);
462   for (f = 0; f < nfields; ++f) {
463     ierr = DataBucketRegisterField(*DB,"DataBucketCreateFromSubset",fields[f]->name,fields[f]->atomic_size,NULL);CHKERRQ(ierr);
464   }
465   ierr = DataBucketFinalize(*DB);CHKERRQ(ierr);
466   ierr = DataBucketSetSizes(*DB,L,buffer);CHKERRQ(ierr);
467   /* now copy the desired guys from DBIn => DB */
468   for (p = 0; p < N; ++p) {
469     ierr = DataBucketCopyPoint(DBIn,list[p], *DB,p);CHKERRQ(ierr);
470   }
471   PetscFunctionReturn(0);
472 }
473 
474 /* insert into an exisitng location */
475 PetscErrorCode DataFieldInsertPoint(const DataField field,const PetscInt index,const void *ctx)
476 {
477   PetscErrorCode ierr;
478 
479   PetscFunctionBegin;
480 #ifdef 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(__DATATFIELD_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 DataBucketRemovePointAtIndex(const DataBucket db,const PetscInt index)
491 {
492   PetscInt       f;
493   PetscBool      any_active_fields;
494   PetscErrorCode ierr;
495 
496   PetscFunctionBegin;
497 #ifdef 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 = DataBucketQueryForActiveFields(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 DataField 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       DataField field = db->field[f];
510 
511       /* copy then remove */
512       ierr = DataFieldCopyPoint(db->L-1, field, index, field);CHKERRQ(ierr);
513       /* DataFieldZeroPoint(field,index); */
514     }
515   }
516   /* decrement size */
517   /* this will zero out an crap at the end of the list */
518   ierr = DataBucketRemovePoint(db);CHKERRQ(ierr);
519   PetscFunctionReturn(0);
520 }
521 
522 /* copy x into y */
523 PetscErrorCode DataFieldCopyPoint(const PetscInt pid_x,const DataField field_x,
524                         const PetscInt pid_y,const DataField field_y )
525 {
526   PetscErrorCode ierr;
527 
528   PetscFunctionBegin;
529 #ifdef 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(__DATATFIELD_point_access(field_y->data,pid_y,field_y->atomic_size),
538                      __DATATFIELD_point_access(field_x->data,pid_x,field_x->atomic_size),
539                      field_y->atomic_size);CHKERRQ(ierr);
540   PetscFunctionReturn(0);
541 }
542 
543 
544 /* zero only the datafield at this point */
545 PetscErrorCode DataFieldZeroPoint(const DataField field,const PetscInt index)
546 {
547   PetscErrorCode ierr;
548 
549   PetscFunctionBegin;
550 #ifdef DATAFIELD_POINT_ACCESS_GUARD
551   /* check point is valid */
552   if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
553   if (index >= field->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",field->L);
554 #endif
555   ierr = PetscMemzero(__DATATFIELD_point_access(field->data,index,field->atomic_size), field->atomic_size);CHKERRQ(ierr);
556   PetscFunctionReturn(0);
557 }
558 
559 /* zero ALL data for this point */
560 PetscErrorCode DataBucketZeroPoint(const DataBucket db,const PetscInt index)
561 {
562   PetscInt f;
563   PetscErrorCode ierr;
564 
565   PetscFunctionBegin;
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->allocated);
569   for (f = 0; f < db->nfields; ++f) {
570     DataField field = db->field[f];
571     ierr = DataFieldZeroPoint(field,index);CHKERRQ(ierr);
572   }
573   PetscFunctionReturn(0);
574 }
575 
576 /* increment */
577 PetscErrorCode DataBucketAddPoint(DataBucket db)
578 {
579   PetscErrorCode ierr;
580 
581   PetscFunctionBegin;
582   ierr = DataBucketSetSizes(db,db->L+1,DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr);
583   PetscFunctionReturn(0);
584 }
585 
586 /* decrement */
587 PetscErrorCode DataBucketRemovePoint(DataBucket db)
588 {
589   PetscErrorCode ierr;
590 
591   PetscFunctionBegin;
592   ierr = DataBucketSetSizes(db,db->L-1,DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr);
593   PetscFunctionReturn(0);
594 }
595 
596 PetscErrorCode DataBucketView_stdout(MPI_Comm comm,DataBucket db)
597 {
598   PetscInt f;
599   double memory_usage_total,memory_usage_total_local = 0.0;
600   PetscErrorCode ierr;
601 
602   PetscFunctionBegin;
603   ierr = PetscPrintf(comm,"DataBucketView: \n");CHKERRQ(ierr);
604   ierr = PetscPrintf(comm,"  L                  = %D \n", db->L);CHKERRQ(ierr);
605   ierr = PetscPrintf(comm,"  buffer             = %D \n", db->buffer);CHKERRQ(ierr);
606   ierr = PetscPrintf(comm,"  allocated          = %D \n", db->allocated);CHKERRQ(ierr);
607   ierr = PetscPrintf(comm,"  nfields registered = %D \n", db->nfields);CHKERRQ(ierr);
608 
609   for (f = 0; f < db->nfields; ++f) {
610     double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
611     memory_usage_total_local += memory_usage_f;
612   }
613   ierr = MPI_Allreduce(&memory_usage_total_local,&memory_usage_total,1,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
614 
615   for (f = 0; f < db->nfields; ++f) {
616     double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
617     ierr = PetscPrintf(comm,"    [%3D] %15s : Mem. usage       = %1.2e (MB) [rank0]\n", f, db->field[f]->name, memory_usage_f );CHKERRQ(ierr);
618     ierr = PetscPrintf(comm,"                            blocksize        = %D \n", db->field[f]->bs);CHKERRQ(ierr);
619     if (db->field[f]->bs != 1) {
620       ierr = PetscPrintf(comm,"                            atomic size      = %zu [full block, bs=%D]\n", db->field[f]->atomic_size,db->field[f]->bs);CHKERRQ(ierr);
621       ierr = PetscPrintf(comm,"                            atomic size/item = %zu \n", db->field[f]->atomic_size/db->field[f]->bs);CHKERRQ(ierr);
622     } else {
623       ierr = PetscPrintf(comm,"                            atomic size      = %zu \n", db->field[f]->atomic_size);CHKERRQ(ierr);
624     }
625   }
626   ierr = PetscPrintf(comm,"  Total mem. usage                           = %1.2e (MB) (collective)\n", memory_usage_total);CHKERRQ(ierr);
627   PetscFunctionReturn(0);
628 }
629 
630 PetscErrorCode DataBucketView_SEQ(MPI_Comm comm,DataBucket db,const char filename[],DataBucketViewType type)
631 {
632   PetscErrorCode ierr;
633 
634   PetscFunctionBegin;
635   switch (type) {
636   case DATABUCKET_VIEW_STDOUT:
637     ierr = DataBucketView_stdout(PETSC_COMM_SELF,db);CHKERRQ(ierr);
638     break;
639   case DATABUCKET_VIEW_ASCII:
640     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for ascii output");
641     break;
642   case DATABUCKET_VIEW_BINARY:
643     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for binary output");
644     break;
645   case DATABUCKET_VIEW_HDF5:
646     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for HDF5 output");
647     break;
648   default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer method requested");
649   }
650   PetscFunctionReturn(0);
651 }
652 
653 PetscErrorCode DataBucketView_MPI(MPI_Comm comm,DataBucket db,const char filename[],DataBucketViewType type)
654 {
655   PetscErrorCode ierr;
656 
657   PetscFunctionBegin;
658   switch (type) {
659   case DATABUCKET_VIEW_STDOUT:
660     ierr = DataBucketView_stdout(comm,db);CHKERRQ(ierr);
661     break;
662   case DATABUCKET_VIEW_ASCII:
663     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for ascii output");
664     break;
665   case DATABUCKET_VIEW_BINARY:
666     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for binary output");
667     break;
668   case DATABUCKET_VIEW_HDF5:
669     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for HDF5 output");
670     break;
671   default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer method requested");
672   }
673   PetscFunctionReturn(0);
674 }
675 
676 PetscErrorCode DataBucketView(MPI_Comm comm,DataBucket db,const char filename[],DataBucketViewType type)
677 {
678   PetscMPIInt nproc;
679   PetscErrorCode ierr;
680 
681   PetscFunctionBegin;
682   ierr = MPI_Comm_size(comm,&nproc);CHKERRQ(ierr);
683   if (nproc == 1) {
684     ierr = DataBucketView_SEQ(comm,db,filename,type);CHKERRQ(ierr);
685   } else {
686     ierr = DataBucketView_MPI(comm,db,filename,type);CHKERRQ(ierr);
687   }
688   PetscFunctionReturn(0);
689 }
690 
691 PetscErrorCode DataBucketDuplicateFields(DataBucket dbA,DataBucket *dbB)
692 {
693   DataBucket db2;
694   PetscInt f;
695   PetscErrorCode ierr;
696 
697   PetscFunctionBegin;
698   ierr = DataBucketCreate(&db2);CHKERRQ(ierr);
699   /* copy contents from dbA into db2 */
700   for (f = 0; f < dbA->nfields; ++f) {
701     DataField field;
702     size_t    atomic_size;
703     char      *name;
704 
705     field = dbA->field[f];
706     atomic_size = field->atomic_size;
707     name        = field->name;
708     ierr = DataBucketRegisterField(db2,"DataBucketDuplicateFields",name,atomic_size,NULL);CHKERRQ(ierr);
709   }
710   ierr = DataBucketFinalize(db2);CHKERRQ(ierr);
711   ierr = DataBucketSetInitialSizes(db2,0,1000);CHKERRQ(ierr);
712   *dbB = db2;
713   PetscFunctionReturn(0);
714 }
715 
716 /*
717  Insert points from db2 into db1
718  db1 <<== db2
719  */
720 PetscErrorCode DataBucketInsertValues(DataBucket db1,DataBucket db2)
721 {
722   PetscInt n_mp_points1,n_mp_points2;
723   PetscInt n_mp_points1_new,p;
724   PetscErrorCode ierr;
725 
726   PetscFunctionBegin;
727   ierr = DataBucketGetSizes(db1,&n_mp_points1,0,0);CHKERRQ(ierr);
728   ierr = DataBucketGetSizes(db2,&n_mp_points2,0,0);CHKERRQ(ierr);
729   n_mp_points1_new = n_mp_points1 + n_mp_points2;
730   ierr = DataBucketSetSizes(db1,n_mp_points1_new,DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr);
731   for (p = 0; p < n_mp_points2; ++p) {
732     /* db1 <<== db2 */
733     ierr = DataBucketCopyPoint(db2,p, db1,(n_mp_points1 + p));CHKERRQ(ierr);
734   }
735   PetscFunctionReturn(0);
736 }
737 
738 /* helpers for parallel send/recv */
739 PetscErrorCode DataBucketCreatePackedArray(DataBucket db,size_t *bytes,void **buf)
740 {
741   PetscInt       f;
742   size_t         sizeof_marker_contents;
743   void          *buffer;
744   PetscErrorCode ierr;
745 
746   PetscFunctionBegin;
747   sizeof_marker_contents = 0;
748   for (f = 0; f < db->nfields; ++f) {
749     DataField df = db->field[f];
750     sizeof_marker_contents += df->atomic_size;
751   }
752   ierr = PetscMalloc(sizeof_marker_contents, &buffer);CHKERRQ(ierr);
753   ierr = PetscMemzero(buffer, sizeof_marker_contents);CHKERRQ(ierr);
754   if (bytes) {*bytes = sizeof_marker_contents;}
755   if (buf)   {*buf   = buffer;}
756   PetscFunctionReturn(0);
757 }
758 
759 PetscErrorCode DataBucketDestroyPackedArray(DataBucket db,void **buf)
760 {
761   PetscErrorCode ierr;
762 
763   PetscFunctionBegin;
764   if (buf) {
765     ierr = PetscFree(*buf);CHKERRQ(ierr);
766     *buf = NULL;
767   }
768   PetscFunctionReturn(0);
769 }
770 
771 PetscErrorCode DataBucketFillPackedArray(DataBucket db,const PetscInt index,void *buf)
772 {
773   PetscInt       f;
774   void          *data, *data_p;
775   size_t         asize, offset;
776   PetscErrorCode ierr;
777 
778   PetscFunctionBegin;
779   offset = 0;
780   for (f = 0; f < db->nfields; ++f) {
781     DataField df = db->field[f];
782 
783     asize = df->atomic_size;
784     data = (void*)( df->data );
785     data_p = (void*)( (char*)data + index*asize );
786     ierr = PetscMemcpy((void*)((char*)buf + offset), data_p, asize);CHKERRQ(ierr);
787     offset = offset + asize;
788   }
789   PetscFunctionReturn(0);
790 }
791 
792 PetscErrorCode DataBucketInsertPackedArray(DataBucket db,const PetscInt idx,void *data)
793 {
794   PetscInt f;
795   void *data_p;
796   size_t offset;
797   PetscErrorCode ierr;
798 
799   PetscFunctionBegin;
800   offset = 0;
801   for (f = 0; f < db->nfields; ++f) {
802     DataField df = db->field[f];
803 
804     data_p = (void*)( (char*)data + offset );
805     ierr = DataFieldInsertPoint(df, idx, (void*)data_p);CHKERRQ(ierr);
806     offset = offset + df->atomic_size;
807   }
808   PetscFunctionReturn(0);
809 }
810