152849c42SDave May #include "data_bucket.h" 252849c42SDave May 352849c42SDave May /* string helpers */ 42eac95f8SDave May PetscErrorCode StringInList(const char name[],const PetscInt N,const DataField gfield[],PetscBool *val) 552849c42SDave May { 65c18a9d6SDave May PetscInt i; 7521f74f9SMatthew G. Knepley PetscErrorCode ierr; 852849c42SDave May 9521f74f9SMatthew G. Knepley PetscFunctionBegin; 1052849c42SDave May *val = PETSC_FALSE; 11521f74f9SMatthew G. Knepley for (i = 0; i < N; ++i) { 12521f74f9SMatthew G. Knepley PetscBool flg; 13521f74f9SMatthew G. Knepley ierr = PetscStrcmp(name, gfield[i]->name, &flg);CHKERRQ(ierr); 14521f74f9SMatthew G. Knepley if (flg) { 1552849c42SDave May *val = PETSC_TRUE; 162eac95f8SDave May PetscFunctionReturn(0); 1752849c42SDave May } 1852849c42SDave May } 192eac95f8SDave May PetscFunctionReturn(0); 2052849c42SDave May } 2152849c42SDave May 222eac95f8SDave May PetscErrorCode StringFindInList(const char name[],const PetscInt N,const DataField gfield[],PetscInt *index) 2352849c42SDave May { 245c18a9d6SDave May PetscInt i; 25521f74f9SMatthew G. Knepley PetscErrorCode ierr; 2652849c42SDave May 27521f74f9SMatthew G. Knepley PetscFunctionBegin; 2852849c42SDave May *index = -1; 29521f74f9SMatthew G. Knepley for (i = 0; i < N; ++i) { 30521f74f9SMatthew G. Knepley PetscBool flg; 31521f74f9SMatthew G. Knepley ierr = PetscStrcmp(name, gfield[i]->name, &flg);CHKERRQ(ierr); 32521f74f9SMatthew G. Knepley if (flg) { 3352849c42SDave May *index = i; 342eac95f8SDave May PetscFunctionReturn(0); 3552849c42SDave May } 3652849c42SDave May } 372eac95f8SDave May PetscFunctionReturn(0); 3852849c42SDave May } 3952849c42SDave May 402eac95f8SDave May PetscErrorCode DataFieldCreate(const char registeration_function[],const char name[],const size_t size,const PetscInt L,DataField *DF) 4152849c42SDave May { 4252849c42SDave May DataField df; 43521f74f9SMatthew G. Knepley PetscErrorCode ierr; 4452849c42SDave May 45521f74f9SMatthew G. Knepley PetscFunctionBegin; 46521f74f9SMatthew G. Knepley ierr = PetscMalloc(sizeof(struct _p_DataField), &df);CHKERRQ(ierr); 47521f74f9SMatthew G. Knepley ierr = PetscMemzero(df, sizeof(struct _p_DataField));CHKERRQ(ierr); 48521f74f9SMatthew G. Knepley ierr = PetscStrallocpy(registeration_function, &df->registeration_function);CHKERRQ(ierr); 49521f74f9SMatthew G. Knepley ierr = PetscStrallocpy(name, &df->name);CHKERRQ(ierr); 5052849c42SDave May df->atomic_size = size; 5152849c42SDave May df->L = L; 5252c3ed93SDave May df->bs = 1; 53521f74f9SMatthew G. Knepley /* allocate something so we don't have to reallocate */ 54521f74f9SMatthew G. Knepley ierr = PetscMalloc(size * L, &df->data);CHKERRQ(ierr); 55521f74f9SMatthew G. Knepley ierr = PetscMemzero(df->data, size * L);CHKERRQ(ierr); 5652849c42SDave May *DF = df; 572eac95f8SDave May PetscFunctionReturn(0); 5852849c42SDave May } 5952849c42SDave May 602eac95f8SDave May PetscErrorCode DataFieldDestroy(DataField *DF) 6152849c42SDave May { 6252849c42SDave May DataField df = *DF; 63521f74f9SMatthew G. Knepley PetscErrorCode ierr; 6452849c42SDave May 65521f74f9SMatthew G. Knepley PetscFunctionBegin; 66521f74f9SMatthew G. Knepley ierr = PetscFree(df->registeration_function);CHKERRQ(ierr); 67521f74f9SMatthew G. Knepley ierr = PetscFree(df->name);CHKERRQ(ierr); 68521f74f9SMatthew G. Knepley ierr = PetscFree(df->data);CHKERRQ(ierr); 69521f74f9SMatthew G. Knepley ierr = PetscFree(df);CHKERRQ(ierr); 7052849c42SDave May *DF = NULL; 712eac95f8SDave May PetscFunctionReturn(0); 7252849c42SDave May } 7352849c42SDave May 7452849c42SDave May /* data bucket */ 752eac95f8SDave May PetscErrorCode DataBucketCreate(DataBucket *DB) 7652849c42SDave May { 7752849c42SDave May DataBucket db; 78521f74f9SMatthew G. Knepley PetscErrorCode ierr; 7952849c42SDave May 80521f74f9SMatthew G. Knepley PetscFunctionBegin; 81521f74f9SMatthew G. Knepley ierr = PetscMalloc(sizeof(struct _p_DataBucket), &db);CHKERRQ(ierr); 82521f74f9SMatthew G. Knepley ierr = PetscMemzero(db, sizeof(struct _p_DataBucket));CHKERRQ(ierr); 8352849c42SDave May 8452849c42SDave May db->finalised = PETSC_FALSE; 8552849c42SDave May /* create empty spaces for fields */ 863454631fSDave May db->L = -1; 8752849c42SDave May db->buffer = 1; 8852849c42SDave May db->allocated = 1; 8952849c42SDave May db->nfields = 0; 90521f74f9SMatthew G. Knepley ierr = PetscMalloc1(1, &db->field);CHKERRQ(ierr); 9152849c42SDave May *DB = db; 922eac95f8SDave May PetscFunctionReturn(0); 9352849c42SDave May } 9452849c42SDave May 952eac95f8SDave May PetscErrorCode DataBucketDestroy(DataBucket *DB) 9652849c42SDave May { 9752849c42SDave May DataBucket db = *DB; 985c18a9d6SDave May PetscInt f; 99dbe06d34SDave May PetscErrorCode ierr; 10052849c42SDave May 101521f74f9SMatthew G. Knepley PetscFunctionBegin; 10252849c42SDave May /* release fields */ 103521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 104dbe06d34SDave May ierr = DataFieldDestroy(&db->field[f]);CHKERRQ(ierr); 10552849c42SDave May } 10652849c42SDave May /* this will catch the initially allocated objects in the event that no fields are registered */ 10752849c42SDave May if (db->field != NULL) { 108521f74f9SMatthew G. Knepley ierr = PetscFree(db->field);CHKERRQ(ierr); 10952849c42SDave May } 110521f74f9SMatthew G. Knepley ierr = PetscFree(db);CHKERRQ(ierr); 11152849c42SDave May *DB = NULL; 1122eac95f8SDave May PetscFunctionReturn(0); 11352849c42SDave May } 11452849c42SDave May 1150cb17e37SDave May PetscErrorCode DataBucketQueryForActiveFields(DataBucket db,PetscBool *any_active_fields) 1160cb17e37SDave May { 1170cb17e37SDave May PetscInt f; 118521f74f9SMatthew G. Knepley 119521f74f9SMatthew G. Knepley PetscFunctionBegin; 1200cb17e37SDave May *any_active_fields = PETSC_FALSE; 121521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 1220cb17e37SDave May if (db->field[f]->active) { 1230cb17e37SDave May *any_active_fields = PETSC_TRUE; 1240cb17e37SDave May PetscFunctionReturn(0); 1250cb17e37SDave May } 1260cb17e37SDave May } 1270cb17e37SDave May PetscFunctionReturn(0); 1280cb17e37SDave May } 1290cb17e37SDave May 1302eac95f8SDave May PetscErrorCode DataBucketRegisterField( 13152849c42SDave May DataBucket db, 13252849c42SDave May const char registeration_function[], 13352849c42SDave May const char field_name[], 13452849c42SDave May size_t atomic_size, DataField *_gfield) 13552849c42SDave May { 13652849c42SDave May PetscBool val; 1374be7464cSMatthew G. Knepley DataField fp; 138dbe06d34SDave May PetscErrorCode ierr; 13952849c42SDave May 140521f74f9SMatthew G. Knepley PetscFunctionBegin; 14152849c42SDave May /* check we haven't finalised the registration of fields */ 14252849c42SDave May /* 14352849c42SDave May if(db->finalised==PETSC_TRUE) { 14452849c42SDave May printf("ERROR: DataBucketFinalize() has been called. Cannot register more fields\n"); 14552849c42SDave May ERROR(); 14652849c42SDave May } 14752849c42SDave May */ 14852849c42SDave May /* check for repeated name */ 1492635f519SDave May ierr = StringInList(field_name, db->nfields, (const DataField*) db->field, &val);CHKERRQ(ierr); 1502eac95f8SDave May if (val == PETSC_TRUE) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field %s already exists. Cannot add same field twice",field_name); 15152849c42SDave May /* create new space for data */ 1524be7464cSMatthew G. Knepley ierr = PetscRealloc(sizeof(DataField)*(db->nfields+1), &db->field);CHKERRQ(ierr); 15352849c42SDave May /* add field */ 154dbe06d34SDave May ierr = DataFieldCreate(registeration_function, field_name, atomic_size, db->allocated, &fp);CHKERRQ(ierr); 15552849c42SDave May db->field[db->nfields] = fp; 15652849c42SDave May db->nfields++; 15752849c42SDave May if (_gfield != NULL) { 15852849c42SDave May *_gfield = fp; 15952849c42SDave May } 1602eac95f8SDave May PetscFunctionReturn(0); 16152849c42SDave May } 16252849c42SDave May 16352849c42SDave May /* 16452849c42SDave May #define DataBucketRegisterField(db,name,size,k) {\ 16552849c42SDave May char *location;\ 16652849c42SDave May asprintf(&location,"Registered by %s() at line %d within file %s", __FUNCTION__, __LINE__, __FILE__);\ 16752849c42SDave May _DataBucketRegisterField( (db), location, (name), (size), (k) );\ 168521f74f9SMatthew G. Knepley ierr = PetscFree(location);\ 16952849c42SDave May } 17052849c42SDave May */ 17152849c42SDave May 1722eac95f8SDave May PetscErrorCode DataBucketGetDataFieldByName(DataBucket db,const char name[],DataField *gfield) 17352849c42SDave May { 1745c18a9d6SDave May PetscInt idx; 17552849c42SDave May PetscBool found; 176521f74f9SMatthew G. Knepley PetscErrorCode ierr; 17752849c42SDave May 178521f74f9SMatthew G. Knepley PetscFunctionBegin; 1792635f519SDave May ierr = StringInList(name,db->nfields,(const DataField*)db->field,&found);CHKERRQ(ierr); 1802eac95f8SDave May if (!found) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot find DataField with name %s",name); 1812635f519SDave May ierr = StringFindInList(name,db->nfields,(const DataField*)db->field,&idx);CHKERRQ(ierr); 18252849c42SDave May *gfield = db->field[idx]; 1832eac95f8SDave May PetscFunctionReturn(0); 18452849c42SDave May } 18552849c42SDave May 1867fbf63aeSDave May PetscErrorCode DataBucketQueryDataFieldByName(DataBucket db,const char name[],PetscBool *found) 18752849c42SDave May { 1882635f519SDave May PetscErrorCode ierr; 189521f74f9SMatthew G. Knepley 190521f74f9SMatthew G. Knepley PetscFunctionBegin; 19152849c42SDave May *found = PETSC_FALSE; 1922635f519SDave May ierr = StringInList(name,db->nfields,(const DataField*)db->field,found);CHKERRQ(ierr); 1937fbf63aeSDave May PetscFunctionReturn(0); 19452849c42SDave May } 19552849c42SDave May 1967fbf63aeSDave May PetscErrorCode DataBucketFinalize(DataBucket db) 19752849c42SDave May { 198521f74f9SMatthew G. Knepley PetscFunctionBegin; 19952849c42SDave May db->finalised = PETSC_TRUE; 2007fbf63aeSDave May PetscFunctionReturn(0); 20152849c42SDave May } 20252849c42SDave May 2037fbf63aeSDave May PetscErrorCode DataFieldGetNumEntries(DataField df,PetscInt *sum) 20452849c42SDave May { 205521f74f9SMatthew G. Knepley PetscFunctionBegin; 20652849c42SDave May *sum = df->L; 2077fbf63aeSDave May PetscFunctionReturn(0); 20852849c42SDave May } 20952849c42SDave May 21052c3ed93SDave May PetscErrorCode DataFieldSetBlockSize(DataField df,PetscInt blocksize) 21152c3ed93SDave May { 212521f74f9SMatthew G. Knepley PetscFunctionBegin; 21352c3ed93SDave May df->bs = blocksize; 21452c3ed93SDave May PetscFunctionReturn(0); 21552c3ed93SDave May } 21652c3ed93SDave May 2177fbf63aeSDave May PetscErrorCode DataFieldSetSize(DataField df,const PetscInt new_L) 21852849c42SDave May { 219521f74f9SMatthew G. Knepley PetscErrorCode ierr; 22052849c42SDave May 221521f74f9SMatthew G. Knepley PetscFunctionBegin; 222c8d15862SDave May if (new_L < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot set size of DataField to be < 0"); 2237fbf63aeSDave May if (new_L == df->L) PetscFunctionReturn(0); 22452849c42SDave May if (new_L > df->L) { 2254be7464cSMatthew G. Knepley ierr = PetscRealloc(df->atomic_size * (new_L), &df->data);CHKERRQ(ierr); 22652849c42SDave May /* init new contents */ 227521f74f9SMatthew G. Knepley ierr = PetscMemzero(( ((char*)df->data)+df->L*df->atomic_size), (new_L-df->L)*df->atomic_size);CHKERRQ(ierr); 2287fbf63aeSDave May } else { 22952849c42SDave May /* reallocate pointer list, add +1 in case new_L = 0 */ 2304be7464cSMatthew G. Knepley ierr = PetscRealloc(df->atomic_size * (new_L+1), &df->data);CHKERRQ(ierr); 23152849c42SDave May } 23252849c42SDave May df->L = new_L; 2337fbf63aeSDave May PetscFunctionReturn(0); 23452849c42SDave May } 23552849c42SDave May 2367fbf63aeSDave May PetscErrorCode DataFieldZeroBlock(DataField df,const PetscInt start,const PetscInt end) 23752849c42SDave May { 238521f74f9SMatthew G. Knepley PetscErrorCode ierr; 239521f74f9SMatthew G. Knepley 240521f74f9SMatthew G. Knepley PetscFunctionBegin; 2417fbf63aeSDave May if (start > end) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if start(%D) > end(%D)",start,end); 2427fbf63aeSDave May if (start < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if start(%D) < 0",start); 243a233d522SDave 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); 244521f74f9SMatthew G. Knepley ierr = PetscMemzero((((char*)df->data)+start*df->atomic_size), (end-start)*df->atomic_size);CHKERRQ(ierr); 2457fbf63aeSDave May PetscFunctionReturn(0); 24652849c42SDave May } 24752849c42SDave May 24852849c42SDave May /* 24952849c42SDave May A negative buffer value will simply be ignored and the old buffer value will be used. 25052849c42SDave May */ 2517fbf63aeSDave May PetscErrorCode DataBucketSetSizes(DataBucket db,const PetscInt L,const PetscInt buffer) 25252849c42SDave May { 2535c18a9d6SDave May PetscInt current_allocated,new_used,new_unused,new_buffer,new_allocated,f; 2540cb17e37SDave May PetscBool any_active_fields; 255dbe06d34SDave May PetscErrorCode ierr; 25652849c42SDave May 257521f74f9SMatthew G. Knepley PetscFunctionBegin; 2587fbf63aeSDave May if (db->finalised == PETSC_FALSE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"You must call DataBucketFinalize() before DataBucketSetSizes()"); 2590cb17e37SDave May ierr = DataBucketQueryForActiveFields(db,&any_active_fields);CHKERRQ(ierr); 2600cb17e37SDave 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"); 2610cb17e37SDave May 26252849c42SDave May current_allocated = db->allocated; 26352849c42SDave May new_used = L; 26452849c42SDave May new_unused = current_allocated - new_used; 26552849c42SDave May new_buffer = db->buffer; 26652849c42SDave May if (buffer >= 0) { /* update the buffer value */ 26752849c42SDave May new_buffer = buffer; 26852849c42SDave May } 26952849c42SDave May new_allocated = new_used + new_buffer; 27052849c42SDave May /* action */ 27152849c42SDave May if (new_allocated > current_allocated) { 27252849c42SDave May /* increase size to new_used + new_buffer */ 27352849c42SDave May for (f=0; f<db->nfields; f++) { 274dbe06d34SDave May ierr = DataFieldSetSize(db->field[f], new_allocated);CHKERRQ(ierr); 27552849c42SDave May } 27652849c42SDave May db->L = new_used; 27752849c42SDave May db->buffer = new_buffer; 27852849c42SDave May db->allocated = new_used + new_buffer; 279521f74f9SMatthew G. Knepley } else { 28052849c42SDave May if (new_unused > 2 * new_buffer) { 28152849c42SDave May /* shrink array to new_used + new_buffer */ 282521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 283dbe06d34SDave May ierr = DataFieldSetSize(db->field[f], new_allocated);CHKERRQ(ierr); 28452849c42SDave May } 28552849c42SDave May db->L = new_used; 28652849c42SDave May db->buffer = new_buffer; 28752849c42SDave May db->allocated = new_used + new_buffer; 288521f74f9SMatthew G. Knepley } else { 28952849c42SDave May db->L = new_used; 29052849c42SDave May db->buffer = new_buffer; 29152849c42SDave May } 29252849c42SDave May } 29352849c42SDave May /* zero all entries from db->L to db->allocated */ 294521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 29552849c42SDave May DataField field = db->field[f]; 296dbe06d34SDave May ierr = DataFieldZeroBlock(field, db->L,db->allocated);CHKERRQ(ierr); 29752849c42SDave May } 2987fbf63aeSDave May PetscFunctionReturn(0); 29952849c42SDave May } 30052849c42SDave May 3017fbf63aeSDave May PetscErrorCode DataBucketSetInitialSizes(DataBucket db,const PetscInt L,const PetscInt buffer) 30252849c42SDave May { 3035c18a9d6SDave May PetscInt f; 304dbe06d34SDave May PetscErrorCode ierr; 305dbe06d34SDave May 306521f74f9SMatthew G. Knepley PetscFunctionBegin; 307dbe06d34SDave May ierr = DataBucketSetSizes(db,L,buffer);CHKERRQ(ierr); 308521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 30952849c42SDave May DataField field = db->field[f]; 310dbe06d34SDave May ierr = DataFieldZeroBlock(field,0,db->allocated);CHKERRQ(ierr); 31152849c42SDave May } 3127fbf63aeSDave May PetscFunctionReturn(0); 31352849c42SDave May } 31452849c42SDave May 3157fbf63aeSDave May PetscErrorCode DataBucketGetSizes(DataBucket db,PetscInt *L,PetscInt *buffer,PetscInt *allocated) 31652849c42SDave May { 317521f74f9SMatthew G. Knepley PetscFunctionBegin; 31852849c42SDave May if (L) {*L = db->L;} 31952849c42SDave May if (buffer) {*buffer = db->buffer;} 32052849c42SDave May if (allocated) {*allocated = db->allocated;} 3217fbf63aeSDave May PetscFunctionReturn(0); 32252849c42SDave May } 32352849c42SDave May 3247fbf63aeSDave May PetscErrorCode DataBucketGetGlobalSizes(MPI_Comm comm,DataBucket db,PetscInt *L,PetscInt *buffer,PetscInt *allocated) 32552849c42SDave May { 3265c18a9d6SDave May PetscInt _L,_buffer,_allocated; 3275c18a9d6SDave May PetscInt ierr; 32852849c42SDave May 329521f74f9SMatthew G. Knepley PetscFunctionBegin; 3305c18a9d6SDave May _L = db->L; 3315c18a9d6SDave May _buffer = db->buffer; 3325c18a9d6SDave May _allocated = db->allocated; 33352849c42SDave May 33433564166SDave May if (L) { ierr = MPI_Allreduce(&_L,L,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); } 33533564166SDave May if (buffer) { ierr = MPI_Allreduce(&_buffer,buffer,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); } 33633564166SDave May if (allocated) { ierr = MPI_Allreduce(&_allocated,allocated,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); } 3377fbf63aeSDave May PetscFunctionReturn(0); 33852849c42SDave May } 33952849c42SDave May 3407fbf63aeSDave May PetscErrorCode DataBucketGetDataFields(DataBucket db,PetscInt *L,DataField *fields[]) 34152849c42SDave May { 342521f74f9SMatthew G. Knepley PetscFunctionBegin; 34352849c42SDave May if (L) {*L = db->nfields;} 34452849c42SDave May if (fields) {*fields = db->field;} 3457fbf63aeSDave May PetscFunctionReturn(0); 34652849c42SDave May } 34752849c42SDave May 3487fbf63aeSDave May PetscErrorCode DataFieldGetAccess(const DataField gfield) 34952849c42SDave May { 350521f74f9SMatthew G. Knepley PetscFunctionBegin; 3517fbf63aeSDave May if (gfield->active) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is already active. You must call DataFieldRestoreAccess()",gfield->name); 35252849c42SDave May gfield->active = PETSC_TRUE; 3537fbf63aeSDave May PetscFunctionReturn(0); 35452849c42SDave May } 35552849c42SDave May 3567fbf63aeSDave May PetscErrorCode DataFieldAccessPoint(const DataField gfield,const PetscInt pid,void **ctx_p) 35752849c42SDave May { 358521f74f9SMatthew G. Knepley PetscFunctionBegin; 35972505a4dSJed Brown *ctx_p = NULL; 36052849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD 36152849c42SDave May /* debug mode */ 36284bcda08SDave May /* check point is valid */ 3637fbf63aeSDave May if (pid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0"); 3647fbf63aeSDave May if (pid >= gfield->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",gfield->L); 36584bcda08SDave 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); 36652849c42SDave May #endif 36752849c42SDave May *ctx_p = __DATATFIELD_point_access(gfield->data,pid,gfield->atomic_size); 3687fbf63aeSDave May PetscFunctionReturn(0); 36952849c42SDave May } 37052849c42SDave May 3717fbf63aeSDave May PetscErrorCode DataFieldAccessPointOffset(const DataField gfield,const size_t offset,const PetscInt pid,void **ctx_p) 37252849c42SDave May { 373521f74f9SMatthew G. Knepley PetscFunctionBegin; 37452849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD 37552849c42SDave May /* debug mode */ 37684bcda08SDave May /* check point is valid */ 377521f74f9SMatthew G. Knepley /* if( offset < 0 ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be >= 0");*/ 378521f74f9SMatthew G. Knepley /* Note compiler realizes this can never happen with an unsigned PetscInt */ 3797fbf63aeSDave May if (offset >= gfield->atomic_size) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be < %zu",gfield->atomic_size); 38084bcda08SDave May /* check point is valid */ 3817fbf63aeSDave May if (pid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0"); 382a233d522SDave May if (pid >= gfield->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",gfield->L); 383a233d522SDave 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); 38452849c42SDave May #endif 38552849c42SDave May *ctx_p = __DATATFIELD_point_access_offset(gfield->data,pid,gfield->atomic_size,offset); 3867fbf63aeSDave May PetscFunctionReturn(0); 38752849c42SDave May } 38852849c42SDave May 3897fbf63aeSDave May PetscErrorCode DataFieldRestoreAccess(DataField gfield) 39052849c42SDave May { 391521f74f9SMatthew G. Knepley PetscFunctionBegin; 3927fbf63aeSDave May if (gfield->active == PETSC_FALSE) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is not active. You must call DataFieldGetAccess()", gfield->name); 39352849c42SDave May gfield->active = PETSC_FALSE; 3947fbf63aeSDave May PetscFunctionReturn(0); 39552849c42SDave May } 39652849c42SDave May 3977fbf63aeSDave May PetscErrorCode DataFieldVerifyAccess(const DataField gfield,const size_t size) 39852849c42SDave May { 399521f74f9SMatthew G. Knepley PetscFunctionBegin; 40052849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD 401521f74f9SMatthew 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 ); 40252849c42SDave May #endif 4037fbf63aeSDave May PetscFunctionReturn(0); 40452849c42SDave May } 40552849c42SDave May 4067fbf63aeSDave May PetscErrorCode DataFieldGetAtomicSize(const DataField gfield,size_t *size) 40752849c42SDave May { 408521f74f9SMatthew G. Knepley PetscFunctionBegin; 40952849c42SDave May if (size) {*size = gfield->atomic_size;} 4107fbf63aeSDave May PetscFunctionReturn(0); 41152849c42SDave May } 41252849c42SDave May 4137fbf63aeSDave May PetscErrorCode DataFieldGetEntries(const DataField gfield,void **data) 41452849c42SDave May { 415521f74f9SMatthew G. Knepley PetscFunctionBegin; 416521f74f9SMatthew G. Knepley if (data) {*data = gfield->data;} 4177fbf63aeSDave May PetscFunctionReturn(0); 41852849c42SDave May } 41952849c42SDave May 4207fbf63aeSDave May PetscErrorCode DataFieldRestoreEntries(const DataField gfield,void **data) 42152849c42SDave May { 422521f74f9SMatthew G. Knepley PetscFunctionBegin; 423521f74f9SMatthew G. Knepley if (data) {*data = NULL;} 4247fbf63aeSDave May PetscFunctionReturn(0); 42552849c42SDave May } 42652849c42SDave May 42752849c42SDave May /* y = x */ 4287fbf63aeSDave May PetscErrorCode DataBucketCopyPoint(const DataBucket xb,const PetscInt pid_x, 4295c18a9d6SDave May const DataBucket yb,const PetscInt pid_y) 43052849c42SDave May { 4315c18a9d6SDave May PetscInt f; 432dbe06d34SDave May PetscErrorCode ierr; 433dbe06d34SDave May 434521f74f9SMatthew G. Knepley PetscFunctionBegin; 435521f74f9SMatthew G. Knepley for (f = 0; f < xb->nfields; ++f) { 43652849c42SDave May void *dest; 43752849c42SDave May void *src; 43852849c42SDave May 439dbe06d34SDave May ierr = DataFieldGetAccess(xb->field[f]);CHKERRQ(ierr); 4406845f8f5SDave May if (xb != yb) { ierr = DataFieldGetAccess( yb->field[f]);CHKERRQ(ierr); } 441dbe06d34SDave May ierr = DataFieldAccessPoint(xb->field[f],pid_x, &src);CHKERRQ(ierr); 442dbe06d34SDave May ierr = DataFieldAccessPoint(yb->field[f],pid_y, &dest);CHKERRQ(ierr); 443521f74f9SMatthew G. Knepley ierr = PetscMemcpy(dest, src, xb->field[f]->atomic_size);CHKERRQ(ierr); 444dbe06d34SDave May ierr = DataFieldRestoreAccess(xb->field[f]);CHKERRQ(ierr); 4456845f8f5SDave May if (xb != yb) {ierr = DataFieldRestoreAccess(yb->field[f]);CHKERRQ(ierr);} 44652849c42SDave May } 4477fbf63aeSDave May PetscFunctionReturn(0); 44852849c42SDave May } 44952849c42SDave May 4507fbf63aeSDave May PetscErrorCode DataBucketCreateFromSubset(DataBucket DBIn,const PetscInt N,const PetscInt list[],DataBucket *DB) 45152849c42SDave May { 4525c18a9d6SDave May PetscInt nfields; 45352849c42SDave May DataField *fields; 4545c18a9d6SDave May PetscInt f,L,buffer,allocated,p; 455dbe06d34SDave May PetscErrorCode ierr; 45652849c42SDave May 457521f74f9SMatthew G. Knepley PetscFunctionBegin; 458521f74f9SMatthew G. Knepley ierr = DataBucketCreate(DB);CHKERRQ(ierr); 45952849c42SDave May /* copy contents of DBIn */ 460dbe06d34SDave May ierr = DataBucketGetDataFields(DBIn,&nfields,&fields);CHKERRQ(ierr); 461dbe06d34SDave May ierr = DataBucketGetSizes(DBIn,&L,&buffer,&allocated);CHKERRQ(ierr); 462521f74f9SMatthew G. Knepley for (f = 0; f < nfields; ++f) { 463dbe06d34SDave May ierr = DataBucketRegisterField(*DB,"DataBucketCreateFromSubset",fields[f]->name,fields[f]->atomic_size,NULL);CHKERRQ(ierr); 46452849c42SDave May } 465dbe06d34SDave May ierr = DataBucketFinalize(*DB);CHKERRQ(ierr); 466dbe06d34SDave May ierr = DataBucketSetSizes(*DB,L,buffer);CHKERRQ(ierr); 46752849c42SDave May /* now copy the desired guys from DBIn => DB */ 468521f74f9SMatthew G. Knepley for (p = 0; p < N; ++p) { 469dbe06d34SDave May ierr = DataBucketCopyPoint(DBIn,list[p], *DB,p);CHKERRQ(ierr); 47052849c42SDave May } 4717fbf63aeSDave May PetscFunctionReturn(0); 47252849c42SDave May } 47352849c42SDave May 474521f74f9SMatthew G. Knepley /* insert into an exisitng location */ 4757fbf63aeSDave May PetscErrorCode DataFieldInsertPoint(const DataField field,const PetscInt index,const void *ctx) 47652849c42SDave May { 477521f74f9SMatthew G. Knepley PetscErrorCode ierr; 47852849c42SDave May 479521f74f9SMatthew G. Knepley PetscFunctionBegin; 48052849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD 48184bcda08SDave May /* check point is valid */ 482a233d522SDave May if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0"); 483a233d522SDave May if (index >= field->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",field->L); 48452849c42SDave May #endif 485521f74f9SMatthew G. Knepley ierr = PetscMemcpy(__DATATFIELD_point_access(field->data,index,field->atomic_size), ctx, field->atomic_size);CHKERRQ(ierr); 4867fbf63aeSDave May PetscFunctionReturn(0); 48752849c42SDave May } 48852849c42SDave May 489521f74f9SMatthew G. Knepley /* remove data at index - replace with last point */ 4907fbf63aeSDave May PetscErrorCode DataBucketRemovePointAtIndex(const DataBucket db,const PetscInt index) 49152849c42SDave May { 4925c18a9d6SDave May PetscInt f; 4930cb17e37SDave May PetscBool any_active_fields; 494dbe06d34SDave May PetscErrorCode ierr; 49552849c42SDave May 496521f74f9SMatthew G. Knepley PetscFunctionBegin; 49752849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD 49884bcda08SDave May /* check point is valid */ 499a233d522SDave May if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0"); 500a233d522SDave May if (index >= db->allocated) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",db->L+db->buffer); 50152849c42SDave May #endif 5020cb17e37SDave May ierr = DataBucketQueryForActiveFields(db,&any_active_fields);CHKERRQ(ierr); 5030cb17e37SDave 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"); 504a233d522SDave May if (index >= db->L) { /* this point is not in the list - no need to error, but I will anyway */ 505a233d522SDave 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); 50652849c42SDave May } 507a233d522SDave May if (index != db->L-1) { /* not last point in list */ 508521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 50952849c42SDave May DataField field = db->field[f]; 51052849c42SDave May 51152849c42SDave May /* copy then remove */ 512dbe06d34SDave May ierr = DataFieldCopyPoint(db->L-1, field, index, field);CHKERRQ(ierr); 513521f74f9SMatthew G. Knepley /* DataFieldZeroPoint(field,index); */ 51452849c42SDave May } 51552849c42SDave May } 51652849c42SDave May /* decrement size */ 51752849c42SDave May /* this will zero out an crap at the end of the list */ 518dbe06d34SDave May ierr = DataBucketRemovePoint(db);CHKERRQ(ierr); 5197fbf63aeSDave May PetscFunctionReturn(0); 52052849c42SDave May } 52152849c42SDave May 52252849c42SDave May /* copy x into y */ 5237fbf63aeSDave May PetscErrorCode DataFieldCopyPoint(const PetscInt pid_x,const DataField field_x, 5245c18a9d6SDave May const PetscInt pid_y,const DataField field_y ) 52552849c42SDave May { 526521f74f9SMatthew G. Knepley PetscErrorCode ierr; 52752849c42SDave May 528521f74f9SMatthew G. Knepley PetscFunctionBegin; 52952849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD 53084bcda08SDave May /* check point is valid */ 531a233d522SDave May if (pid_x < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"(IN) index must be >= 0"); 532a233d522SDave May if (pid_x >= field_x->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"(IN) index must be < %D",field_x->L); 533a233d522SDave May if (pid_y < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"(OUT) index must be >= 0"); 534a233d522SDave May if (pid_y >= field_y->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"(OUT) index must be < %D",field_y->L); 535521f74f9SMatthew G. Knepley if( field_y->atomic_size != field_x->atomic_size ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"atomic size must match"); 53652849c42SDave May #endif 537*298827fbSBarry Smith ierr = PetscMemcpy(__DATATFIELD_point_access(field_y->data,pid_y,field_y->atomic_size),__DATATFIELD_point_access(field_x->data,pid_x,field_x->atomic_size),field_y->atomic_size);CHKERRQ(ierr); 5387fbf63aeSDave May PetscFunctionReturn(0); 53952849c42SDave May } 54052849c42SDave May 54152849c42SDave May 542521f74f9SMatthew G. Knepley /* zero only the datafield at this point */ 5437fbf63aeSDave May PetscErrorCode DataFieldZeroPoint(const DataField field,const PetscInt index) 54452849c42SDave May { 545521f74f9SMatthew G. Knepley PetscErrorCode ierr; 546521f74f9SMatthew G. Knepley 547521f74f9SMatthew G. Knepley PetscFunctionBegin; 54852849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD 54984bcda08SDave May /* check point is valid */ 550a233d522SDave May if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0"); 551a233d522SDave May if (index >= field->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",field->L); 55252849c42SDave May #endif 553521f74f9SMatthew G. Knepley ierr = PetscMemzero(__DATATFIELD_point_access(field->data,index,field->atomic_size), field->atomic_size);CHKERRQ(ierr); 5547fbf63aeSDave May PetscFunctionReturn(0); 55552849c42SDave May } 55652849c42SDave May 557521f74f9SMatthew G. Knepley /* zero ALL data for this point */ 5587fbf63aeSDave May PetscErrorCode DataBucketZeroPoint(const DataBucket db,const PetscInt index) 55952849c42SDave May { 5605c18a9d6SDave May PetscInt f; 561dbe06d34SDave May PetscErrorCode ierr; 56252849c42SDave May 563521f74f9SMatthew G. Knepley PetscFunctionBegin; 56484bcda08SDave May /* check point is valid */ 565a233d522SDave May if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0"); 566a233d522SDave May if (index >= db->allocated) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",db->allocated); 567521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 56852849c42SDave May DataField field = db->field[f]; 569dbe06d34SDave May ierr = DataFieldZeroPoint(field,index);CHKERRQ(ierr); 57052849c42SDave May } 5717fbf63aeSDave May PetscFunctionReturn(0); 57252849c42SDave May } 57352849c42SDave May 57452849c42SDave May /* increment */ 5757fbf63aeSDave May PetscErrorCode DataBucketAddPoint(DataBucket db) 57652849c42SDave May { 577dbe06d34SDave May PetscErrorCode ierr; 578dbe06d34SDave May 579521f74f9SMatthew G. Knepley PetscFunctionBegin; 58065141ba8SDave May ierr = DataBucketSetSizes(db,db->L+1,DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr); 5817fbf63aeSDave May PetscFunctionReturn(0); 58252849c42SDave May } 58352849c42SDave May 5847fbf63aeSDave May /* decrement */ 5857fbf63aeSDave May PetscErrorCode DataBucketRemovePoint(DataBucket db) 5867fbf63aeSDave May { 587dbe06d34SDave May PetscErrorCode ierr; 588dbe06d34SDave May 589521f74f9SMatthew G. Knepley PetscFunctionBegin; 59065141ba8SDave May ierr = DataBucketSetSizes(db,db->L-1,DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr); 5917fbf63aeSDave May PetscFunctionReturn(0); 5927fbf63aeSDave May } 5937fbf63aeSDave May 594ab2be9a3SDave May PetscErrorCode DataBucketView_stdout(MPI_Comm comm,DataBucket db) 595ab2be9a3SDave May { 596ab2be9a3SDave May PetscInt f; 597ab2be9a3SDave May double memory_usage_total,memory_usage_total_local = 0.0; 598ab2be9a3SDave May PetscErrorCode ierr; 599ab2be9a3SDave May 600ab2be9a3SDave May PetscFunctionBegin; 601ab2be9a3SDave May ierr = PetscPrintf(comm,"DataBucketView: \n");CHKERRQ(ierr); 602ab2be9a3SDave May ierr = PetscPrintf(comm," L = %D \n", db->L);CHKERRQ(ierr); 603ab2be9a3SDave May ierr = PetscPrintf(comm," buffer = %D \n", db->buffer);CHKERRQ(ierr); 604ab2be9a3SDave May ierr = PetscPrintf(comm," allocated = %D \n", db->allocated);CHKERRQ(ierr); 605ab2be9a3SDave May ierr = PetscPrintf(comm," nfields registered = %D \n", db->nfields);CHKERRQ(ierr); 606ab2be9a3SDave May 607ab2be9a3SDave May for (f = 0; f < db->nfields; ++f) { 608ab2be9a3SDave May double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6; 609ab2be9a3SDave May memory_usage_total_local += memory_usage_f; 610ab2be9a3SDave May } 611ab2be9a3SDave May ierr = MPI_Allreduce(&memory_usage_total_local,&memory_usage_total,1,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr); 612ab2be9a3SDave May 613ab2be9a3SDave May for (f = 0; f < db->nfields; ++f) { 614ab2be9a3SDave May double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6; 615ab2be9a3SDave May ierr = PetscPrintf(comm," [%3D] %15s : Mem. usage = %1.2e (MB) [rank0]\n", f, db->field[f]->name, memory_usage_f );CHKERRQ(ierr); 616ab2be9a3SDave May ierr = PetscPrintf(comm," blocksize = %D \n", db->field[f]->bs);CHKERRQ(ierr); 617ab2be9a3SDave May if (db->field[f]->bs != 1) { 618ab2be9a3SDave May ierr = PetscPrintf(comm," atomic size = %zu [full block, bs=%D]\n", db->field[f]->atomic_size,db->field[f]->bs);CHKERRQ(ierr); 619ab2be9a3SDave May ierr = PetscPrintf(comm," atomic size/item = %zu \n", db->field[f]->atomic_size/db->field[f]->bs);CHKERRQ(ierr); 620ab2be9a3SDave May } else { 621ab2be9a3SDave May ierr = PetscPrintf(comm," atomic size = %zu \n", db->field[f]->atomic_size);CHKERRQ(ierr); 622ab2be9a3SDave May } 623ab2be9a3SDave May } 624ab2be9a3SDave May ierr = PetscPrintf(comm," Total mem. usage = %1.2e (MB) (collective)\n", memory_usage_total);CHKERRQ(ierr); 625ab2be9a3SDave May PetscFunctionReturn(0); 626ab2be9a3SDave May } 627ab2be9a3SDave May 628ab2be9a3SDave May PetscErrorCode DataBucketView_SEQ(MPI_Comm comm,DataBucket db,const char filename[],DataBucketViewType type) 62952849c42SDave May { 630dbe06d34SDave May PetscErrorCode ierr; 631dbe06d34SDave May 632521f74f9SMatthew G. Knepley PetscFunctionBegin; 63352849c42SDave May switch (type) { 63452849c42SDave May case DATABUCKET_VIEW_STDOUT: 635ab2be9a3SDave May ierr = DataBucketView_stdout(PETSC_COMM_SELF,db);CHKERRQ(ierr); 63652849c42SDave May break; 63752849c42SDave May case DATABUCKET_VIEW_ASCII: 638071ce369SDave May SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for ascii output"); 63952849c42SDave May break; 64052849c42SDave May case DATABUCKET_VIEW_BINARY: 641071ce369SDave May SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for binary output"); 64252849c42SDave May break; 64352849c42SDave May case DATABUCKET_VIEW_HDF5: 644071ce369SDave May SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for HDF5 output"); 64552849c42SDave May break; 646521f74f9SMatthew G. Knepley default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer method requested"); 64752849c42SDave May } 6487fbf63aeSDave May PetscFunctionReturn(0); 64952849c42SDave May } 65052849c42SDave May 6517fbf63aeSDave May PetscErrorCode DataBucketView_MPI(MPI_Comm comm,DataBucket db,const char filename[],DataBucketViewType type) 65252849c42SDave May { 653997fa542SDave May PetscErrorCode ierr; 654997fa542SDave May 655521f74f9SMatthew G. Knepley PetscFunctionBegin; 65652849c42SDave May switch (type) { 65752849c42SDave May case DATABUCKET_VIEW_STDOUT: 658ab2be9a3SDave May ierr = DataBucketView_stdout(comm,db);CHKERRQ(ierr); 65952849c42SDave May break; 66052849c42SDave May case DATABUCKET_VIEW_ASCII: 661071ce369SDave May SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for ascii output"); 66252849c42SDave May break; 66352849c42SDave May case DATABUCKET_VIEW_BINARY: 664071ce369SDave May SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for binary output"); 66552849c42SDave May break; 66652849c42SDave May case DATABUCKET_VIEW_HDF5: 667071ce369SDave May SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for HDF5 output"); 66852849c42SDave May break; 669521f74f9SMatthew G. Knepley default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer method requested"); 67052849c42SDave May } 6717fbf63aeSDave May PetscFunctionReturn(0); 67252849c42SDave May } 67352849c42SDave May 6747fbf63aeSDave May PetscErrorCode DataBucketView(MPI_Comm comm,DataBucket db,const char filename[],DataBucketViewType type) 67552849c42SDave May { 6765c18a9d6SDave May PetscMPIInt nproc; 677997fa542SDave May PetscErrorCode ierr; 67852849c42SDave May 679521f74f9SMatthew G. Knepley PetscFunctionBegin; 680997fa542SDave May ierr = MPI_Comm_size(comm,&nproc);CHKERRQ(ierr); 68152849c42SDave May if (nproc == 1) { 682ab2be9a3SDave May ierr = DataBucketView_SEQ(comm,db,filename,type);CHKERRQ(ierr); 68352849c42SDave May } else { 684dbe06d34SDave May ierr = DataBucketView_MPI(comm,db,filename,type);CHKERRQ(ierr); 68552849c42SDave May } 6867fbf63aeSDave May PetscFunctionReturn(0); 68752849c42SDave May } 68852849c42SDave May 6897fbf63aeSDave May PetscErrorCode DataBucketDuplicateFields(DataBucket dbA,DataBucket *dbB) 69052849c42SDave May { 69152849c42SDave May DataBucket db2; 6925c18a9d6SDave May PetscInt f; 693dbe06d34SDave May PetscErrorCode ierr; 69452849c42SDave May 695521f74f9SMatthew G. Knepley PetscFunctionBegin; 696dbe06d34SDave May ierr = DataBucketCreate(&db2);CHKERRQ(ierr); 69752849c42SDave May /* copy contents from dbA into db2 */ 698521f74f9SMatthew G. Knepley for (f = 0; f < dbA->nfields; ++f) { 69952849c42SDave May DataField field; 70052849c42SDave May size_t atomic_size; 70152849c42SDave May char *name; 70252849c42SDave May 70352849c42SDave May field = dbA->field[f]; 70452849c42SDave May atomic_size = field->atomic_size; 70552849c42SDave May name = field->name; 706dbe06d34SDave May ierr = DataBucketRegisterField(db2,"DataBucketDuplicateFields",name,atomic_size,NULL);CHKERRQ(ierr); 70752849c42SDave May } 708dbe06d34SDave May ierr = DataBucketFinalize(db2);CHKERRQ(ierr); 709dbe06d34SDave May ierr = DataBucketSetInitialSizes(db2,0,1000);CHKERRQ(ierr); 71052849c42SDave May *dbB = db2; 7117fbf63aeSDave May PetscFunctionReturn(0); 71252849c42SDave May } 71352849c42SDave May 71452849c42SDave May /* 71552849c42SDave May Insert points from db2 into db1 71652849c42SDave May db1 <<== db2 71752849c42SDave May */ 7187fbf63aeSDave May PetscErrorCode DataBucketInsertValues(DataBucket db1,DataBucket db2) 71952849c42SDave May { 7205c18a9d6SDave May PetscInt n_mp_points1,n_mp_points2; 7215c18a9d6SDave May PetscInt n_mp_points1_new,p; 722dbe06d34SDave May PetscErrorCode ierr; 72352849c42SDave May 724521f74f9SMatthew G. Knepley PetscFunctionBegin; 725dbe06d34SDave May ierr = DataBucketGetSizes(db1,&n_mp_points1,0,0);CHKERRQ(ierr); 726dbe06d34SDave May ierr = DataBucketGetSizes(db2,&n_mp_points2,0,0);CHKERRQ(ierr); 72752849c42SDave May n_mp_points1_new = n_mp_points1 + n_mp_points2; 72865141ba8SDave May ierr = DataBucketSetSizes(db1,n_mp_points1_new,DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr); 729521f74f9SMatthew G. Knepley for (p = 0; p < n_mp_points2; ++p) { 730521f74f9SMatthew G. Knepley /* db1 <<== db2 */ 731dbe06d34SDave May ierr = DataBucketCopyPoint(db2,p, db1,(n_mp_points1 + p));CHKERRQ(ierr); 73252849c42SDave May } 7337fbf63aeSDave May PetscFunctionReturn(0); 73452849c42SDave May } 73552849c42SDave May 73652849c42SDave May /* helpers for parallel send/recv */ 7377fbf63aeSDave May PetscErrorCode DataBucketCreatePackedArray(DataBucket db,size_t *bytes,void **buf) 73852849c42SDave May { 7395c18a9d6SDave May PetscInt f; 74052849c42SDave May size_t sizeof_marker_contents; 74152849c42SDave May void *buffer; 742521f74f9SMatthew G. Knepley PetscErrorCode ierr; 74352849c42SDave May 744521f74f9SMatthew G. Knepley PetscFunctionBegin; 74552849c42SDave May sizeof_marker_contents = 0; 746521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 74752849c42SDave May DataField df = db->field[f]; 74852849c42SDave May sizeof_marker_contents += df->atomic_size; 74952849c42SDave May } 750521f74f9SMatthew G. Knepley ierr = PetscMalloc(sizeof_marker_contents, &buffer);CHKERRQ(ierr); 751521f74f9SMatthew G. Knepley ierr = PetscMemzero(buffer, sizeof_marker_contents);CHKERRQ(ierr); 75252849c42SDave May if (bytes) {*bytes = sizeof_marker_contents;} 75352849c42SDave May if (buf) {*buf = buffer;} 7547fbf63aeSDave May PetscFunctionReturn(0); 75552849c42SDave May } 75652849c42SDave May 7577fbf63aeSDave May PetscErrorCode DataBucketDestroyPackedArray(DataBucket db,void **buf) 75852849c42SDave May { 759521f74f9SMatthew G. Knepley PetscErrorCode ierr; 760521f74f9SMatthew G. Knepley 761521f74f9SMatthew G. Knepley PetscFunctionBegin; 76252849c42SDave May if (buf) { 763521f74f9SMatthew G. Knepley ierr = PetscFree(*buf);CHKERRQ(ierr); 76452849c42SDave May *buf = NULL; 76552849c42SDave May } 7667fbf63aeSDave May PetscFunctionReturn(0); 76752849c42SDave May } 76852849c42SDave May 7697fbf63aeSDave May PetscErrorCode DataBucketFillPackedArray(DataBucket db,const PetscInt index,void *buf) 77052849c42SDave May { 7715c18a9d6SDave May PetscInt f; 77252849c42SDave May void *data, *data_p; 77352849c42SDave May size_t asize, offset; 774521f74f9SMatthew G. Knepley PetscErrorCode ierr; 77552849c42SDave May 776521f74f9SMatthew G. Knepley PetscFunctionBegin; 77752849c42SDave May offset = 0; 778521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 77952849c42SDave May DataField df = db->field[f]; 78052849c42SDave May 78152849c42SDave May asize = df->atomic_size; 78252849c42SDave May data = (void*)( df->data ); 78352849c42SDave May data_p = (void*)( (char*)data + index*asize ); 784521f74f9SMatthew G. Knepley ierr = PetscMemcpy((void*)((char*)buf + offset), data_p, asize);CHKERRQ(ierr); 78552849c42SDave May offset = offset + asize; 78652849c42SDave May } 7877fbf63aeSDave May PetscFunctionReturn(0); 78852849c42SDave May } 78952849c42SDave May 7907fbf63aeSDave May PetscErrorCode DataBucketInsertPackedArray(DataBucket db,const PetscInt idx,void *data) 79152849c42SDave May { 7925c18a9d6SDave May PetscInt f; 79352849c42SDave May void *data_p; 79452849c42SDave May size_t offset; 795dbe06d34SDave May PetscErrorCode ierr; 79652849c42SDave May 797521f74f9SMatthew G. Knepley PetscFunctionBegin; 79852849c42SDave May offset = 0; 799521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 80052849c42SDave May DataField df = db->field[f]; 80152849c42SDave May 80252849c42SDave May data_p = (void*)( (char*)data + offset ); 803dbe06d34SDave May ierr = DataFieldInsertPoint(df, idx, (void*)data_p);CHKERRQ(ierr); 80452849c42SDave May offset = offset + df->atomic_size; 80552849c42SDave May } 8067fbf63aeSDave May PetscFunctionReturn(0); 80752849c42SDave May } 808