152849c42SDave May 252849c42SDave May #include "data_bucket.h" 352849c42SDave May 452849c42SDave May /* string helpers */ 52eac95f8SDave May PetscErrorCode StringInList(const char name[],const PetscInt N,const DataField gfield[],PetscBool *val) 652849c42SDave May { 75c18a9d6SDave May PetscInt i; 8521f74f9SMatthew G. Knepley PetscErrorCode ierr; 952849c42SDave May 10521f74f9SMatthew G. Knepley PetscFunctionBegin; 1152849c42SDave May *val = PETSC_FALSE; 12521f74f9SMatthew G. Knepley for (i = 0; i < N; ++i) { 13521f74f9SMatthew G. Knepley PetscBool flg; 14521f74f9SMatthew G. Knepley ierr = PetscStrcmp(name, gfield[i]->name, &flg);CHKERRQ(ierr); 15521f74f9SMatthew G. Knepley if (flg) { 1652849c42SDave May *val = PETSC_TRUE; 172eac95f8SDave May PetscFunctionReturn(0); 1852849c42SDave May } 1952849c42SDave May } 202eac95f8SDave May PetscFunctionReturn(0); 2152849c42SDave May } 2252849c42SDave May 232eac95f8SDave May PetscErrorCode StringFindInList(const char name[],const PetscInt N,const DataField gfield[],PetscInt *index) 2452849c42SDave May { 255c18a9d6SDave May PetscInt i; 26521f74f9SMatthew G. Knepley PetscErrorCode ierr; 2752849c42SDave May 28521f74f9SMatthew G. Knepley PetscFunctionBegin; 2952849c42SDave May *index = -1; 30521f74f9SMatthew G. Knepley for (i = 0; i < N; ++i) { 31521f74f9SMatthew G. Knepley PetscBool flg; 32521f74f9SMatthew G. Knepley ierr = PetscStrcmp(name, gfield[i]->name, &flg);CHKERRQ(ierr); 33521f74f9SMatthew G. Knepley if (flg) { 3452849c42SDave May *index = i; 352eac95f8SDave May PetscFunctionReturn(0); 3652849c42SDave May } 3752849c42SDave May } 382eac95f8SDave May PetscFunctionReturn(0); 3952849c42SDave May } 4052849c42SDave May 412eac95f8SDave May PetscErrorCode DataFieldCreate(const char registeration_function[],const char name[],const size_t size,const PetscInt L,DataField *DF) 4252849c42SDave May { 4352849c42SDave May DataField df; 44521f74f9SMatthew G. Knepley PetscErrorCode ierr; 4552849c42SDave May 46521f74f9SMatthew G. Knepley PetscFunctionBegin; 47521f74f9SMatthew G. Knepley ierr = PetscMalloc(sizeof(struct _p_DataField), &df);CHKERRQ(ierr); 48521f74f9SMatthew G. Knepley ierr = PetscMemzero(df, sizeof(struct _p_DataField));CHKERRQ(ierr); 49521f74f9SMatthew G. Knepley ierr = PetscStrallocpy(registeration_function, &df->registeration_function);CHKERRQ(ierr); 50521f74f9SMatthew G. Knepley ierr = PetscStrallocpy(name, &df->name);CHKERRQ(ierr); 5152849c42SDave May df->atomic_size = size; 5252849c42SDave May df->L = L; 5352c3ed93SDave May df->bs = 1; 54521f74f9SMatthew G. Knepley /* allocate something so we don't have to reallocate */ 55521f74f9SMatthew G. Knepley ierr = PetscMalloc(size * L, &df->data);CHKERRQ(ierr); 56521f74f9SMatthew G. Knepley ierr = PetscMemzero(df->data, size * L);CHKERRQ(ierr); 5752849c42SDave May *DF = df; 582eac95f8SDave May PetscFunctionReturn(0); 5952849c42SDave May } 6052849c42SDave May 612eac95f8SDave May PetscErrorCode DataFieldDestroy(DataField *DF) 6252849c42SDave May { 6352849c42SDave May DataField df = *DF; 64521f74f9SMatthew G. Knepley PetscErrorCode ierr; 6552849c42SDave May 66521f74f9SMatthew G. Knepley PetscFunctionBegin; 67521f74f9SMatthew G. Knepley ierr = PetscFree(df->registeration_function);CHKERRQ(ierr); 68521f74f9SMatthew G. Knepley ierr = PetscFree(df->name);CHKERRQ(ierr); 69521f74f9SMatthew G. Knepley ierr = PetscFree(df->data);CHKERRQ(ierr); 70521f74f9SMatthew G. Knepley ierr = PetscFree(df);CHKERRQ(ierr); 7152849c42SDave May *DF = NULL; 722eac95f8SDave May PetscFunctionReturn(0); 7352849c42SDave May } 7452849c42SDave May 7552849c42SDave May /* data bucket */ 762eac95f8SDave May PetscErrorCode DataBucketCreate(DataBucket *DB) 7752849c42SDave May { 7852849c42SDave May DataBucket db; 79521f74f9SMatthew G. Knepley PetscErrorCode ierr; 8052849c42SDave May 81521f74f9SMatthew G. Knepley PetscFunctionBegin; 82521f74f9SMatthew G. Knepley ierr = PetscMalloc(sizeof(struct _p_DataBucket), &db);CHKERRQ(ierr); 83521f74f9SMatthew G. Knepley ierr = PetscMemzero(db, sizeof(struct _p_DataBucket));CHKERRQ(ierr); 8452849c42SDave May 8552849c42SDave May db->finalised = PETSC_FALSE; 8652849c42SDave May /* create empty spaces for fields */ 873454631fSDave May db->L = -1; 8852849c42SDave May db->buffer = 1; 8952849c42SDave May db->allocated = 1; 9052849c42SDave May db->nfields = 0; 91521f74f9SMatthew G. Knepley ierr = PetscMalloc1(1, &db->field);CHKERRQ(ierr); 9252849c42SDave May *DB = db; 932eac95f8SDave May PetscFunctionReturn(0); 9452849c42SDave May } 9552849c42SDave May 962eac95f8SDave May PetscErrorCode DataBucketDestroy(DataBucket *DB) 9752849c42SDave May { 9852849c42SDave May DataBucket db = *DB; 995c18a9d6SDave May PetscInt f; 100dbe06d34SDave May PetscErrorCode ierr; 10152849c42SDave May 102521f74f9SMatthew G. Knepley PetscFunctionBegin; 10352849c42SDave May /* release fields */ 104521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 105dbe06d34SDave May ierr = DataFieldDestroy(&db->field[f]);CHKERRQ(ierr); 10652849c42SDave May } 10752849c42SDave May /* this will catch the initially allocated objects in the event that no fields are registered */ 10852849c42SDave May if (db->field != NULL) { 109521f74f9SMatthew G. Knepley ierr = PetscFree(db->field);CHKERRQ(ierr); 11052849c42SDave May } 111521f74f9SMatthew G. Knepley ierr = PetscFree(db);CHKERRQ(ierr); 11252849c42SDave May *DB = NULL; 1132eac95f8SDave May PetscFunctionReturn(0); 11452849c42SDave May } 11552849c42SDave May 1160cb17e37SDave May PetscErrorCode DataBucketQueryForActiveFields(DataBucket db,PetscBool *any_active_fields) 1170cb17e37SDave May { 1180cb17e37SDave May PetscInt f; 119521f74f9SMatthew G. Knepley 120521f74f9SMatthew G. Knepley PetscFunctionBegin; 1210cb17e37SDave May *any_active_fields = PETSC_FALSE; 122521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 1230cb17e37SDave May if (db->field[f]->active) { 1240cb17e37SDave May *any_active_fields = PETSC_TRUE; 1250cb17e37SDave May PetscFunctionReturn(0); 1260cb17e37SDave May } 1270cb17e37SDave May } 1280cb17e37SDave May PetscFunctionReturn(0); 1290cb17e37SDave May } 1300cb17e37SDave May 1312eac95f8SDave May PetscErrorCode DataBucketRegisterField( 13252849c42SDave May DataBucket db, 13352849c42SDave May const char registeration_function[], 13452849c42SDave May const char field_name[], 13552849c42SDave May size_t atomic_size, DataField *_gfield) 13652849c42SDave May { 13752849c42SDave May PetscBool val; 1384be7464cSMatthew G. Knepley DataField fp; 139dbe06d34SDave May PetscErrorCode ierr; 14052849c42SDave May 141521f74f9SMatthew G. Knepley PetscFunctionBegin; 14252849c42SDave May /* check we haven't finalised the registration of fields */ 14352849c42SDave May /* 14452849c42SDave May if(db->finalised==PETSC_TRUE) { 14552849c42SDave May printf("ERROR: DataBucketFinalize() has been called. Cannot register more fields\n"); 14652849c42SDave May ERROR(); 14752849c42SDave May } 14852849c42SDave May */ 14952849c42SDave May /* check for repeated name */ 1502635f519SDave May ierr = StringInList(field_name, db->nfields, (const DataField*) db->field, &val);CHKERRQ(ierr); 1512eac95f8SDave May if (val == PETSC_TRUE) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field %s already exists. Cannot add same field twice",field_name); 15252849c42SDave May /* create new space for data */ 1534be7464cSMatthew G. Knepley ierr = PetscRealloc(sizeof(DataField)*(db->nfields+1), &db->field);CHKERRQ(ierr); 15452849c42SDave May /* add field */ 155dbe06d34SDave May ierr = DataFieldCreate(registeration_function, field_name, atomic_size, db->allocated, &fp);CHKERRQ(ierr); 15652849c42SDave May db->field[db->nfields] = fp; 15752849c42SDave May db->nfields++; 15852849c42SDave May if (_gfield != NULL) { 15952849c42SDave May *_gfield = fp; 16052849c42SDave May } 1612eac95f8SDave May PetscFunctionReturn(0); 16252849c42SDave May } 16352849c42SDave May 16452849c42SDave May /* 16552849c42SDave May #define DataBucketRegisterField(db,name,size,k) {\ 16652849c42SDave May char *location;\ 16752849c42SDave May asprintf(&location,"Registered by %s() at line %d within file %s", __FUNCTION__, __LINE__, __FILE__);\ 16852849c42SDave May _DataBucketRegisterField( (db), location, (name), (size), (k) );\ 169521f74f9SMatthew G. Knepley ierr = PetscFree(location);\ 17052849c42SDave May } 17152849c42SDave May */ 17252849c42SDave May 1732eac95f8SDave May PetscErrorCode DataBucketGetDataFieldByName(DataBucket db,const char name[],DataField *gfield) 17452849c42SDave May { 1755c18a9d6SDave May PetscInt idx; 17652849c42SDave May PetscBool found; 177521f74f9SMatthew G. Knepley PetscErrorCode ierr; 17852849c42SDave May 179521f74f9SMatthew G. Knepley PetscFunctionBegin; 1802635f519SDave May ierr = StringInList(name,db->nfields,(const DataField*)db->field,&found);CHKERRQ(ierr); 1812eac95f8SDave May if (!found) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot find DataField with name %s",name); 1822635f519SDave May ierr = StringFindInList(name,db->nfields,(const DataField*)db->field,&idx);CHKERRQ(ierr); 18352849c42SDave May *gfield = db->field[idx]; 1842eac95f8SDave May PetscFunctionReturn(0); 18552849c42SDave May } 18652849c42SDave May 1877fbf63aeSDave May PetscErrorCode DataBucketQueryDataFieldByName(DataBucket db,const char name[],PetscBool *found) 18852849c42SDave May { 1892635f519SDave May PetscErrorCode ierr; 190521f74f9SMatthew G. Knepley 191521f74f9SMatthew G. Knepley PetscFunctionBegin; 19252849c42SDave May *found = PETSC_FALSE; 1932635f519SDave May ierr = StringInList(name,db->nfields,(const DataField*)db->field,found);CHKERRQ(ierr); 1947fbf63aeSDave May PetscFunctionReturn(0); 19552849c42SDave May } 19652849c42SDave May 1977fbf63aeSDave May PetscErrorCode DataBucketFinalize(DataBucket db) 19852849c42SDave May { 199521f74f9SMatthew G. Knepley PetscFunctionBegin; 20052849c42SDave May db->finalised = PETSC_TRUE; 2017fbf63aeSDave May PetscFunctionReturn(0); 20252849c42SDave May } 20352849c42SDave May 2047fbf63aeSDave May PetscErrorCode DataFieldGetNumEntries(DataField df,PetscInt *sum) 20552849c42SDave May { 206521f74f9SMatthew G. Knepley PetscFunctionBegin; 20752849c42SDave May *sum = df->L; 2087fbf63aeSDave May PetscFunctionReturn(0); 20952849c42SDave May } 21052849c42SDave May 21152c3ed93SDave May PetscErrorCode DataFieldSetBlockSize(DataField df,PetscInt blocksize) 21252c3ed93SDave May { 213521f74f9SMatthew G. Knepley PetscFunctionBegin; 21452c3ed93SDave May df->bs = blocksize; 21552c3ed93SDave May PetscFunctionReturn(0); 21652c3ed93SDave May } 21752c3ed93SDave May 2187fbf63aeSDave May PetscErrorCode DataFieldSetSize(DataField df,const PetscInt new_L) 21952849c42SDave May { 220521f74f9SMatthew G. Knepley PetscErrorCode ierr; 22152849c42SDave May 222521f74f9SMatthew G. Knepley PetscFunctionBegin; 223*c8d15862SDave May if (new_L < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot set size of DataField to be < 0"); 2247fbf63aeSDave May if (new_L == df->L) PetscFunctionReturn(0); 22552849c42SDave May if (new_L > df->L) { 2264be7464cSMatthew G. Knepley ierr = PetscRealloc(df->atomic_size * (new_L), &df->data);CHKERRQ(ierr); 22752849c42SDave May /* init new contents */ 228521f74f9SMatthew G. Knepley ierr = PetscMemzero(( ((char*)df->data)+df->L*df->atomic_size), (new_L-df->L)*df->atomic_size);CHKERRQ(ierr); 2297fbf63aeSDave May } else { 23052849c42SDave May /* reallocate pointer list, add +1 in case new_L = 0 */ 2314be7464cSMatthew G. Knepley ierr = PetscRealloc(df->atomic_size * (new_L+1), &df->data);CHKERRQ(ierr); 23252849c42SDave May } 23352849c42SDave May df->L = new_L; 2347fbf63aeSDave May PetscFunctionReturn(0); 23552849c42SDave May } 23652849c42SDave May 2377fbf63aeSDave May PetscErrorCode DataFieldZeroBlock(DataField df,const PetscInt start,const PetscInt end) 23852849c42SDave May { 239521f74f9SMatthew G. Knepley PetscErrorCode ierr; 240521f74f9SMatthew G. Knepley 241521f74f9SMatthew G. Knepley PetscFunctionBegin; 2427fbf63aeSDave May if (start > end) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if start(%D) > end(%D)",start,end); 2437fbf63aeSDave May if (start < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if start(%D) < 0",start); 244a233d522SDave 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); 245521f74f9SMatthew G. Knepley ierr = PetscMemzero((((char*)df->data)+start*df->atomic_size), (end-start)*df->atomic_size);CHKERRQ(ierr); 2467fbf63aeSDave May PetscFunctionReturn(0); 24752849c42SDave May } 24852849c42SDave May 24952849c42SDave May /* 25052849c42SDave May A negative buffer value will simply be ignored and the old buffer value will be used. 25152849c42SDave May */ 2527fbf63aeSDave May PetscErrorCode DataBucketSetSizes(DataBucket db,const PetscInt L,const PetscInt buffer) 25352849c42SDave May { 2545c18a9d6SDave May PetscInt current_allocated,new_used,new_unused,new_buffer,new_allocated,f; 2550cb17e37SDave May PetscBool any_active_fields; 256dbe06d34SDave May PetscErrorCode ierr; 25752849c42SDave May 258521f74f9SMatthew G. Knepley PetscFunctionBegin; 2597fbf63aeSDave May if (db->finalised == PETSC_FALSE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"You must call DataBucketFinalize() before DataBucketSetSizes()"); 2600cb17e37SDave May ierr = DataBucketQueryForActiveFields(db,&any_active_fields);CHKERRQ(ierr); 2610cb17e37SDave 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"); 2620cb17e37SDave May 26352849c42SDave May current_allocated = db->allocated; 26452849c42SDave May new_used = L; 26552849c42SDave May new_unused = current_allocated - new_used; 26652849c42SDave May new_buffer = db->buffer; 26752849c42SDave May if (buffer >= 0) { /* update the buffer value */ 26852849c42SDave May new_buffer = buffer; 26952849c42SDave May } 27052849c42SDave May new_allocated = new_used + new_buffer; 27152849c42SDave May /* action */ 27252849c42SDave May if (new_allocated > current_allocated) { 27352849c42SDave May /* increase size to new_used + new_buffer */ 27452849c42SDave May for (f=0; f<db->nfields; f++) { 275dbe06d34SDave May ierr = DataFieldSetSize(db->field[f], new_allocated);CHKERRQ(ierr); 27652849c42SDave May } 27752849c42SDave May db->L = new_used; 27852849c42SDave May db->buffer = new_buffer; 27952849c42SDave May db->allocated = new_used + new_buffer; 280521f74f9SMatthew G. Knepley } else { 28152849c42SDave May if (new_unused > 2 * new_buffer) { 28252849c42SDave May /* shrink array to new_used + new_buffer */ 283521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 284dbe06d34SDave May ierr = DataFieldSetSize(db->field[f], new_allocated);CHKERRQ(ierr); 28552849c42SDave May } 28652849c42SDave May db->L = new_used; 28752849c42SDave May db->buffer = new_buffer; 28852849c42SDave May db->allocated = new_used + new_buffer; 289521f74f9SMatthew G. Knepley } else { 29052849c42SDave May db->L = new_used; 29152849c42SDave May db->buffer = new_buffer; 29252849c42SDave May } 29352849c42SDave May } 29452849c42SDave May /* zero all entries from db->L to db->allocated */ 295521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 29652849c42SDave May DataField field = db->field[f]; 297dbe06d34SDave May ierr = DataFieldZeroBlock(field, db->L,db->allocated);CHKERRQ(ierr); 29852849c42SDave May } 2997fbf63aeSDave May PetscFunctionReturn(0); 30052849c42SDave May } 30152849c42SDave May 3027fbf63aeSDave May PetscErrorCode DataBucketSetInitialSizes(DataBucket db,const PetscInt L,const PetscInt buffer) 30352849c42SDave May { 3045c18a9d6SDave May PetscInt f; 305dbe06d34SDave May PetscErrorCode ierr; 306dbe06d34SDave May 307521f74f9SMatthew G. Knepley PetscFunctionBegin; 308dbe06d34SDave May ierr = DataBucketSetSizes(db,L,buffer);CHKERRQ(ierr); 309521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 31052849c42SDave May DataField field = db->field[f]; 311dbe06d34SDave May ierr = DataFieldZeroBlock(field,0,db->allocated);CHKERRQ(ierr); 31252849c42SDave May } 3137fbf63aeSDave May PetscFunctionReturn(0); 31452849c42SDave May } 31552849c42SDave May 3167fbf63aeSDave May PetscErrorCode DataBucketGetSizes(DataBucket db,PetscInt *L,PetscInt *buffer,PetscInt *allocated) 31752849c42SDave May { 318521f74f9SMatthew G. Knepley PetscFunctionBegin; 31952849c42SDave May if (L) {*L = db->L;} 32052849c42SDave May if (buffer) {*buffer = db->buffer;} 32152849c42SDave May if (allocated) {*allocated = db->allocated;} 3227fbf63aeSDave May PetscFunctionReturn(0); 32352849c42SDave May } 32452849c42SDave May 3257fbf63aeSDave May PetscErrorCode DataBucketGetGlobalSizes(MPI_Comm comm,DataBucket db,PetscInt *L,PetscInt *buffer,PetscInt *allocated) 32652849c42SDave May { 3275c18a9d6SDave May PetscInt _L,_buffer,_allocated; 3285c18a9d6SDave May PetscInt ierr; 32952849c42SDave May 330521f74f9SMatthew G. Knepley PetscFunctionBegin; 3315c18a9d6SDave May _L = db->L; 3325c18a9d6SDave May _buffer = db->buffer; 3335c18a9d6SDave May _allocated = db->allocated; 33452849c42SDave May 33533564166SDave May if (L) { ierr = MPI_Allreduce(&_L,L,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); } 33633564166SDave May if (buffer) { ierr = MPI_Allreduce(&_buffer,buffer,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); } 33733564166SDave May if (allocated) { ierr = MPI_Allreduce(&_allocated,allocated,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); } 3387fbf63aeSDave May PetscFunctionReturn(0); 33952849c42SDave May } 34052849c42SDave May 3417fbf63aeSDave May PetscErrorCode DataBucketGetDataFields(DataBucket db,PetscInt *L,DataField *fields[]) 34252849c42SDave May { 343521f74f9SMatthew G. Knepley PetscFunctionBegin; 34452849c42SDave May if (L) {*L = db->nfields;} 34552849c42SDave May if (fields) {*fields = db->field;} 3467fbf63aeSDave May PetscFunctionReturn(0); 34752849c42SDave May } 34852849c42SDave May 3497fbf63aeSDave May PetscErrorCode DataFieldGetAccess(const DataField gfield) 35052849c42SDave May { 351521f74f9SMatthew G. Knepley PetscFunctionBegin; 3527fbf63aeSDave May if (gfield->active) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is already active. You must call DataFieldRestoreAccess()",gfield->name); 35352849c42SDave May gfield->active = PETSC_TRUE; 3547fbf63aeSDave May PetscFunctionReturn(0); 35552849c42SDave May } 35652849c42SDave May 3577fbf63aeSDave May PetscErrorCode DataFieldAccessPoint(const DataField gfield,const PetscInt pid,void **ctx_p) 35852849c42SDave May { 359521f74f9SMatthew G. Knepley PetscFunctionBegin; 36072505a4dSJed Brown *ctx_p = NULL; 36152849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD 36252849c42SDave May /* debug mode */ 36384bcda08SDave May /* check point is valid */ 3647fbf63aeSDave May if (pid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0"); 3657fbf63aeSDave May if (pid >= gfield->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",gfield->L); 36684bcda08SDave 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); 36752849c42SDave May #endif 36852849c42SDave May *ctx_p = __DATATFIELD_point_access(gfield->data,pid,gfield->atomic_size); 3697fbf63aeSDave May PetscFunctionReturn(0); 37052849c42SDave May } 37152849c42SDave May 3727fbf63aeSDave May PetscErrorCode DataFieldAccessPointOffset(const DataField gfield,const size_t offset,const PetscInt pid,void **ctx_p) 37352849c42SDave May { 374521f74f9SMatthew G. Knepley PetscFunctionBegin; 37552849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD 37652849c42SDave May /* debug mode */ 37784bcda08SDave May /* check point is valid */ 378521f74f9SMatthew G. Knepley /* if( offset < 0 ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be >= 0");*/ 379521f74f9SMatthew G. Knepley /* Note compiler realizes this can never happen with an unsigned PetscInt */ 3807fbf63aeSDave May if (offset >= gfield->atomic_size) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be < %zu",gfield->atomic_size); 38184bcda08SDave May /* check point is valid */ 3827fbf63aeSDave May if (pid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0"); 383a233d522SDave May if (pid >= gfield->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",gfield->L); 384a233d522SDave 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); 38552849c42SDave May #endif 38652849c42SDave May *ctx_p = __DATATFIELD_point_access_offset(gfield->data,pid,gfield->atomic_size,offset); 3877fbf63aeSDave May PetscFunctionReturn(0); 38852849c42SDave May } 38952849c42SDave May 3907fbf63aeSDave May PetscErrorCode DataFieldRestoreAccess(DataField gfield) 39152849c42SDave May { 392521f74f9SMatthew G. Knepley PetscFunctionBegin; 3937fbf63aeSDave May if (gfield->active == PETSC_FALSE) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is not active. You must call DataFieldGetAccess()", gfield->name); 39452849c42SDave May gfield->active = PETSC_FALSE; 3957fbf63aeSDave May PetscFunctionReturn(0); 39652849c42SDave May } 39752849c42SDave May 3987fbf63aeSDave May PetscErrorCode DataFieldVerifyAccess(const DataField gfield,const size_t size) 39952849c42SDave May { 400521f74f9SMatthew G. Knepley PetscFunctionBegin; 40152849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD 402521f74f9SMatthew 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 ); 40352849c42SDave May #endif 4047fbf63aeSDave May PetscFunctionReturn(0); 40552849c42SDave May } 40652849c42SDave May 4077fbf63aeSDave May PetscErrorCode DataFieldGetAtomicSize(const DataField gfield,size_t *size) 40852849c42SDave May { 409521f74f9SMatthew G. Knepley PetscFunctionBegin; 41052849c42SDave May if (size) {*size = gfield->atomic_size;} 4117fbf63aeSDave May PetscFunctionReturn(0); 41252849c42SDave May } 41352849c42SDave May 4147fbf63aeSDave May PetscErrorCode DataFieldGetEntries(const DataField gfield,void **data) 41552849c42SDave May { 416521f74f9SMatthew G. Knepley PetscFunctionBegin; 417521f74f9SMatthew G. Knepley if (data) {*data = gfield->data;} 4187fbf63aeSDave May PetscFunctionReturn(0); 41952849c42SDave May } 42052849c42SDave May 4217fbf63aeSDave May PetscErrorCode DataFieldRestoreEntries(const DataField gfield,void **data) 42252849c42SDave May { 423521f74f9SMatthew G. Knepley PetscFunctionBegin; 424521f74f9SMatthew G. Knepley if (data) {*data = NULL;} 4257fbf63aeSDave May PetscFunctionReturn(0); 42652849c42SDave May } 42752849c42SDave May 42852849c42SDave May /* y = x */ 4297fbf63aeSDave May PetscErrorCode DataBucketCopyPoint(const DataBucket xb,const PetscInt pid_x, 4305c18a9d6SDave May const DataBucket yb,const PetscInt pid_y) 43152849c42SDave May { 4325c18a9d6SDave May PetscInt f; 433dbe06d34SDave May PetscErrorCode ierr; 434dbe06d34SDave May 435521f74f9SMatthew G. Knepley PetscFunctionBegin; 436521f74f9SMatthew G. Knepley for (f = 0; f < xb->nfields; ++f) { 43752849c42SDave May void *dest; 43852849c42SDave May void *src; 43952849c42SDave May 440dbe06d34SDave May ierr = DataFieldGetAccess(xb->field[f]);CHKERRQ(ierr); 4416845f8f5SDave May if (xb != yb) { ierr = DataFieldGetAccess( yb->field[f]);CHKERRQ(ierr); } 442dbe06d34SDave May ierr = DataFieldAccessPoint(xb->field[f],pid_x, &src);CHKERRQ(ierr); 443dbe06d34SDave May ierr = DataFieldAccessPoint(yb->field[f],pid_y, &dest);CHKERRQ(ierr); 444521f74f9SMatthew G. Knepley ierr = PetscMemcpy(dest, src, xb->field[f]->atomic_size);CHKERRQ(ierr); 445dbe06d34SDave May ierr = DataFieldRestoreAccess(xb->field[f]);CHKERRQ(ierr); 4466845f8f5SDave May if (xb != yb) {ierr = DataFieldRestoreAccess(yb->field[f]);CHKERRQ(ierr);} 44752849c42SDave May } 4487fbf63aeSDave May PetscFunctionReturn(0); 44952849c42SDave May } 45052849c42SDave May 4517fbf63aeSDave May PetscErrorCode DataBucketCreateFromSubset(DataBucket DBIn,const PetscInt N,const PetscInt list[],DataBucket *DB) 45252849c42SDave May { 4535c18a9d6SDave May PetscInt nfields; 45452849c42SDave May DataField *fields; 4555c18a9d6SDave May PetscInt f,L,buffer,allocated,p; 456dbe06d34SDave May PetscErrorCode ierr; 45752849c42SDave May 458521f74f9SMatthew G. Knepley PetscFunctionBegin; 459521f74f9SMatthew G. Knepley ierr = DataBucketCreate(DB);CHKERRQ(ierr); 46052849c42SDave May /* copy contents of DBIn */ 461dbe06d34SDave May ierr = DataBucketGetDataFields(DBIn,&nfields,&fields);CHKERRQ(ierr); 462dbe06d34SDave May ierr = DataBucketGetSizes(DBIn,&L,&buffer,&allocated);CHKERRQ(ierr); 463521f74f9SMatthew G. Knepley for (f = 0; f < nfields; ++f) { 464dbe06d34SDave May ierr = DataBucketRegisterField(*DB,"DataBucketCreateFromSubset",fields[f]->name,fields[f]->atomic_size,NULL);CHKERRQ(ierr); 46552849c42SDave May } 466dbe06d34SDave May ierr = DataBucketFinalize(*DB);CHKERRQ(ierr); 467dbe06d34SDave May ierr = DataBucketSetSizes(*DB,L,buffer);CHKERRQ(ierr); 46852849c42SDave May /* now copy the desired guys from DBIn => DB */ 469521f74f9SMatthew G. Knepley for (p = 0; p < N; ++p) { 470dbe06d34SDave May ierr = DataBucketCopyPoint(DBIn,list[p], *DB,p);CHKERRQ(ierr); 47152849c42SDave May } 4727fbf63aeSDave May PetscFunctionReturn(0); 47352849c42SDave May } 47452849c42SDave May 475521f74f9SMatthew G. Knepley /* insert into an exisitng location */ 4767fbf63aeSDave May PetscErrorCode DataFieldInsertPoint(const DataField field,const PetscInt index,const void *ctx) 47752849c42SDave May { 478521f74f9SMatthew G. Knepley PetscErrorCode ierr; 47952849c42SDave May 480521f74f9SMatthew G. Knepley PetscFunctionBegin; 48152849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD 48284bcda08SDave May /* check point is valid */ 483a233d522SDave May if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0"); 484a233d522SDave May if (index >= field->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",field->L); 48552849c42SDave May #endif 486521f74f9SMatthew G. Knepley ierr = PetscMemcpy(__DATATFIELD_point_access(field->data,index,field->atomic_size), ctx, field->atomic_size);CHKERRQ(ierr); 4877fbf63aeSDave May PetscFunctionReturn(0); 48852849c42SDave May } 48952849c42SDave May 490521f74f9SMatthew G. Knepley /* remove data at index - replace with last point */ 4917fbf63aeSDave May PetscErrorCode DataBucketRemovePointAtIndex(const DataBucket db,const PetscInt index) 49252849c42SDave May { 4935c18a9d6SDave May PetscInt f; 4940cb17e37SDave May PetscBool any_active_fields; 495dbe06d34SDave May PetscErrorCode ierr; 49652849c42SDave May 497521f74f9SMatthew G. Knepley PetscFunctionBegin; 49852849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD 49984bcda08SDave May /* check point is valid */ 500a233d522SDave May if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0"); 501a233d522SDave May if (index >= db->allocated) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",db->L+db->buffer); 50252849c42SDave May #endif 5030cb17e37SDave May ierr = DataBucketQueryForActiveFields(db,&any_active_fields);CHKERRQ(ierr); 5040cb17e37SDave 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"); 505a233d522SDave May if (index >= db->L) { /* this point is not in the list - no need to error, but I will anyway */ 506a233d522SDave 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); 50752849c42SDave May } 508a233d522SDave May if (index != db->L-1) { /* not last point in list */ 509521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 51052849c42SDave May DataField field = db->field[f]; 51152849c42SDave May 51252849c42SDave May /* copy then remove */ 513dbe06d34SDave May ierr = DataFieldCopyPoint(db->L-1, field, index, field);CHKERRQ(ierr); 514521f74f9SMatthew G. Knepley /* DataFieldZeroPoint(field,index); */ 51552849c42SDave May } 51652849c42SDave May } 51752849c42SDave May /* decrement size */ 51852849c42SDave May /* this will zero out an crap at the end of the list */ 519dbe06d34SDave May ierr = DataBucketRemovePoint(db);CHKERRQ(ierr); 5207fbf63aeSDave May PetscFunctionReturn(0); 52152849c42SDave May } 52252849c42SDave May 52352849c42SDave May /* copy x into y */ 5247fbf63aeSDave May PetscErrorCode DataFieldCopyPoint(const PetscInt pid_x,const DataField field_x, 5255c18a9d6SDave May const PetscInt pid_y,const DataField field_y ) 52652849c42SDave May { 527521f74f9SMatthew G. Knepley PetscErrorCode ierr; 52852849c42SDave May 529521f74f9SMatthew G. Knepley PetscFunctionBegin; 53052849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD 53184bcda08SDave May /* check point is valid */ 532a233d522SDave May if (pid_x < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"(IN) index must be >= 0"); 533a233d522SDave May if (pid_x >= field_x->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"(IN) index must be < %D",field_x->L); 534a233d522SDave May if (pid_y < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"(OUT) index must be >= 0"); 535a233d522SDave May if (pid_y >= field_y->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"(OUT) index must be < %D",field_y->L); 536521f74f9SMatthew G. Knepley if( field_y->atomic_size != field_x->atomic_size ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"atomic size must match"); 53752849c42SDave May #endif 538521f74f9SMatthew G. Knepley ierr = PetscMemcpy(__DATATFIELD_point_access(field_y->data,pid_y,field_y->atomic_size), 53952849c42SDave May __DATATFIELD_point_access(field_x->data,pid_x,field_x->atomic_size), 540521f74f9SMatthew G. Knepley field_y->atomic_size);CHKERRQ(ierr); 5417fbf63aeSDave May PetscFunctionReturn(0); 54252849c42SDave May } 54352849c42SDave May 54452849c42SDave May 545521f74f9SMatthew G. Knepley /* zero only the datafield at this point */ 5467fbf63aeSDave May PetscErrorCode DataFieldZeroPoint(const DataField field,const PetscInt index) 54752849c42SDave May { 548521f74f9SMatthew G. Knepley PetscErrorCode ierr; 549521f74f9SMatthew G. Knepley 550521f74f9SMatthew G. Knepley PetscFunctionBegin; 55152849c42SDave May #ifdef DATAFIELD_POINT_ACCESS_GUARD 55284bcda08SDave May /* check point is valid */ 553a233d522SDave May if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0"); 554a233d522SDave May if (index >= field->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",field->L); 55552849c42SDave May #endif 556521f74f9SMatthew G. Knepley ierr = PetscMemzero(__DATATFIELD_point_access(field->data,index,field->atomic_size), field->atomic_size);CHKERRQ(ierr); 5577fbf63aeSDave May PetscFunctionReturn(0); 55852849c42SDave May } 55952849c42SDave May 560521f74f9SMatthew G. Knepley /* zero ALL data for this point */ 5617fbf63aeSDave May PetscErrorCode DataBucketZeroPoint(const DataBucket db,const PetscInt index) 56252849c42SDave May { 5635c18a9d6SDave May PetscInt f; 564dbe06d34SDave May PetscErrorCode ierr; 56552849c42SDave May 566521f74f9SMatthew G. Knepley PetscFunctionBegin; 56784bcda08SDave May /* check point is valid */ 568a233d522SDave May if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0"); 569a233d522SDave May if (index >= db->allocated) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",db->allocated); 570521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 57152849c42SDave May DataField field = db->field[f]; 572dbe06d34SDave May ierr = DataFieldZeroPoint(field,index);CHKERRQ(ierr); 57352849c42SDave May } 5747fbf63aeSDave May PetscFunctionReturn(0); 57552849c42SDave May } 57652849c42SDave May 57752849c42SDave May /* increment */ 5787fbf63aeSDave May PetscErrorCode DataBucketAddPoint(DataBucket db) 57952849c42SDave May { 580dbe06d34SDave May PetscErrorCode ierr; 581dbe06d34SDave May 582521f74f9SMatthew G. Knepley PetscFunctionBegin; 583dbe06d34SDave May ierr = DataBucketSetSizes(db,db->L+1,-1);CHKERRQ(ierr); 5847fbf63aeSDave May PetscFunctionReturn(0); 58552849c42SDave May } 58652849c42SDave May 5877fbf63aeSDave May /* decrement */ 5887fbf63aeSDave May PetscErrorCode DataBucketRemovePoint(DataBucket db) 5897fbf63aeSDave May { 590dbe06d34SDave May PetscErrorCode ierr; 591dbe06d34SDave May 592521f74f9SMatthew G. Knepley PetscFunctionBegin; 593dbe06d34SDave May ierr = DataBucketSetSizes(db,db->L-1,-1);CHKERRQ(ierr); 5947fbf63aeSDave May PetscFunctionReturn(0); 5957fbf63aeSDave May } 5967fbf63aeSDave May 597ab2be9a3SDave May PetscErrorCode DataBucketView_stdout(MPI_Comm comm,DataBucket db) 598ab2be9a3SDave May { 599ab2be9a3SDave May PetscInt f; 600ab2be9a3SDave May double memory_usage_total,memory_usage_total_local = 0.0; 601ab2be9a3SDave May PetscErrorCode ierr; 602ab2be9a3SDave May 603ab2be9a3SDave May PetscFunctionBegin; 604ab2be9a3SDave May ierr = PetscPrintf(comm,"DataBucketView: \n");CHKERRQ(ierr); 605ab2be9a3SDave May ierr = PetscPrintf(comm," L = %D \n", db->L);CHKERRQ(ierr); 606ab2be9a3SDave May ierr = PetscPrintf(comm," buffer = %D \n", db->buffer);CHKERRQ(ierr); 607ab2be9a3SDave May ierr = PetscPrintf(comm," allocated = %D \n", db->allocated);CHKERRQ(ierr); 608ab2be9a3SDave May ierr = PetscPrintf(comm," nfields registered = %D \n", db->nfields);CHKERRQ(ierr); 609ab2be9a3SDave May 610ab2be9a3SDave May for (f = 0; f < db->nfields; ++f) { 611ab2be9a3SDave May double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6; 612ab2be9a3SDave May memory_usage_total_local += memory_usage_f; 613ab2be9a3SDave May } 614ab2be9a3SDave May ierr = MPI_Allreduce(&memory_usage_total_local,&memory_usage_total,1,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr); 615ab2be9a3SDave May 616ab2be9a3SDave May for (f = 0; f < db->nfields; ++f) { 617ab2be9a3SDave May double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6; 618ab2be9a3SDave May ierr = PetscPrintf(comm," [%3D] %15s : Mem. usage = %1.2e (MB) [rank0]\n", f, db->field[f]->name, memory_usage_f );CHKERRQ(ierr); 619ab2be9a3SDave May ierr = PetscPrintf(comm," blocksize = %D \n", db->field[f]->bs);CHKERRQ(ierr); 620ab2be9a3SDave May if (db->field[f]->bs != 1) { 621ab2be9a3SDave May ierr = PetscPrintf(comm," atomic size = %zu [full block, bs=%D]\n", db->field[f]->atomic_size,db->field[f]->bs);CHKERRQ(ierr); 622ab2be9a3SDave May ierr = PetscPrintf(comm," atomic size/item = %zu \n", db->field[f]->atomic_size/db->field[f]->bs);CHKERRQ(ierr); 623ab2be9a3SDave May } else { 624ab2be9a3SDave May ierr = PetscPrintf(comm," atomic size = %zu \n", db->field[f]->atomic_size);CHKERRQ(ierr); 625ab2be9a3SDave May } 626ab2be9a3SDave May } 627ab2be9a3SDave May ierr = PetscPrintf(comm," Total mem. usage = %1.2e (MB) (collective)\n", memory_usage_total);CHKERRQ(ierr); 628ab2be9a3SDave May PetscFunctionReturn(0); 629ab2be9a3SDave May } 630ab2be9a3SDave May 631ab2be9a3SDave May PetscErrorCode DataBucketView_SEQ(MPI_Comm comm,DataBucket db,const char filename[],DataBucketViewType type) 63252849c42SDave May { 633dbe06d34SDave May PetscErrorCode ierr; 634dbe06d34SDave May 635521f74f9SMatthew G. Knepley PetscFunctionBegin; 63652849c42SDave May switch (type) { 63752849c42SDave May case DATABUCKET_VIEW_STDOUT: 638ab2be9a3SDave May ierr = DataBucketView_stdout(PETSC_COMM_SELF,db);CHKERRQ(ierr); 63952849c42SDave May break; 64052849c42SDave May case DATABUCKET_VIEW_ASCII: 641071ce369SDave May SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for ascii output"); 64252849c42SDave May break; 64352849c42SDave May case DATABUCKET_VIEW_BINARY: 644071ce369SDave May SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for binary output"); 64552849c42SDave May break; 64652849c42SDave May case DATABUCKET_VIEW_HDF5: 647071ce369SDave May SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for HDF5 output"); 64852849c42SDave May break; 649521f74f9SMatthew G. Knepley default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer method requested"); 65052849c42SDave May } 6517fbf63aeSDave May PetscFunctionReturn(0); 65252849c42SDave May } 65352849c42SDave May 6547fbf63aeSDave May PetscErrorCode DataBucketView_MPI(MPI_Comm comm,DataBucket db,const char filename[],DataBucketViewType type) 65552849c42SDave May { 656997fa542SDave May PetscErrorCode ierr; 657997fa542SDave May 658521f74f9SMatthew G. Knepley PetscFunctionBegin; 65952849c42SDave May switch (type) { 66052849c42SDave May case DATABUCKET_VIEW_STDOUT: 661ab2be9a3SDave May ierr = DataBucketView_stdout(comm,db);CHKERRQ(ierr); 66252849c42SDave May break; 66352849c42SDave May case DATABUCKET_VIEW_ASCII: 664071ce369SDave May SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for ascii output"); 66552849c42SDave May break; 66652849c42SDave May case DATABUCKET_VIEW_BINARY: 667071ce369SDave May SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for binary output"); 66852849c42SDave May break; 66952849c42SDave May case DATABUCKET_VIEW_HDF5: 670071ce369SDave May SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for HDF5 output"); 67152849c42SDave May break; 672521f74f9SMatthew G. Knepley default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer method requested"); 67352849c42SDave May } 6747fbf63aeSDave May PetscFunctionReturn(0); 67552849c42SDave May } 67652849c42SDave May 6777fbf63aeSDave May PetscErrorCode DataBucketView(MPI_Comm comm,DataBucket db,const char filename[],DataBucketViewType type) 67852849c42SDave May { 6795c18a9d6SDave May PetscMPIInt nproc; 680997fa542SDave May PetscErrorCode ierr; 68152849c42SDave May 682521f74f9SMatthew G. Knepley PetscFunctionBegin; 683997fa542SDave May ierr = MPI_Comm_size(comm,&nproc);CHKERRQ(ierr); 68452849c42SDave May if (nproc == 1) { 685ab2be9a3SDave May ierr = DataBucketView_SEQ(comm,db,filename,type);CHKERRQ(ierr); 68652849c42SDave May } else { 687dbe06d34SDave May ierr = DataBucketView_MPI(comm,db,filename,type);CHKERRQ(ierr); 68852849c42SDave May } 6897fbf63aeSDave May PetscFunctionReturn(0); 69052849c42SDave May } 69152849c42SDave May 6927fbf63aeSDave May PetscErrorCode DataBucketDuplicateFields(DataBucket dbA,DataBucket *dbB) 69352849c42SDave May { 69452849c42SDave May DataBucket db2; 6955c18a9d6SDave May PetscInt f; 696dbe06d34SDave May PetscErrorCode ierr; 69752849c42SDave May 698521f74f9SMatthew G. Knepley PetscFunctionBegin; 699dbe06d34SDave May ierr = DataBucketCreate(&db2);CHKERRQ(ierr); 70052849c42SDave May /* copy contents from dbA into db2 */ 701521f74f9SMatthew G. Knepley for (f = 0; f < dbA->nfields; ++f) { 70252849c42SDave May DataField field; 70352849c42SDave May size_t atomic_size; 70452849c42SDave May char *name; 70552849c42SDave May 70652849c42SDave May field = dbA->field[f]; 70752849c42SDave May atomic_size = field->atomic_size; 70852849c42SDave May name = field->name; 709dbe06d34SDave May ierr = DataBucketRegisterField(db2,"DataBucketDuplicateFields",name,atomic_size,NULL);CHKERRQ(ierr); 71052849c42SDave May } 711dbe06d34SDave May ierr = DataBucketFinalize(db2);CHKERRQ(ierr); 712dbe06d34SDave May ierr = DataBucketSetInitialSizes(db2,0,1000);CHKERRQ(ierr); 71352849c42SDave May *dbB = db2; 7147fbf63aeSDave May PetscFunctionReturn(0); 71552849c42SDave May } 71652849c42SDave May 71752849c42SDave May /* 71852849c42SDave May Insert points from db2 into db1 71952849c42SDave May db1 <<== db2 72052849c42SDave May */ 7217fbf63aeSDave May PetscErrorCode DataBucketInsertValues(DataBucket db1,DataBucket db2) 72252849c42SDave May { 7235c18a9d6SDave May PetscInt n_mp_points1,n_mp_points2; 7245c18a9d6SDave May PetscInt n_mp_points1_new,p; 725dbe06d34SDave May PetscErrorCode ierr; 72652849c42SDave May 727521f74f9SMatthew G. Knepley PetscFunctionBegin; 728dbe06d34SDave May ierr = DataBucketGetSizes(db1,&n_mp_points1,0,0);CHKERRQ(ierr); 729dbe06d34SDave May ierr = DataBucketGetSizes(db2,&n_mp_points2,0,0);CHKERRQ(ierr); 73052849c42SDave May n_mp_points1_new = n_mp_points1 + n_mp_points2; 731dbe06d34SDave May ierr = DataBucketSetSizes(db1,n_mp_points1_new,-1);CHKERRQ(ierr); 732521f74f9SMatthew G. Knepley for (p = 0; p < n_mp_points2; ++p) { 733521f74f9SMatthew G. Knepley /* db1 <<== db2 */ 734dbe06d34SDave May ierr = DataBucketCopyPoint(db2,p, db1,(n_mp_points1 + p));CHKERRQ(ierr); 73552849c42SDave May } 7367fbf63aeSDave May PetscFunctionReturn(0); 73752849c42SDave May } 73852849c42SDave May 73952849c42SDave May /* helpers for parallel send/recv */ 7407fbf63aeSDave May PetscErrorCode DataBucketCreatePackedArray(DataBucket db,size_t *bytes,void **buf) 74152849c42SDave May { 7425c18a9d6SDave May PetscInt f; 74352849c42SDave May size_t sizeof_marker_contents; 74452849c42SDave May void *buffer; 745521f74f9SMatthew G. Knepley PetscErrorCode ierr; 74652849c42SDave May 747521f74f9SMatthew G. Knepley PetscFunctionBegin; 74852849c42SDave May sizeof_marker_contents = 0; 749521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 75052849c42SDave May DataField df = db->field[f]; 75152849c42SDave May sizeof_marker_contents += df->atomic_size; 75252849c42SDave May } 753521f74f9SMatthew G. Knepley ierr = PetscMalloc(sizeof_marker_contents, &buffer);CHKERRQ(ierr); 754521f74f9SMatthew G. Knepley ierr = PetscMemzero(buffer, sizeof_marker_contents);CHKERRQ(ierr); 75552849c42SDave May if (bytes) {*bytes = sizeof_marker_contents;} 75652849c42SDave May if (buf) {*buf = buffer;} 7577fbf63aeSDave May PetscFunctionReturn(0); 75852849c42SDave May } 75952849c42SDave May 7607fbf63aeSDave May PetscErrorCode DataBucketDestroyPackedArray(DataBucket db,void **buf) 76152849c42SDave May { 762521f74f9SMatthew G. Knepley PetscErrorCode ierr; 763521f74f9SMatthew G. Knepley 764521f74f9SMatthew G. Knepley PetscFunctionBegin; 76552849c42SDave May if (buf) { 766521f74f9SMatthew G. Knepley ierr = PetscFree(*buf);CHKERRQ(ierr); 76752849c42SDave May *buf = NULL; 76852849c42SDave May } 7697fbf63aeSDave May PetscFunctionReturn(0); 77052849c42SDave May } 77152849c42SDave May 7727fbf63aeSDave May PetscErrorCode DataBucketFillPackedArray(DataBucket db,const PetscInt index,void *buf) 77352849c42SDave May { 7745c18a9d6SDave May PetscInt f; 77552849c42SDave May void *data, *data_p; 77652849c42SDave May size_t asize, offset; 777521f74f9SMatthew G. Knepley PetscErrorCode ierr; 77852849c42SDave May 779521f74f9SMatthew G. Knepley PetscFunctionBegin; 78052849c42SDave May offset = 0; 781521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 78252849c42SDave May DataField df = db->field[f]; 78352849c42SDave May 78452849c42SDave May asize = df->atomic_size; 78552849c42SDave May data = (void*)( df->data ); 78652849c42SDave May data_p = (void*)( (char*)data + index*asize ); 787521f74f9SMatthew G. Knepley ierr = PetscMemcpy((void*)((char*)buf + offset), data_p, asize);CHKERRQ(ierr); 78852849c42SDave May offset = offset + asize; 78952849c42SDave May } 7907fbf63aeSDave May PetscFunctionReturn(0); 79152849c42SDave May } 79252849c42SDave May 7937fbf63aeSDave May PetscErrorCode DataBucketInsertPackedArray(DataBucket db,const PetscInt idx,void *data) 79452849c42SDave May { 7955c18a9d6SDave May PetscInt f; 79652849c42SDave May void *data_p; 79752849c42SDave May size_t offset; 798dbe06d34SDave May PetscErrorCode ierr; 79952849c42SDave May 800521f74f9SMatthew G. Knepley PetscFunctionBegin; 80152849c42SDave May offset = 0; 802521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 80352849c42SDave May DataField df = db->field[f]; 80452849c42SDave May 80552849c42SDave May data_p = (void*)( (char*)data + offset ); 806dbe06d34SDave May ierr = DataFieldInsertPoint(df, idx, (void*)data_p);CHKERRQ(ierr); 80752849c42SDave May offset = offset + df->atomic_size; 80852849c42SDave May } 8097fbf63aeSDave May PetscFunctionReturn(0); 81052849c42SDave May } 811