xref: /petsc/src/dm/impls/swarm/data_bucket.c (revision 72505a4dc06f22260039f76cc0a3da7876ed95ae)
152849c42SDave May 
252849c42SDave May #include "data_bucket.h"
352849c42SDave May 
452849c42SDave May /* string helpers */
52635f519SDave May #undef __FUNCT__
62635f519SDave May #define __FUNCT__ "StringInList"
72eac95f8SDave May PetscErrorCode StringInList(const char name[],const PetscInt N,const DataField gfield[],PetscBool *val)
852849c42SDave May {
95c18a9d6SDave May   PetscInt       i;
10521f74f9SMatthew G. Knepley   PetscErrorCode ierr;
1152849c42SDave May 
12521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1352849c42SDave May   *val = PETSC_FALSE;
14521f74f9SMatthew G. Knepley   for (i = 0; i < N; ++i) {
15521f74f9SMatthew G. Knepley     PetscBool flg;
16521f74f9SMatthew G. Knepley     ierr = PetscStrcmp(name, gfield[i]->name, &flg);CHKERRQ(ierr);
17521f74f9SMatthew G. Knepley     if (flg) {
1852849c42SDave May       *val = PETSC_TRUE;
192eac95f8SDave May       PetscFunctionReturn(0);
2052849c42SDave May     }
2152849c42SDave May   }
222eac95f8SDave May   PetscFunctionReturn(0);
2352849c42SDave May }
2452849c42SDave May 
252635f519SDave May #undef __FUNCT__
262635f519SDave May #define __FUNCT__ "StringFindInList"
272eac95f8SDave May PetscErrorCode StringFindInList(const char name[],const PetscInt N,const DataField gfield[],PetscInt *index)
2852849c42SDave May {
295c18a9d6SDave May   PetscInt       i;
30521f74f9SMatthew G. Knepley   PetscErrorCode ierr;
3152849c42SDave May 
32521f74f9SMatthew G. Knepley   PetscFunctionBegin;
3352849c42SDave May   *index = -1;
34521f74f9SMatthew G. Knepley   for (i = 0; i < N; ++i) {
35521f74f9SMatthew G. Knepley     PetscBool flg;
36521f74f9SMatthew G. Knepley     ierr = PetscStrcmp(name, gfield[i]->name, &flg);CHKERRQ(ierr);
37521f74f9SMatthew G. Knepley     if (flg) {
3852849c42SDave May       *index = i;
392eac95f8SDave May       PetscFunctionReturn(0);
4052849c42SDave May     }
4152849c42SDave May   }
422eac95f8SDave May   PetscFunctionReturn(0);
4352849c42SDave May }
4452849c42SDave May 
452eac95f8SDave May #undef __FUNCT__
462eac95f8SDave May #define __FUNCT__ "DataFieldCreate"
472eac95f8SDave May PetscErrorCode DataFieldCreate(const char registeration_function[],const char name[],const size_t size,const PetscInt L,DataField *DF)
4852849c42SDave May {
4952849c42SDave May   DataField      df;
50521f74f9SMatthew G. Knepley   PetscErrorCode ierr;
5152849c42SDave May 
52521f74f9SMatthew G. Knepley   PetscFunctionBegin;
53521f74f9SMatthew G. Knepley   ierr = PetscMalloc(sizeof(struct _p_DataField), &df);CHKERRQ(ierr);
54521f74f9SMatthew G. Knepley   ierr = PetscMemzero(df, sizeof(struct _p_DataField));CHKERRQ(ierr);
55521f74f9SMatthew G. Knepley   ierr = PetscStrallocpy(registeration_function, &df->registeration_function);CHKERRQ(ierr);
56521f74f9SMatthew G. Knepley   ierr = PetscStrallocpy(name, &df->name);CHKERRQ(ierr);
5752849c42SDave May   df->atomic_size = size;
5852849c42SDave May   df->L  = L;
5952c3ed93SDave May   df->bs = 1;
60521f74f9SMatthew G. Knepley   /* allocate something so we don't have to reallocate */
61521f74f9SMatthew G. Knepley   ierr = PetscMalloc(size * L, &df->data);CHKERRQ(ierr);
62521f74f9SMatthew G. Knepley   ierr = PetscMemzero(df->data, size * L);CHKERRQ(ierr);
6352849c42SDave May   *DF = df;
642eac95f8SDave May   PetscFunctionReturn(0);
6552849c42SDave May }
6652849c42SDave May 
672eac95f8SDave May #undef __FUNCT__
682eac95f8SDave May #define __FUNCT__ "DataFieldDestroy"
692eac95f8SDave May PetscErrorCode DataFieldDestroy(DataField *DF)
7052849c42SDave May {
7152849c42SDave May   DataField      df = *DF;
72521f74f9SMatthew G. Knepley   PetscErrorCode ierr;
7352849c42SDave May 
74521f74f9SMatthew G. Knepley   PetscFunctionBegin;
75521f74f9SMatthew G. Knepley   ierr = PetscFree(df->registeration_function);CHKERRQ(ierr);
76521f74f9SMatthew G. Knepley   ierr = PetscFree(df->name);CHKERRQ(ierr);
77521f74f9SMatthew G. Knepley   ierr = PetscFree(df->data);CHKERRQ(ierr);
78521f74f9SMatthew G. Knepley   ierr = PetscFree(df);CHKERRQ(ierr);
7952849c42SDave May   *DF  = NULL;
802eac95f8SDave May   PetscFunctionReturn(0);
8152849c42SDave May }
8252849c42SDave May 
8352849c42SDave May /* data bucket */
842eac95f8SDave May #undef __FUNCT__
852eac95f8SDave May #define __FUNCT__ "DataBucketCreate"
862eac95f8SDave May PetscErrorCode DataBucketCreate(DataBucket *DB)
8752849c42SDave May {
8852849c42SDave May   DataBucket     db;
89521f74f9SMatthew G. Knepley   PetscErrorCode ierr;
9052849c42SDave May 
91521f74f9SMatthew G. Knepley   PetscFunctionBegin;
92521f74f9SMatthew G. Knepley   ierr = PetscMalloc(sizeof(struct _p_DataBucket), &db);CHKERRQ(ierr);
93521f74f9SMatthew G. Knepley   ierr = PetscMemzero(db, sizeof(struct _p_DataBucket));CHKERRQ(ierr);
9452849c42SDave May 
9552849c42SDave May   db->finalised = PETSC_FALSE;
9652849c42SDave May   /* create empty spaces for fields */
973454631fSDave May   db->L         = -1;
9852849c42SDave May   db->buffer    = 1;
9952849c42SDave May   db->allocated = 1;
10052849c42SDave May   db->nfields   = 0;
101521f74f9SMatthew G. Knepley   ierr = PetscMalloc1(1, &db->field);CHKERRQ(ierr);
10252849c42SDave May   *DB  = db;
1032eac95f8SDave May   PetscFunctionReturn(0);
10452849c42SDave May }
10552849c42SDave May 
1062eac95f8SDave May #undef __FUNCT__
1072eac95f8SDave May #define __FUNCT__ "DataBucketDestroy"
1082eac95f8SDave May PetscErrorCode DataBucketDestroy(DataBucket *DB)
10952849c42SDave May {
11052849c42SDave May   DataBucket     db = *DB;
1115c18a9d6SDave May   PetscInt       f;
112dbe06d34SDave May   PetscErrorCode ierr;
11352849c42SDave May 
114521f74f9SMatthew G. Knepley   PetscFunctionBegin;
11552849c42SDave May   /* release fields */
116521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
117dbe06d34SDave May     ierr = DataFieldDestroy(&db->field[f]);CHKERRQ(ierr);
11852849c42SDave May   }
11952849c42SDave May   /* this will catch the initially allocated objects in the event that no fields are registered */
12052849c42SDave May   if (db->field != NULL) {
121521f74f9SMatthew G. Knepley     ierr = PetscFree(db->field);CHKERRQ(ierr);
12252849c42SDave May   }
123521f74f9SMatthew G. Knepley   ierr = PetscFree(db);CHKERRQ(ierr);
12452849c42SDave May   *DB = NULL;
1252eac95f8SDave May   PetscFunctionReturn(0);
12652849c42SDave May }
12752849c42SDave May 
1282eac95f8SDave May #undef __FUNCT__
1290cb17e37SDave May #define __FUNCT__ "DataBucketQueryForActiveFields"
1300cb17e37SDave May PetscErrorCode DataBucketQueryForActiveFields(DataBucket db,PetscBool *any_active_fields)
1310cb17e37SDave May {
1320cb17e37SDave May   PetscInt f;
133521f74f9SMatthew G. Knepley 
134521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1350cb17e37SDave May   *any_active_fields = PETSC_FALSE;
136521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
1370cb17e37SDave May     if (db->field[f]->active) {
1380cb17e37SDave May       *any_active_fields = PETSC_TRUE;
1390cb17e37SDave May       PetscFunctionReturn(0);
1400cb17e37SDave May     }
1410cb17e37SDave May   }
1420cb17e37SDave May   PetscFunctionReturn(0);
1430cb17e37SDave May }
1440cb17e37SDave May 
1450cb17e37SDave May #undef __FUNCT__
1462eac95f8SDave May #define __FUNCT__ "DataBucketRegisterField"
1472eac95f8SDave May PetscErrorCode DataBucketRegisterField(
14852849c42SDave May                               DataBucket db,
14952849c42SDave May                               const char registeration_function[],
15052849c42SDave May                               const char field_name[],
15152849c42SDave May                               size_t atomic_size, DataField *_gfield)
15252849c42SDave May {
15352849c42SDave May   PetscBool val;
1544be7464cSMatthew G. Knepley   DataField fp;
155dbe06d34SDave May   PetscErrorCode ierr;
15652849c42SDave May 
157521f74f9SMatthew G. Knepley   PetscFunctionBegin;
15852849c42SDave May 	/* check we haven't finalised the registration of fields */
15952849c42SDave May 	/*
16052849c42SDave May    if(db->finalised==PETSC_TRUE) {
16152849c42SDave May    printf("ERROR: DataBucketFinalize() has been called. Cannot register more fields\n");
16252849c42SDave May    ERROR();
16352849c42SDave May    }
16452849c42SDave May    */
16552849c42SDave May   /* check for repeated name */
1662635f519SDave May   ierr = StringInList(field_name, db->nfields, (const DataField*) db->field, &val);CHKERRQ(ierr);
1672eac95f8SDave May   if (val == PETSC_TRUE) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field %s already exists. Cannot add same field twice",field_name);
16852849c42SDave May   /* create new space for data */
1694be7464cSMatthew G. Knepley   ierr = PetscRealloc(sizeof(DataField)*(db->nfields+1), &db->field);CHKERRQ(ierr);
17052849c42SDave May   /* add field */
171dbe06d34SDave May   ierr = DataFieldCreate(registeration_function, field_name, atomic_size, db->allocated, &fp);CHKERRQ(ierr);
17252849c42SDave May   db->field[db->nfields] = fp;
17352849c42SDave May   db->nfields++;
17452849c42SDave May   if (_gfield != NULL) {
17552849c42SDave May     *_gfield = fp;
17652849c42SDave May   }
1772eac95f8SDave May   PetscFunctionReturn(0);
17852849c42SDave May }
17952849c42SDave May 
18052849c42SDave May /*
18152849c42SDave May  #define DataBucketRegisterField(db,name,size,k) {\
18252849c42SDave May  char *location;\
18352849c42SDave May  asprintf(&location,"Registered by %s() at line %d within file %s", __FUNCTION__, __LINE__, __FILE__);\
18452849c42SDave May  _DataBucketRegisterField( (db), location, (name), (size), (k) );\
185521f74f9SMatthew G. Knepley  ierr = PetscFree(location);\
18652849c42SDave May  }
18752849c42SDave May  */
18852849c42SDave May 
1892eac95f8SDave May #undef __FUNCT__
1902eac95f8SDave May #define __FUNCT__ "DataBucketGetDataFieldByName"
1912eac95f8SDave May PetscErrorCode DataBucketGetDataFieldByName(DataBucket db,const char name[],DataField *gfield)
19252849c42SDave May {
1935c18a9d6SDave May   PetscInt       idx;
19452849c42SDave May   PetscBool      found;
195521f74f9SMatthew G. Knepley   PetscErrorCode ierr;
19652849c42SDave May 
197521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1982635f519SDave May   ierr = StringInList(name,db->nfields,(const DataField*)db->field,&found);CHKERRQ(ierr);
1992eac95f8SDave May   if (!found) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot find DataField with name %s",name);
2002635f519SDave May   ierr = StringFindInList(name,db->nfields,(const DataField*)db->field,&idx);CHKERRQ(ierr);
20152849c42SDave May   *gfield = db->field[idx];
2022eac95f8SDave May   PetscFunctionReturn(0);
20352849c42SDave May }
20452849c42SDave May 
2057fbf63aeSDave May #undef __FUNCT__
2067fbf63aeSDave May #define __FUNCT__ "DataBucketQueryDataFieldByName"
2077fbf63aeSDave May PetscErrorCode DataBucketQueryDataFieldByName(DataBucket db,const char name[],PetscBool *found)
20852849c42SDave May {
2092635f519SDave May   PetscErrorCode ierr;
210521f74f9SMatthew G. Knepley 
211521f74f9SMatthew G. Knepley   PetscFunctionBegin;
21252849c42SDave May   *found = PETSC_FALSE;
2132635f519SDave May   ierr = StringInList(name,db->nfields,(const DataField*)db->field,found);CHKERRQ(ierr);
2147fbf63aeSDave May   PetscFunctionReturn(0);
21552849c42SDave May }
21652849c42SDave May 
2177fbf63aeSDave May #undef __FUNCT__
2187fbf63aeSDave May #define __FUNCT__ "DataBucketFinalize"
2197fbf63aeSDave May PetscErrorCode DataBucketFinalize(DataBucket db)
22052849c42SDave May {
221521f74f9SMatthew G. Knepley   PetscFunctionBegin;
22252849c42SDave May   db->finalised = PETSC_TRUE;
2237fbf63aeSDave May   PetscFunctionReturn(0);
22452849c42SDave May }
22552849c42SDave May 
2267fbf63aeSDave May #undef __FUNCT__
2277fbf63aeSDave May #define __FUNCT__ "DataFieldGetNumEntries"
2287fbf63aeSDave May PetscErrorCode DataFieldGetNumEntries(DataField df,PetscInt *sum)
22952849c42SDave May {
230521f74f9SMatthew G. Knepley   PetscFunctionBegin;
23152849c42SDave May   *sum = df->L;
2327fbf63aeSDave May   PetscFunctionReturn(0);
23352849c42SDave May }
23452849c42SDave May 
2357fbf63aeSDave May #undef __FUNCT__
23652c3ed93SDave May #define __FUNCT__ "DataFieldSetBlockSize"
23752c3ed93SDave May PetscErrorCode DataFieldSetBlockSize(DataField df,PetscInt blocksize)
23852c3ed93SDave May {
239521f74f9SMatthew G. Knepley   PetscFunctionBegin;
24052c3ed93SDave May   df->bs = blocksize;
24152c3ed93SDave May   PetscFunctionReturn(0);
24252c3ed93SDave May }
24352c3ed93SDave May 
24452c3ed93SDave May #undef __FUNCT__
2457fbf63aeSDave May #define __FUNCT__ "DataFieldSetSize"
2467fbf63aeSDave May PetscErrorCode DataFieldSetSize(DataField df,const PetscInt new_L)
24752849c42SDave May {
248521f74f9SMatthew G. Knepley   PetscErrorCode ierr;
24952849c42SDave May 
250521f74f9SMatthew G. Knepley   PetscFunctionBegin;
251a233d522SDave May   if (new_L <= 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot set size of DataField to be <= 0");
2527fbf63aeSDave May   if (new_L == df->L) PetscFunctionReturn(0);
25352849c42SDave May   if (new_L > df->L) {
2544be7464cSMatthew G. Knepley     ierr = PetscRealloc(df->atomic_size * (new_L), &df->data);CHKERRQ(ierr);
25552849c42SDave May     /* init new contents */
256521f74f9SMatthew G. Knepley     ierr = PetscMemzero(( ((char*)df->data)+df->L*df->atomic_size), (new_L-df->L)*df->atomic_size);CHKERRQ(ierr);
2577fbf63aeSDave May   } else {
25852849c42SDave May     /* reallocate pointer list, add +1 in case new_L = 0 */
2594be7464cSMatthew G. Knepley     ierr = PetscRealloc(df->atomic_size * (new_L+1), &df->data);CHKERRQ(ierr);
26052849c42SDave May   }
26152849c42SDave May   df->L = new_L;
2627fbf63aeSDave May   PetscFunctionReturn(0);
26352849c42SDave May }
26452849c42SDave May 
2657fbf63aeSDave May #undef __FUNCT__
2667fbf63aeSDave May #define __FUNCT__ "DataFieldZeroBlock"
2677fbf63aeSDave May PetscErrorCode DataFieldZeroBlock(DataField df,const PetscInt start,const PetscInt end)
26852849c42SDave May {
269521f74f9SMatthew G. Knepley   PetscErrorCode ierr;
270521f74f9SMatthew G. Knepley 
271521f74f9SMatthew G. Knepley   PetscFunctionBegin;
2727fbf63aeSDave May   if (start > end) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if start(%D) > end(%D)",start,end);
2737fbf63aeSDave May   if (start < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if start(%D) < 0",start);
274a233d522SDave May   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);
275521f74f9SMatthew G. Knepley   ierr = PetscMemzero((((char*)df->data)+start*df->atomic_size), (end-start)*df->atomic_size);CHKERRQ(ierr);
2767fbf63aeSDave May   PetscFunctionReturn(0);
27752849c42SDave May }
27852849c42SDave May 
27952849c42SDave May /*
28052849c42SDave May  A negative buffer value will simply be ignored and the old buffer value will be used.
28152849c42SDave May  */
2827fbf63aeSDave May #undef __FUNCT__
2837fbf63aeSDave May #define __FUNCT__ "DataBucketSetSizes"
2847fbf63aeSDave May PetscErrorCode DataBucketSetSizes(DataBucket db,const PetscInt L,const PetscInt buffer)
28552849c42SDave May {
2865c18a9d6SDave May   PetscInt       current_allocated,new_used,new_unused,new_buffer,new_allocated,f;
2870cb17e37SDave May   PetscBool      any_active_fields;
288dbe06d34SDave May   PetscErrorCode ierr;
28952849c42SDave May 
290521f74f9SMatthew G. Knepley   PetscFunctionBegin;
2917fbf63aeSDave May   if (db->finalised == PETSC_FALSE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"You must call DataBucketFinalize() before DataBucketSetSizes()");
2920cb17e37SDave May   ierr = DataBucketQueryForActiveFields(db,&any_active_fields);CHKERRQ(ierr);
2930cb17e37SDave May   if (any_active_fields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot safely re-size as at least one DataField is currently being accessed");
2940cb17e37SDave May 
29552849c42SDave May   current_allocated = db->allocated;
29652849c42SDave May   new_used   = L;
29752849c42SDave May   new_unused = current_allocated - new_used;
29852849c42SDave May   new_buffer = db->buffer;
29952849c42SDave May   if (buffer >= 0) { /* update the buffer value */
30052849c42SDave May     new_buffer = buffer;
30152849c42SDave May   }
30252849c42SDave May   new_allocated = new_used + new_buffer;
30352849c42SDave May   /* action */
30452849c42SDave May   if (new_allocated > current_allocated) {
30552849c42SDave May     /* increase size to new_used + new_buffer */
30652849c42SDave May     for (f=0; f<db->nfields; f++) {
307dbe06d34SDave May       ierr = DataFieldSetSize(db->field[f], new_allocated);CHKERRQ(ierr);
30852849c42SDave May     }
30952849c42SDave May     db->L         = new_used;
31052849c42SDave May     db->buffer    = new_buffer;
31152849c42SDave May     db->allocated = new_used + new_buffer;
312521f74f9SMatthew G. Knepley   } else {
31352849c42SDave May     if (new_unused > 2 * new_buffer) {
31452849c42SDave May       /* shrink array to new_used + new_buffer */
315521f74f9SMatthew G. Knepley       for (f = 0; f < db->nfields; ++f) {
316dbe06d34SDave May         ierr = DataFieldSetSize(db->field[f], new_allocated);CHKERRQ(ierr);
31752849c42SDave May       }
31852849c42SDave May       db->L         = new_used;
31952849c42SDave May       db->buffer    = new_buffer;
32052849c42SDave May       db->allocated = new_used + new_buffer;
321521f74f9SMatthew G. Knepley     } else {
32252849c42SDave May       db->L      = new_used;
32352849c42SDave May       db->buffer = new_buffer;
32452849c42SDave May     }
32552849c42SDave May   }
32652849c42SDave May   /* zero all entries from db->L to db->allocated */
327521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
32852849c42SDave May     DataField field = db->field[f];
329dbe06d34SDave May     ierr = DataFieldZeroBlock(field, db->L,db->allocated);CHKERRQ(ierr);
33052849c42SDave May   }
3317fbf63aeSDave May   PetscFunctionReturn(0);
33252849c42SDave May }
33352849c42SDave May 
3347fbf63aeSDave May #undef __FUNCT__
3357fbf63aeSDave May #define __FUNCT__ "DataBucketSetInitialSizes"
3367fbf63aeSDave May PetscErrorCode DataBucketSetInitialSizes(DataBucket db,const PetscInt L,const PetscInt buffer)
33752849c42SDave May {
3385c18a9d6SDave May   PetscInt       f;
339dbe06d34SDave May   PetscErrorCode ierr;
340dbe06d34SDave May 
341521f74f9SMatthew G. Knepley   PetscFunctionBegin;
342dbe06d34SDave May   ierr = DataBucketSetSizes(db,L,buffer);CHKERRQ(ierr);
343521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
34452849c42SDave May     DataField field = db->field[f];
345dbe06d34SDave May     ierr = DataFieldZeroBlock(field,0,db->allocated);CHKERRQ(ierr);
34652849c42SDave May   }
3477fbf63aeSDave May   PetscFunctionReturn(0);
34852849c42SDave May }
34952849c42SDave May 
3507fbf63aeSDave May #undef __FUNCT__
3517fbf63aeSDave May #define __FUNCT__ "DataBucketGetSizes"
3527fbf63aeSDave May PetscErrorCode DataBucketGetSizes(DataBucket db,PetscInt *L,PetscInt *buffer,PetscInt *allocated)
35352849c42SDave May {
354521f74f9SMatthew G. Knepley   PetscFunctionBegin;
35552849c42SDave May   if (L) {*L = db->L;}
35652849c42SDave May   if (buffer) {*buffer = db->buffer;}
35752849c42SDave May   if (allocated) {*allocated = db->allocated;}
3587fbf63aeSDave May   PetscFunctionReturn(0);
35952849c42SDave May }
36052849c42SDave May 
3617fbf63aeSDave May #undef __FUNCT__
3627fbf63aeSDave May #define __FUNCT__ "DataBucketGetGlobalSizes"
3637fbf63aeSDave May PetscErrorCode DataBucketGetGlobalSizes(MPI_Comm comm,DataBucket db,PetscInt *L,PetscInt *buffer,PetscInt *allocated)
36452849c42SDave May {
3655c18a9d6SDave May   PetscInt _L,_buffer,_allocated;
3665c18a9d6SDave May   PetscInt ierr;
36752849c42SDave May 
368521f74f9SMatthew G. Knepley   PetscFunctionBegin;
3695c18a9d6SDave May   _L = db->L;
3705c18a9d6SDave May   _buffer = db->buffer;
3715c18a9d6SDave May   _allocated = db->allocated;
37252849c42SDave May 
37333564166SDave May   if (L) {         ierr = MPI_Allreduce(&_L,L,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); }
37433564166SDave May   if (buffer) {    ierr = MPI_Allreduce(&_buffer,buffer,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); }
37533564166SDave May   if (allocated) { ierr = MPI_Allreduce(&_allocated,allocated,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); }
3767fbf63aeSDave May   PetscFunctionReturn(0);
37752849c42SDave May }
37852849c42SDave May 
3797fbf63aeSDave May #undef __FUNCT__
38089884300SDave May #define __FUNCT__ "DataBucketGetDataFields"
3817fbf63aeSDave May PetscErrorCode DataBucketGetDataFields(DataBucket db,PetscInt *L,DataField *fields[])
38252849c42SDave May {
383521f74f9SMatthew G. Knepley   PetscFunctionBegin;
38452849c42SDave May   if (L)      {*L      = db->nfields;}
38552849c42SDave May   if (fields) {*fields = db->field;}
3867fbf63aeSDave May   PetscFunctionReturn(0);
38752849c42SDave May }
38852849c42SDave May 
3897fbf63aeSDave May #undef __FUNCT__
3907fbf63aeSDave May #define __FUNCT__ "DataFieldGetAccess"
3917fbf63aeSDave May PetscErrorCode DataFieldGetAccess(const DataField gfield)
39252849c42SDave May {
393521f74f9SMatthew G. Knepley   PetscFunctionBegin;
3947fbf63aeSDave May   if (gfield->active) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is already active. You must call DataFieldRestoreAccess()",gfield->name);
39552849c42SDave May   gfield->active = PETSC_TRUE;
3967fbf63aeSDave May   PetscFunctionReturn(0);
39752849c42SDave May }
39852849c42SDave May 
3997fbf63aeSDave May #undef __FUNCT__
4007fbf63aeSDave May #define __FUNCT__ "DataFieldAccessPoint"
4017fbf63aeSDave May PetscErrorCode DataFieldAccessPoint(const DataField gfield,const PetscInt pid,void **ctx_p)
40252849c42SDave May {
403521f74f9SMatthew G. Knepley   PetscFunctionBegin;
404*72505a4dSJed Brown   *ctx_p = NULL;
40552849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
40652849c42SDave May   /* debug mode */
40784bcda08SDave May   /* check point is valid */
4087fbf63aeSDave May   if (pid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
4097fbf63aeSDave May   if (pid >= gfield->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",gfield->L);
41084bcda08SDave May   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);
41152849c42SDave May #endif
41252849c42SDave May   *ctx_p = __DATATFIELD_point_access(gfield->data,pid,gfield->atomic_size);
4137fbf63aeSDave May   PetscFunctionReturn(0);
41452849c42SDave May }
41552849c42SDave May 
4167fbf63aeSDave May #undef __FUNCT__
4177fbf63aeSDave May #define __FUNCT__ "DataFieldAccessPointOffset"
4187fbf63aeSDave May PetscErrorCode DataFieldAccessPointOffset(const DataField gfield,const size_t offset,const PetscInt pid,void **ctx_p)
41952849c42SDave May {
420521f74f9SMatthew G. Knepley   PetscFunctionBegin;
42152849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
42252849c42SDave May   /* debug mode */
42384bcda08SDave May   /* check point is valid */
424521f74f9SMatthew G. Knepley   /* if( offset < 0 ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be >= 0");*/
425521f74f9SMatthew G. Knepley   /* Note compiler realizes this can never happen with an unsigned PetscInt */
4267fbf63aeSDave May   if (offset >= gfield->atomic_size) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be < %zu",gfield->atomic_size);
42784bcda08SDave May   /* check point is valid */
4287fbf63aeSDave May   if (pid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
429a233d522SDave May   if (pid >= gfield->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",gfield->L);
430a233d522SDave May   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);
43152849c42SDave May #endif
43252849c42SDave May   *ctx_p = __DATATFIELD_point_access_offset(gfield->data,pid,gfield->atomic_size,offset);
4337fbf63aeSDave May   PetscFunctionReturn(0);
43452849c42SDave May }
43552849c42SDave May 
4367fbf63aeSDave May #undef __FUNCT__
43789884300SDave May #define __FUNCT__ "DataFieldRestoreAccess"
4387fbf63aeSDave May PetscErrorCode DataFieldRestoreAccess(DataField gfield)
43952849c42SDave May {
440521f74f9SMatthew G. Knepley   PetscFunctionBegin;
4417fbf63aeSDave May   if (gfield->active == PETSC_FALSE) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is not active. You must call DataFieldGetAccess()", gfield->name);
44252849c42SDave May   gfield->active = PETSC_FALSE;
4437fbf63aeSDave May   PetscFunctionReturn(0);
44452849c42SDave May }
44552849c42SDave May 
4467fbf63aeSDave May #undef __FUNCT__
4477fbf63aeSDave May #define __FUNCT__ "DataFieldVerifyAccess"
4487fbf63aeSDave May PetscErrorCode DataFieldVerifyAccess(const DataField gfield,const size_t size)
44952849c42SDave May {
450521f74f9SMatthew G. Knepley   PetscFunctionBegin;
45152849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
452521f74f9SMatthew G. Knepley   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 );
45352849c42SDave May #endif
4547fbf63aeSDave May   PetscFunctionReturn(0);
45552849c42SDave May }
45652849c42SDave May 
4577fbf63aeSDave May #undef __FUNCT__
4587fbf63aeSDave May #define __FUNCT__ "DataFieldGetAtomicSize"
4597fbf63aeSDave May PetscErrorCode DataFieldGetAtomicSize(const DataField gfield,size_t *size)
46052849c42SDave May {
461521f74f9SMatthew G. Knepley   PetscFunctionBegin;
46252849c42SDave May   if (size) {*size = gfield->atomic_size;}
4637fbf63aeSDave May   PetscFunctionReturn(0);
46452849c42SDave May }
46552849c42SDave May 
4667fbf63aeSDave May #undef __FUNCT__
4677fbf63aeSDave May #define __FUNCT__ "DataFieldGetEntries"
4687fbf63aeSDave May PetscErrorCode DataFieldGetEntries(const DataField gfield,void **data)
46952849c42SDave May {
470521f74f9SMatthew G. Knepley   PetscFunctionBegin;
471521f74f9SMatthew G. Knepley   if (data) {*data = gfield->data;}
4727fbf63aeSDave May   PetscFunctionReturn(0);
47352849c42SDave May }
47452849c42SDave May 
4757fbf63aeSDave May #undef __FUNCT__
4767fbf63aeSDave May #define __FUNCT__ "DataFieldRestoreEntries"
4777fbf63aeSDave May PetscErrorCode DataFieldRestoreEntries(const DataField gfield,void **data)
47852849c42SDave May {
479521f74f9SMatthew G. Knepley   PetscFunctionBegin;
480521f74f9SMatthew G. Knepley   if (data) {*data = NULL;}
4817fbf63aeSDave May   PetscFunctionReturn(0);
48252849c42SDave May }
48352849c42SDave May 
48452849c42SDave May /* y = x */
4857fbf63aeSDave May #undef __FUNCT__
4867fbf63aeSDave May #define __FUNCT__ "DataBucketCopyPoint"
4877fbf63aeSDave May PetscErrorCode DataBucketCopyPoint(const DataBucket xb,const PetscInt pid_x,
4885c18a9d6SDave May                          const DataBucket yb,const PetscInt pid_y)
48952849c42SDave May {
4905c18a9d6SDave May   PetscInt f;
491dbe06d34SDave May   PetscErrorCode ierr;
492dbe06d34SDave May 
493521f74f9SMatthew G. Knepley   PetscFunctionBegin;
494521f74f9SMatthew G. Knepley   for (f = 0; f < xb->nfields; ++f) {
49552849c42SDave May     void *dest;
49652849c42SDave May     void *src;
49752849c42SDave May 
498dbe06d34SDave May     ierr = DataFieldGetAccess(xb->field[f]);CHKERRQ(ierr);
4996845f8f5SDave May     if (xb != yb) { ierr = DataFieldGetAccess( yb->field[f]);CHKERRQ(ierr); }
500dbe06d34SDave May     ierr = DataFieldAccessPoint(xb->field[f],pid_x, &src);CHKERRQ(ierr);
501dbe06d34SDave May     ierr = DataFieldAccessPoint(yb->field[f],pid_y, &dest);CHKERRQ(ierr);
502521f74f9SMatthew G. Knepley     ierr = PetscMemcpy(dest, src, xb->field[f]->atomic_size);CHKERRQ(ierr);
503dbe06d34SDave May     ierr = DataFieldRestoreAccess(xb->field[f]);CHKERRQ(ierr);
5046845f8f5SDave May     if (xb != yb) {ierr = DataFieldRestoreAccess(yb->field[f]);CHKERRQ(ierr);}
50552849c42SDave May   }
5067fbf63aeSDave May   PetscFunctionReturn(0);
50752849c42SDave May }
50852849c42SDave May 
5097fbf63aeSDave May #undef __FUNCT__
5107fbf63aeSDave May #define __FUNCT__ "DataBucketCreateFromSubset"
5117fbf63aeSDave May PetscErrorCode DataBucketCreateFromSubset(DataBucket DBIn,const PetscInt N,const PetscInt list[],DataBucket *DB)
51252849c42SDave May {
5135c18a9d6SDave May   PetscInt nfields;
51452849c42SDave May   DataField *fields;
5155c18a9d6SDave May   PetscInt f,L,buffer,allocated,p;
516dbe06d34SDave May   PetscErrorCode ierr;
51752849c42SDave May 
518521f74f9SMatthew G. Knepley   PetscFunctionBegin;
519521f74f9SMatthew G. Knepley   ierr = DataBucketCreate(DB);CHKERRQ(ierr);
52052849c42SDave May   /* copy contents of DBIn */
521dbe06d34SDave May   ierr = DataBucketGetDataFields(DBIn,&nfields,&fields);CHKERRQ(ierr);
522dbe06d34SDave May   ierr = DataBucketGetSizes(DBIn,&L,&buffer,&allocated);CHKERRQ(ierr);
523521f74f9SMatthew G. Knepley   for (f = 0; f < nfields; ++f) {
524dbe06d34SDave May     ierr = DataBucketRegisterField(*DB,"DataBucketCreateFromSubset",fields[f]->name,fields[f]->atomic_size,NULL);CHKERRQ(ierr);
52552849c42SDave May   }
526dbe06d34SDave May   ierr = DataBucketFinalize(*DB);CHKERRQ(ierr);
527dbe06d34SDave May   ierr = DataBucketSetSizes(*DB,L,buffer);CHKERRQ(ierr);
52852849c42SDave May   /* now copy the desired guys from DBIn => DB */
529521f74f9SMatthew G. Knepley   for (p = 0; p < N; ++p) {
530dbe06d34SDave May     ierr = DataBucketCopyPoint(DBIn,list[p], *DB,p);CHKERRQ(ierr);
53152849c42SDave May   }
5327fbf63aeSDave May   PetscFunctionReturn(0);
53352849c42SDave May }
53452849c42SDave May 
535521f74f9SMatthew G. Knepley /* insert into an exisitng location */
5367fbf63aeSDave May #undef __FUNCT__
5377fbf63aeSDave May #define __FUNCT__ "DataFieldInsertPoint"
5387fbf63aeSDave May PetscErrorCode DataFieldInsertPoint(const DataField field,const PetscInt index,const void *ctx)
53952849c42SDave May {
540521f74f9SMatthew G. Knepley   PetscErrorCode ierr;
54152849c42SDave May 
542521f74f9SMatthew G. Knepley   PetscFunctionBegin;
54352849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
54484bcda08SDave May   /* check point is valid */
545a233d522SDave May   if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
546a233d522SDave May   if (index >= field->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",field->L);
54752849c42SDave May #endif
548521f74f9SMatthew G. Knepley   ierr = PetscMemcpy(__DATATFIELD_point_access(field->data,index,field->atomic_size), ctx, field->atomic_size);CHKERRQ(ierr);
5497fbf63aeSDave May   PetscFunctionReturn(0);
55052849c42SDave May }
55152849c42SDave May 
552521f74f9SMatthew G. Knepley /* remove data at index - replace with last point */
553a233d522SDave May #undef __FUNCT__
554a233d522SDave May #define __FUNCT__ "DataBucketRemovePointAtIndex"
5557fbf63aeSDave May PetscErrorCode DataBucketRemovePointAtIndex(const DataBucket db,const PetscInt index)
55652849c42SDave May {
5575c18a9d6SDave May   PetscInt       f;
5580cb17e37SDave May   PetscBool      any_active_fields;
559dbe06d34SDave May   PetscErrorCode ierr;
56052849c42SDave May 
561521f74f9SMatthew G. Knepley   PetscFunctionBegin;
56252849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
56384bcda08SDave May   /* check point is valid */
564a233d522SDave May   if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
565a233d522SDave May   if (index >= db->allocated) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",db->L+db->buffer);
56652849c42SDave May #endif
5670cb17e37SDave May   ierr = DataBucketQueryForActiveFields(db,&any_active_fields);CHKERRQ(ierr);
5680cb17e37SDave May   if (any_active_fields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot safely remove point as at least one DataField is currently being accessed");
569a233d522SDave May   if (index >= db->L) { /* this point is not in the list - no need to error, but I will anyway */
570a233d522SDave May     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);
57152849c42SDave May   }
572a233d522SDave May   if (index != db->L-1) { /* not last point in list */
573521f74f9SMatthew G. Knepley     for (f = 0; f < db->nfields; ++f) {
57452849c42SDave May       DataField field = db->field[f];
57552849c42SDave May 
57652849c42SDave May       /* copy then remove */
577dbe06d34SDave May       ierr = DataFieldCopyPoint(db->L-1, field, index, field);CHKERRQ(ierr);
578521f74f9SMatthew G. Knepley       /* DataFieldZeroPoint(field,index); */
57952849c42SDave May     }
58052849c42SDave May   }
58152849c42SDave May   /* decrement size */
58252849c42SDave May   /* this will zero out an crap at the end of the list */
583dbe06d34SDave May   ierr = DataBucketRemovePoint(db);CHKERRQ(ierr);
5847fbf63aeSDave May   PetscFunctionReturn(0);
58552849c42SDave May }
58652849c42SDave May 
58752849c42SDave May /* copy x into y */
5887fbf63aeSDave May #undef __FUNCT__
5897fbf63aeSDave May #define __FUNCT__ "DataFieldCopyPoint"
5907fbf63aeSDave May PetscErrorCode DataFieldCopyPoint(const PetscInt pid_x,const DataField field_x,
5915c18a9d6SDave May                         const PetscInt pid_y,const DataField field_y )
59252849c42SDave May {
593521f74f9SMatthew G. Knepley   PetscErrorCode ierr;
59452849c42SDave May 
595521f74f9SMatthew G. Knepley   PetscFunctionBegin;
59652849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
59784bcda08SDave May   /* check point is valid */
598a233d522SDave May   if (pid_x < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"(IN) index must be >= 0");
599a233d522SDave May   if (pid_x >= field_x->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"(IN) index must be < %D",field_x->L);
600a233d522SDave May   if (pid_y < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"(OUT) index must be >= 0");
601a233d522SDave May   if (pid_y >= field_y->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"(OUT) index must be < %D",field_y->L);
602521f74f9SMatthew G. Knepley   if( field_y->atomic_size != field_x->atomic_size ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"atomic size must match");
60352849c42SDave May #endif
604521f74f9SMatthew G. Knepley   ierr = PetscMemcpy(__DATATFIELD_point_access(field_y->data,pid_y,field_y->atomic_size),
60552849c42SDave May                      __DATATFIELD_point_access(field_x->data,pid_x,field_x->atomic_size),
606521f74f9SMatthew G. Knepley                      field_y->atomic_size);CHKERRQ(ierr);
6077fbf63aeSDave May   PetscFunctionReturn(0);
60852849c42SDave May }
60952849c42SDave May 
61052849c42SDave May 
611521f74f9SMatthew G. Knepley /* zero only the datafield at this point */
6127fbf63aeSDave May #undef __FUNCT__
6137fbf63aeSDave May #define __FUNCT__ "DataFieldZeroPoint"
6147fbf63aeSDave May PetscErrorCode DataFieldZeroPoint(const DataField field,const PetscInt index)
61552849c42SDave May {
616521f74f9SMatthew G. Knepley   PetscErrorCode ierr;
617521f74f9SMatthew G. Knepley 
618521f74f9SMatthew G. Knepley   PetscFunctionBegin;
61952849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD
62084bcda08SDave May   /* check point is valid */
621a233d522SDave May   if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
622a233d522SDave May   if (index >= field->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",field->L);
62352849c42SDave May #endif
624521f74f9SMatthew G. Knepley   ierr = PetscMemzero(__DATATFIELD_point_access(field->data,index,field->atomic_size), field->atomic_size);CHKERRQ(ierr);
6257fbf63aeSDave May   PetscFunctionReturn(0);
62652849c42SDave May }
62752849c42SDave May 
628521f74f9SMatthew G. Knepley /* zero ALL data for this point */
6297fbf63aeSDave May #undef __FUNCT__
6307fbf63aeSDave May #define __FUNCT__ "DataBucketZeroPoint"
6317fbf63aeSDave May PetscErrorCode DataBucketZeroPoint(const DataBucket db,const PetscInt index)
63252849c42SDave May {
6335c18a9d6SDave May   PetscInt f;
634dbe06d34SDave May   PetscErrorCode ierr;
63552849c42SDave May 
636521f74f9SMatthew G. Knepley   PetscFunctionBegin;
63784bcda08SDave May   /* check point is valid */
638a233d522SDave May   if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
639a233d522SDave May   if (index >= db->allocated) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",db->allocated);
640521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
64152849c42SDave May     DataField field = db->field[f];
642dbe06d34SDave May     ierr = DataFieldZeroPoint(field,index);CHKERRQ(ierr);
64352849c42SDave May   }
6447fbf63aeSDave May   PetscFunctionReturn(0);
64552849c42SDave May }
64652849c42SDave May 
64752849c42SDave May /* increment */
6487fbf63aeSDave May #undef __FUNCT__
6497fbf63aeSDave May #define __FUNCT__ "DataBucketAddPoint"
6507fbf63aeSDave May PetscErrorCode DataBucketAddPoint(DataBucket db)
65152849c42SDave May {
652dbe06d34SDave May   PetscErrorCode ierr;
653dbe06d34SDave May 
654521f74f9SMatthew G. Knepley   PetscFunctionBegin;
655dbe06d34SDave May   ierr = DataBucketSetSizes(db,db->L+1,-1);CHKERRQ(ierr);
6567fbf63aeSDave May   PetscFunctionReturn(0);
65752849c42SDave May }
65852849c42SDave May 
6597fbf63aeSDave May /* decrement */
6607fbf63aeSDave May #undef __FUNCT__
6617fbf63aeSDave May #define __FUNCT__ "DataBucketRemovePoint"
6627fbf63aeSDave May PetscErrorCode DataBucketRemovePoint(DataBucket db)
6637fbf63aeSDave May {
664dbe06d34SDave May   PetscErrorCode ierr;
665dbe06d34SDave May 
666521f74f9SMatthew G. Knepley   PetscFunctionBegin;
667dbe06d34SDave May   ierr = DataBucketSetSizes(db,db->L-1,-1);CHKERRQ(ierr);
6687fbf63aeSDave May   PetscFunctionReturn(0);
6697fbf63aeSDave May }
6707fbf63aeSDave May 
6717fbf63aeSDave May #undef __FUNCT__
6727fbf63aeSDave May #define __FUNCT__ "_DataFieldViewBinary"
6737fbf63aeSDave May PetscErrorCode _DataFieldViewBinary(DataField field,FILE *fp)
67452849c42SDave May {
675521f74f9SMatthew G. Knepley   PetscErrorCode ierr;
67652849c42SDave May 
677521f74f9SMatthew G. Knepley   PetscFunctionBegin;
678521f74f9SMatthew G. Knepley   ierr = PetscFPrintf(PETSC_COMM_SELF, fp,"<DataField>\n");CHKERRQ(ierr);
679521f74f9SMatthew G. Knepley   ierr = PetscFPrintf(PETSC_COMM_SELF, fp,"%D\n", field->L);CHKERRQ(ierr);
680521f74f9SMatthew G. Knepley   ierr = PetscFPrintf(PETSC_COMM_SELF, fp,"%zu\n",field->atomic_size);CHKERRQ(ierr);
681521f74f9SMatthew G. Knepley   ierr = PetscFPrintf(PETSC_COMM_SELF, fp,"%s\n", field->registeration_function);CHKERRQ(ierr);
682521f74f9SMatthew G. Knepley   ierr = PetscFPrintf(PETSC_COMM_SELF, fp,"%s\n", field->name);CHKERRQ(ierr);
68352849c42SDave May   fwrite(field->data, field->atomic_size, field->L, fp);
684521f74f9SMatthew G. Knepley   ierr = PetscFPrintf(PETSC_COMM_SELF, fp,"\n</DataField>\n");CHKERRQ(ierr);
6857fbf63aeSDave May   PetscFunctionReturn(0);
68652849c42SDave May }
68752849c42SDave May 
6887fbf63aeSDave May #undef __FUNCT__
6897fbf63aeSDave May #define __FUNCT__ "_DataBucketRegisterFieldFromFile"
6907fbf63aeSDave May PetscErrorCode _DataBucketRegisterFieldFromFile(FILE *fp,DataBucket db)
69152849c42SDave May {
69252849c42SDave May   PetscBool      val;
69352849c42SDave May   DataField      gfield;
69452849c42SDave May   char           dummy[100];
69552849c42SDave May   char           registeration_function[5000];
69652849c42SDave May   char           field_name[5000];
6975c18a9d6SDave May   PetscInt       L;
69852849c42SDave May   size_t         atomic_size,strL;
699dbe06d34SDave May   PetscErrorCode ierr;
70052849c42SDave May 
701521f74f9SMatthew G. Knepley   PetscFunctionBegin;
70252849c42SDave May   /* check we haven't finalised the registration of fields */
70352849c42SDave May   /*
70452849c42SDave May    if(db->finalised==PETSC_TRUE) {
70552849c42SDave May    printf("ERROR: DataBucketFinalize() has been called. Cannot register more fields\n");
70652849c42SDave May    ERROR();
70752849c42SDave May    }
70852849c42SDave May    */
70952849c42SDave May   /* read file contents */
71099d15df1SBarry Smith   if (!fgets(dummy,99,fp))  SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Incorrect file format");
71199d15df1SBarry Smith   if (fscanf(fp, "%" PetscInt_FMT "\n",&L) != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Incorrect file format");
71299d15df1SBarry Smith   if (fscanf(fp, "%zu\n",&atomic_size) != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Incorrect file format");
71399d15df1SBarry Smith   if (!fgets(registeration_function,4999,fp))  SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Incorrect file format");
71452849c42SDave May   strL = strlen(registeration_function);
71552849c42SDave May   if (strL > 1) {
71652849c42SDave May     registeration_function[strL-1] = 0;
71752849c42SDave May   }
71899d15df1SBarry Smith   if (!fgets(field_name,4999,fp)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Incorrect file format");
71952849c42SDave May   strL = strlen(field_name);
72052849c42SDave May   if (strL > 1) {
72152849c42SDave May     field_name[strL-1] = 0;
72252849c42SDave May   }
72352849c42SDave May 
7244b46c5e1SDave May #ifdef DATA_BUCKET_LOG
725521f74f9SMatthew G. Knepley   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);
72652849c42SDave May #endif
72752849c42SDave May   /* check for repeated name */
7282635f519SDave May   ierr = StringInList( field_name, db->nfields, (const DataField*)db->field, &val );CHKERRQ(ierr);
729a233d522SDave May   if (val == PETSC_TRUE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot add same field twice");
73052849c42SDave May   /* create new space for data */
7314be7464cSMatthew G. Knepley   ierr = PetscRealloc(sizeof(DataField)*(db->nfields+1), &db->field);CHKERRQ(ierr);
73252849c42SDave May   /* add field */
733dbe06d34SDave May   ierr = DataFieldCreate( registeration_function, field_name, atomic_size, L, &gfield );CHKERRQ(ierr);
73452849c42SDave May   /* copy contents of file */
735fbe18ef9SBarry Smith   if (fread(gfield->data, gfield->atomic_size, gfield->L, fp) != (size_t) gfield->L) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Incorrect file format");
7364b46c5e1SDave May #ifdef DATA_BUCKET_LOG
737521f74f9SMatthew G. Knepley   ierr = PetscPrintf(PETSC_COMM_SELF,"  ** read %zu bytes for DataField \"%s\" \n", gfield->atomic_size * gfield->L, field_name);CHKERRQ(ierr);
73852849c42SDave May #endif
73952849c42SDave May   /* finish reading meta data */
74099d15df1SBarry Smith   if (!fgets(dummy,99,fp))  SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Incorrect file format");
74199d15df1SBarry Smith   if (!fgets(dummy,99,fp)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Incorrect file format");
74252849c42SDave May   db->field[db->nfields] = gfield;
74352849c42SDave May   db->nfields++;
7447fbf63aeSDave May   PetscFunctionReturn(0);
74552849c42SDave May }
74652849c42SDave May 
7477fbf63aeSDave May #undef __FUNCT__
7487fbf63aeSDave May #define __FUNCT__ "_DataBucketViewAscii_HeaderWrite_v00"
7497fbf63aeSDave May PetscErrorCode _DataBucketViewAscii_HeaderWrite_v00(FILE *fp)
75052849c42SDave May {
751521f74f9SMatthew G. Knepley   PetscErrorCode ierr;
752521f74f9SMatthew G. Knepley 
753521f74f9SMatthew G. Knepley   PetscFunctionBegin;
754521f74f9SMatthew G. Knepley   ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"<DataBucketHeader>\n");CHKERRQ(ierr);
755521f74f9SMatthew G. Knepley   ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"type=DataBucket\n");CHKERRQ(ierr);
756521f74f9SMatthew G. Knepley   ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"format=ascii\n");CHKERRQ(ierr);
757521f74f9SMatthew G. Knepley   ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"version=0.0\n");CHKERRQ(ierr);
758521f74f9SMatthew G. Knepley   ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"options=\n");CHKERRQ(ierr);
759521f74f9SMatthew G. Knepley   ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"</DataBucketHeader>\n");CHKERRQ(ierr);
7607fbf63aeSDave May   PetscFunctionReturn(0);
76152849c42SDave May }
7627fbf63aeSDave May 
7637fbf63aeSDave May #undef __FUNCT__
7647fbf63aeSDave May #define __FUNCT__ "_DataBucketViewAscii_HeaderRead_v00"
7657fbf63aeSDave May PetscErrorCode _DataBucketViewAscii_HeaderRead_v00(FILE *fp)
76652849c42SDave May {
76752849c42SDave May   char           dummy[100];
76852849c42SDave May   size_t         strL;
769521f74f9SMatthew G. Knepley   PetscBool      flg;
770521f74f9SMatthew G. Knepley   PetscErrorCode ierr;
77152849c42SDave May 
772521f74f9SMatthew G. Knepley   PetscFunctionBegin;
773521f74f9SMatthew G. Knepley   /* header open */
77499d15df1SBarry Smith   if (!fgets(dummy,99,fp))  SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Incorrect file format");
77552849c42SDave May 
776521f74f9SMatthew G. Knepley   /* type */
77799d15df1SBarry Smith   if (!fgets(dummy,99,fp)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Incorrect file format");
77852849c42SDave May   strL = strlen(dummy);
77952849c42SDave May   if (strL > 1) {dummy[strL-1] = 0;}
780521f74f9SMatthew G. Knepley   ierr = PetscStrcmp(dummy, "type=DataBucket", &flg);CHKERRQ(ierr);
781521f74f9SMatthew G. Knepley   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Data file doesn't contain a DataBucket type");
782521f74f9SMatthew G. Knepley   /* format */
78399d15df1SBarry Smith   if (!fgets(dummy,99,fp))  SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Incorrect file format");
784521f74f9SMatthew G. Knepley   /* version */
78599d15df1SBarry Smith   if (!fgets(dummy,99,fp))  SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Incorrect file format");
78652849c42SDave May   strL = strlen(dummy);
78752849c42SDave May   if (strL > 1) { dummy[strL-1] = 0; }
788521f74f9SMatthew G. Knepley   ierr = PetscStrcmp(dummy, "version=0.0", &flg);CHKERRQ(ierr);
789521f74f9SMatthew G. Knepley   if (!flg) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"DataBucket file must be parsed with version=0.0 : You tried %s", dummy);
790521f74f9SMatthew G. Knepley   /* options */
79199d15df1SBarry Smith   if (!fgets(dummy,99,fp))  SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Incorrect file format");
792521f74f9SMatthew G. Knepley   /* header close */
79399d15df1SBarry Smith   if (!fgets(dummy,99,fp))  SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Incorrect file format");
7947fbf63aeSDave May   PetscFunctionReturn(0);
79552849c42SDave May }
79652849c42SDave May 
7977fbf63aeSDave May #undef __FUNCT__
7987fbf63aeSDave May #define __FUNCT__ "_DataBucketLoadFromFileBinary_SEQ"
7997fbf63aeSDave May PetscErrorCode _DataBucketLoadFromFileBinary_SEQ(const char filename[],DataBucket *_db)
80052849c42SDave May {
80152849c42SDave May   DataBucket db;
80252849c42SDave May   FILE *fp;
8035c18a9d6SDave May   PetscInt L,buffer,f,nfields;
804dbe06d34SDave May   PetscErrorCode ierr;
80552849c42SDave May 
806521f74f9SMatthew G. Knepley   PetscFunctionBegin;
8074b46c5e1SDave May #ifdef DATA_BUCKET_LOG
808521f74f9SMatthew G. Knepley   ierr = PetscPrintf(PETSC_COMM_SELF,"** DataBucketLoadFromFile **\n");CHKERRQ(ierr);
80952849c42SDave May #endif
81052849c42SDave May   /* open file */
81152849c42SDave May   fp = fopen(filename,"rb");
812521f74f9SMatthew G. Knepley   if (fp == NULL) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file with name %s", filename);
81352849c42SDave May   /* read header */
814dbe06d34SDave May   ierr = _DataBucketViewAscii_HeaderRead_v00(fp);CHKERRQ(ierr);
81599d15df1SBarry Smith   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");
816dbe06d34SDave May   ierr = DataBucketCreate(&db);CHKERRQ(ierr);
817521f74f9SMatthew G. Knepley   for (f = 0; f < nfields; ++f) {
818dbe06d34SDave May     ierr = _DataBucketRegisterFieldFromFile(fp,db);CHKERRQ(ierr);
81952849c42SDave May   }
82052849c42SDave May   fclose(fp);
821dbe06d34SDave May   ierr = DataBucketFinalize(db);CHKERRQ(ierr);
82252849c42SDave May   /*
82352849c42SDave May    DataBucketSetSizes(db,L,buffer);
82452849c42SDave May    */
82552849c42SDave May   db->L = L;
82652849c42SDave May   db->buffer = buffer;
82752849c42SDave May   db->allocated = L + buffer;
82852849c42SDave May   *_db = db;
8297fbf63aeSDave May   PetscFunctionReturn(0);
83052849c42SDave May }
83152849c42SDave May 
8327fbf63aeSDave May #undef __FUNCT__
8337fbf63aeSDave May #define __FUNCT__ "DataBucketLoadFromFile"
8347fbf63aeSDave May PetscErrorCode DataBucketLoadFromFile(MPI_Comm comm,const char filename[],DataBucketViewType type,DataBucket *db)
83552849c42SDave May {
8365c18a9d6SDave May   PetscMPIInt    nproc,rank;
837997fa542SDave May   PetscErrorCode ierr;
83852849c42SDave May 
839521f74f9SMatthew G. Knepley   PetscFunctionBegin;
840997fa542SDave May   ierr = MPI_Comm_size(comm,&nproc);CHKERRQ(ierr);
841997fa542SDave May   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
8424b46c5e1SDave May #ifdef DATA_BUCKET_LOG
843521f74f9SMatthew G. Knepley   ierr = PetscPrintf(PETSC_COMM_SELF,"** DataBucketLoadFromFile **\n");CHKERRQ(ierr);
84452849c42SDave May #endif
84552849c42SDave May   if (type == DATABUCKET_VIEW_STDOUT) {
84652849c42SDave May   } else if (type == DATABUCKET_VIEW_ASCII) {
847a233d522SDave May     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot be implemented as we don't know the underlying particle data structure");
84852849c42SDave May   } else if (type == DATABUCKET_VIEW_BINARY) {
84952849c42SDave May     if (nproc == 1) {
850dbe06d34SDave May       ierr = _DataBucketLoadFromFileBinary_SEQ(filename,db);CHKERRQ(ierr);
85152849c42SDave May     } else {
852521f74f9SMatthew G. Knepley       char name[PETSC_MAX_PATH_LEN];
85352849c42SDave May 
854521f74f9SMatthew G. Knepley       ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "%s_p%1.5d", filename, rank);CHKERRQ(ierr);
855dbe06d34SDave May       ierr = _DataBucketLoadFromFileBinary_SEQ(name, db);CHKERRQ(ierr);
85652849c42SDave May     }
857521f74f9SMatthew G. Knepley   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown viewer requested");
8587fbf63aeSDave May   PetscFunctionReturn(0);
85952849c42SDave May }
86052849c42SDave May 
8617fbf63aeSDave May #undef __FUNCT__
8627fbf63aeSDave May #define __FUNCT__ "_DataBucketViewBinary"
8637fbf63aeSDave May PetscErrorCode _DataBucketViewBinary(DataBucket db,const char filename[])
86452849c42SDave May {
86552849c42SDave May   FILE          *fp = NULL;
8665c18a9d6SDave May   PetscInt       f;
867dbe06d34SDave May   PetscErrorCode ierr;
86852849c42SDave May 
869521f74f9SMatthew G. Knepley   PetscFunctionBegin;
87052849c42SDave May   fp = fopen(filename,"wb");
871a233d522SDave May   if (fp == NULL) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot open file with name %s", filename);
87252849c42SDave May   /* db header */
873dbe06d34SDave May   ierr =_DataBucketViewAscii_HeaderWrite_v00(fp);CHKERRQ(ierr);
87452849c42SDave May   /* meta-data */
875521f74f9SMatthew G. Knepley   ierr = PetscFPrintf(PETSC_COMM_SELF, fp, "%D\n%D\n%D\n", db->L,db->buffer,db->nfields);CHKERRQ(ierr);
87652849c42SDave May   /* load datafields */
877521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
878dbe06d34SDave May     ierr = _DataFieldViewBinary(db->field[f],fp);CHKERRQ(ierr);
87952849c42SDave May   }
88052849c42SDave May   fclose(fp);
8817fbf63aeSDave May   PetscFunctionReturn(0);
88252849c42SDave May }
88352849c42SDave May 
8847fbf63aeSDave May #undef __FUNCT__
8857fbf63aeSDave May #define __FUNCT__ "DataBucketView_SEQ"
8867fbf63aeSDave May PetscErrorCode DataBucketView_SEQ(DataBucket db,const char filename[],DataBucketViewType type)
88752849c42SDave May {
888dbe06d34SDave May   PetscErrorCode ierr;
889dbe06d34SDave May 
890521f74f9SMatthew G. Knepley   PetscFunctionBegin;
89152849c42SDave May   switch (type) {
89252849c42SDave May   case DATABUCKET_VIEW_STDOUT:
89352849c42SDave May   {
8945c18a9d6SDave May     PetscInt f;
89552849c42SDave May     double memory_usage_total = 0.0;
89652849c42SDave May 
897521f74f9SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_SELF,"DataBucketView(SEQ): (\"%s\")\n",filename);CHKERRQ(ierr);
898521f74f9SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_SELF,"  L                  = %D \n", db->L);CHKERRQ(ierr);
899521f74f9SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_SELF,"  buffer             = %D \n", db->buffer);CHKERRQ(ierr);
900521f74f9SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_SELF,"  allocated          = %D \n", db->allocated);CHKERRQ(ierr);
901521f74f9SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_SELF,"  nfields registered = %D \n", db->nfields);CHKERRQ(ierr);
902521f74f9SMatthew G. Knepley     for (f = 0; f < db->nfields; ++f) {
90352849c42SDave May       double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
904521f74f9SMatthew G. Knepley       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);
905521f74f9SMatthew G. Knepley       ierr = PetscPrintf(PETSC_COMM_SELF,"           blocksize          = %D \n", db->field[f]->bs);CHKERRQ(ierr);
90652c3ed93SDave May       if (db->field[f]->bs != 1) {
907521f74f9SMatthew G. Knepley         ierr = PetscPrintf(PETSC_COMM_SELF,"           atomic size        = %zu [full block, bs=%D]\n", db->field[f]->atomic_size,db->field[f]->bs);CHKERRQ(ierr);
908521f74f9SMatthew G. Knepley         ierr = PetscPrintf(PETSC_COMM_SELF,"           atomic size/item   = %zu \n", db->field[f]->atomic_size/db->field[f]->bs);CHKERRQ(ierr);
90952c3ed93SDave May       } else {
910521f74f9SMatthew G. Knepley         ierr = PetscPrintf(PETSC_COMM_SELF,"           atomic size        = %zu \n", db->field[f]->atomic_size);CHKERRQ(ierr);
91152c3ed93SDave May       }
91252849c42SDave May       memory_usage_total += memory_usage_f;
91352849c42SDave May     }
914521f74f9SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_SELF,"  Total mem. usage                                                      = %1.2e (MB) \n", memory_usage_total);CHKERRQ(ierr);
91552849c42SDave May   }
91652849c42SDave May   break;
91752849c42SDave May   case DATABUCKET_VIEW_ASCII:
91852849c42SDave May   {
919a233d522SDave May     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot be implemented as we don't know the underlying particle data structure");
92052849c42SDave May   }
92152849c42SDave May   break;
92252849c42SDave May   case DATABUCKET_VIEW_BINARY:
92352849c42SDave May   {
924dbe06d34SDave May     ierr = _DataBucketViewBinary(db,filename);CHKERRQ(ierr);
92552849c42SDave May   }
92652849c42SDave May   break;
92752849c42SDave May   case DATABUCKET_VIEW_HDF5:
92852849c42SDave May   {
929a233d522SDave May     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No HDF5 support");
93052849c42SDave May   }
93152849c42SDave May   break;
932521f74f9SMatthew G. Knepley   default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer method requested");
93352849c42SDave May   }
9347fbf63aeSDave May   PetscFunctionReturn(0);
93552849c42SDave May }
93652849c42SDave May 
9377fbf63aeSDave May #undef __FUNCT__
9387fbf63aeSDave May #define __FUNCT__ "DataBucketView_MPI"
9397fbf63aeSDave May PetscErrorCode DataBucketView_MPI(MPI_Comm comm,DataBucket db,const char filename[],DataBucketViewType type)
94052849c42SDave May {
941997fa542SDave May   PetscErrorCode ierr;
942997fa542SDave May 
943521f74f9SMatthew G. Knepley   PetscFunctionBegin;
94452849c42SDave May   switch (type) {
94552849c42SDave May   case DATABUCKET_VIEW_STDOUT:
94652849c42SDave May   {
9475c18a9d6SDave May     PetscInt f;
9485c18a9d6SDave May     PetscInt L,buffer,allocated;
94952849c42SDave May     double memory_usage_total,memory_usage_total_local = 0.0;
9505c18a9d6SDave May     PetscMPIInt rank;
95152849c42SDave May 
952997fa542SDave May     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
953dbe06d34SDave May     ierr = DataBucketGetGlobalSizes(comm,db,&L,&buffer,&allocated);CHKERRQ(ierr);
954521f74f9SMatthew G. Knepley     for (f = 0; f < db->nfields; ++f) {
95552849c42SDave May       double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
95652849c42SDave May       memory_usage_total_local += memory_usage_f;
95752849c42SDave May     }
958dbe06d34SDave May     ierr = MPI_Allreduce(&memory_usage_total_local,&memory_usage_total,1,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
959521f74f9SMatthew G. Knepley     ierr = PetscPrintf(comm,"DataBucketView(MPI): (\"%s\")\n",filename);CHKERRQ(ierr);
960521f74f9SMatthew G. Knepley     ierr = PetscPrintf(comm,"  L                  = %D \n", L);CHKERRQ(ierr);
961521f74f9SMatthew G. Knepley     ierr = PetscPrintf(comm,"  buffer (max)       = %D \n", buffer);CHKERRQ(ierr);
962521f74f9SMatthew G. Knepley     ierr = PetscPrintf(comm,"  allocated          = %D \n", allocated);CHKERRQ(ierr);
963521f74f9SMatthew G. Knepley     ierr = PetscPrintf(comm,"  nfields registered = %D \n", db->nfields);CHKERRQ(ierr);
96452849c42SDave May     for (f=0; f<db->nfields; f++) {
96552849c42SDave May       double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
966521f74f9SMatthew G. Knepley       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);
96752849c42SDave May     }
968521f74f9SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_SELF,"  Total mem. usage                                                      = %1.2e (MB) : collective\n", memory_usage_total);CHKERRQ(ierr);
96952849c42SDave May   }
97052849c42SDave May   break;
97152849c42SDave May   case DATABUCKET_VIEW_ASCII:
97252849c42SDave May   {
973a233d522SDave May     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot be implemented as we don't know the underlying data structure");
97452849c42SDave May   }
97552849c42SDave May   break;
97652849c42SDave May   case DATABUCKET_VIEW_BINARY:
97752849c42SDave May   {
978521f74f9SMatthew G. Knepley     char        name[PETSC_MAX_PATH_LEN];
9795c18a9d6SDave May     PetscMPIInt rank;
98052849c42SDave May 
98152849c42SDave May     /* create correct extension */
982997fa542SDave May     ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
983521f74f9SMatthew G. Knepley     ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "%s_p%1.5d", filename, rank);CHKERRQ(ierr);
984dbe06d34SDave May     ierr = _DataBucketViewBinary(db, name);CHKERRQ(ierr);
98552849c42SDave May   }
98652849c42SDave May   break;
98752849c42SDave May   case DATABUCKET_VIEW_HDF5:
98852849c42SDave May   {
989a233d522SDave May     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for HDF5");
99052849c42SDave May   }
99152849c42SDave May   break;
992521f74f9SMatthew G. Knepley   default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer method requested");
99352849c42SDave May   }
9947fbf63aeSDave May   PetscFunctionReturn(0);
99552849c42SDave May }
99652849c42SDave May 
9977fbf63aeSDave May #undef __FUNCT__
9987fbf63aeSDave May #define __FUNCT__ "DataBucketView"
9997fbf63aeSDave May PetscErrorCode DataBucketView(MPI_Comm comm,DataBucket db,const char filename[],DataBucketViewType type)
100052849c42SDave May {
10015c18a9d6SDave May   PetscMPIInt nproc;
1002997fa542SDave May   PetscErrorCode ierr;
100352849c42SDave May 
1004521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1005997fa542SDave May   ierr = MPI_Comm_size(comm,&nproc);CHKERRQ(ierr);
100652849c42SDave May   if (nproc == 1) {
1007dbe06d34SDave May     ierr = DataBucketView_SEQ(db,filename,type);CHKERRQ(ierr);
100852849c42SDave May   } else {
1009dbe06d34SDave May     ierr = DataBucketView_MPI(comm,db,filename,type);CHKERRQ(ierr);
101052849c42SDave May   }
10117fbf63aeSDave May   PetscFunctionReturn(0);
101252849c42SDave May }
101352849c42SDave May 
10147fbf63aeSDave May #undef __FUNCT__
10157fbf63aeSDave May #define __FUNCT__ "DataBucketDuplicateFields"
10167fbf63aeSDave May PetscErrorCode DataBucketDuplicateFields(DataBucket dbA,DataBucket *dbB)
101752849c42SDave May {
101852849c42SDave May   DataBucket db2;
10195c18a9d6SDave May   PetscInt f;
1020dbe06d34SDave May   PetscErrorCode ierr;
102152849c42SDave May 
1022521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1023dbe06d34SDave May   ierr = DataBucketCreate(&db2);CHKERRQ(ierr);
102452849c42SDave May   /* copy contents from dbA into db2 */
1025521f74f9SMatthew G. Knepley   for (f = 0; f < dbA->nfields; ++f) {
102652849c42SDave May     DataField field;
102752849c42SDave May     size_t    atomic_size;
102852849c42SDave May     char      *name;
102952849c42SDave May 
103052849c42SDave May     field = dbA->field[f];
103152849c42SDave May     atomic_size = field->atomic_size;
103252849c42SDave May     name        = field->name;
1033dbe06d34SDave May     ierr = DataBucketRegisterField(db2,"DataBucketDuplicateFields",name,atomic_size,NULL);CHKERRQ(ierr);
103452849c42SDave May   }
1035dbe06d34SDave May   ierr = DataBucketFinalize(db2);CHKERRQ(ierr);
1036dbe06d34SDave May   ierr = DataBucketSetInitialSizes(db2,0,1000);CHKERRQ(ierr);
103752849c42SDave May   *dbB = db2;
10387fbf63aeSDave May   PetscFunctionReturn(0);
103952849c42SDave May }
104052849c42SDave May 
104152849c42SDave May /*
104252849c42SDave May  Insert points from db2 into db1
104352849c42SDave May  db1 <<== db2
104452849c42SDave May  */
10457fbf63aeSDave May #undef __FUNCT__
10467fbf63aeSDave May #define __FUNCT__ "DataBucketInsertValues"
10477fbf63aeSDave May PetscErrorCode DataBucketInsertValues(DataBucket db1,DataBucket db2)
104852849c42SDave May {
10495c18a9d6SDave May   PetscInt n_mp_points1,n_mp_points2;
10505c18a9d6SDave May   PetscInt n_mp_points1_new,p;
1051dbe06d34SDave May   PetscErrorCode ierr;
105252849c42SDave May 
1053521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1054dbe06d34SDave May   ierr = DataBucketGetSizes(db1,&n_mp_points1,0,0);CHKERRQ(ierr);
1055dbe06d34SDave May   ierr = DataBucketGetSizes(db2,&n_mp_points2,0,0);CHKERRQ(ierr);
105652849c42SDave May   n_mp_points1_new = n_mp_points1 + n_mp_points2;
1057dbe06d34SDave May   ierr = DataBucketSetSizes(db1,n_mp_points1_new,-1);CHKERRQ(ierr);
1058521f74f9SMatthew G. Knepley   for (p = 0; p < n_mp_points2; ++p) {
1059521f74f9SMatthew G. Knepley     /* db1 <<== db2 */
1060dbe06d34SDave May     ierr = DataBucketCopyPoint(db2,p, db1,(n_mp_points1 + p));CHKERRQ(ierr);
106152849c42SDave May   }
10627fbf63aeSDave May   PetscFunctionReturn(0);
106352849c42SDave May }
106452849c42SDave May 
106552849c42SDave May /* helpers for parallel send/recv */
10667fbf63aeSDave May #undef __FUNCT__
10677fbf63aeSDave May #define __FUNCT__ "DataBucketCreatePackedArray"
10687fbf63aeSDave May PetscErrorCode DataBucketCreatePackedArray(DataBucket db,size_t *bytes,void **buf)
106952849c42SDave May {
10705c18a9d6SDave May   PetscInt       f;
107152849c42SDave May   size_t         sizeof_marker_contents;
107252849c42SDave May   void          *buffer;
1073521f74f9SMatthew G. Knepley   PetscErrorCode ierr;
107452849c42SDave May 
1075521f74f9SMatthew G. Knepley   PetscFunctionBegin;
107652849c42SDave May   sizeof_marker_contents = 0;
1077521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
107852849c42SDave May     DataField df = db->field[f];
107952849c42SDave May     sizeof_marker_contents += df->atomic_size;
108052849c42SDave May   }
1081521f74f9SMatthew G. Knepley   ierr = PetscMalloc(sizeof_marker_contents, &buffer);CHKERRQ(ierr);
1082521f74f9SMatthew G. Knepley   ierr = PetscMemzero(buffer, sizeof_marker_contents);CHKERRQ(ierr);
108352849c42SDave May   if (bytes) {*bytes = sizeof_marker_contents;}
108452849c42SDave May   if (buf)   {*buf   = buffer;}
10857fbf63aeSDave May   PetscFunctionReturn(0);
108652849c42SDave May }
108752849c42SDave May 
10887fbf63aeSDave May #undef __FUNCT__
10897fbf63aeSDave May #define __FUNCT__ "DataBucketDestroyPackedArray"
10907fbf63aeSDave May PetscErrorCode DataBucketDestroyPackedArray(DataBucket db,void **buf)
109152849c42SDave May {
1092521f74f9SMatthew G. Knepley   PetscErrorCode ierr;
1093521f74f9SMatthew G. Knepley 
1094521f74f9SMatthew G. Knepley   PetscFunctionBegin;
109552849c42SDave May   if (buf) {
1096521f74f9SMatthew G. Knepley     ierr = PetscFree(*buf);CHKERRQ(ierr);
109752849c42SDave May     *buf = NULL;
109852849c42SDave May   }
10997fbf63aeSDave May   PetscFunctionReturn(0);
110052849c42SDave May }
110152849c42SDave May 
11027fbf63aeSDave May #undef __FUNCT__
11037fbf63aeSDave May #define __FUNCT__ "DataBucketFillPackedArray"
11047fbf63aeSDave May PetscErrorCode DataBucketFillPackedArray(DataBucket db,const PetscInt index,void *buf)
110552849c42SDave May {
11065c18a9d6SDave May   PetscInt       f;
110752849c42SDave May   void          *data, *data_p;
110852849c42SDave May   size_t         asize, offset;
1109521f74f9SMatthew G. Knepley   PetscErrorCode ierr;
111052849c42SDave May 
1111521f74f9SMatthew G. Knepley   PetscFunctionBegin;
111252849c42SDave May   offset = 0;
1113521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
111452849c42SDave May     DataField df = db->field[f];
111552849c42SDave May 
111652849c42SDave May     asize = df->atomic_size;
111752849c42SDave May     data = (void*)( df->data );
111852849c42SDave May     data_p = (void*)( (char*)data + index*asize );
1119521f74f9SMatthew G. Knepley     ierr = PetscMemcpy((void*)((char*)buf + offset), data_p, asize);CHKERRQ(ierr);
112052849c42SDave May     offset = offset + asize;
112152849c42SDave May   }
11227fbf63aeSDave May   PetscFunctionReturn(0);
112352849c42SDave May }
112452849c42SDave May 
11257fbf63aeSDave May #undef __FUNCT__
11267fbf63aeSDave May #define __FUNCT__ "DataBucketInsertPackedArray"
11277fbf63aeSDave May PetscErrorCode DataBucketInsertPackedArray(DataBucket db,const PetscInt idx,void *data)
112852849c42SDave May {
11295c18a9d6SDave May   PetscInt f;
113052849c42SDave May   void *data_p;
113152849c42SDave May   size_t offset;
1132dbe06d34SDave May   PetscErrorCode ierr;
113352849c42SDave May 
1134521f74f9SMatthew G. Knepley   PetscFunctionBegin;
113552849c42SDave May   offset = 0;
1136521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
113752849c42SDave May     DataField df = db->field[f];
113852849c42SDave May 
113952849c42SDave May     data_p = (void*)( (char*)data + offset );
1140dbe06d34SDave May     ierr = DataFieldInsertPoint(df, idx, (void*)data_p);CHKERRQ(ierr);
114152849c42SDave May     offset = offset + df->atomic_size;
114252849c42SDave May   }
11437fbf63aeSDave May   PetscFunctionReturn(0);
114452849c42SDave May }
1145