1279f676cSBarry Smith #include "../src/dm/impls/swarm/data_bucket.h" 252849c42SDave May 352849c42SDave May /* string helpers */ 477048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldStringInList(const char name[],const PetscInt N,const DMSwarmDataField gfield[],PetscBool *val) 552849c42SDave May { 65c18a9d6SDave May PetscInt i; 752849c42SDave May 8521f74f9SMatthew G. Knepley PetscFunctionBegin; 952849c42SDave May *val = PETSC_FALSE; 10521f74f9SMatthew G. Knepley for (i = 0; i < N; ++i) { 11521f74f9SMatthew G. Knepley PetscBool flg; 129566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(name, gfield[i]->name, &flg)); 13521f74f9SMatthew G. Knepley if (flg) { 1452849c42SDave May *val = PETSC_TRUE; 152eac95f8SDave May PetscFunctionReturn(0); 1652849c42SDave May } 1752849c42SDave May } 182eac95f8SDave May PetscFunctionReturn(0); 1952849c42SDave May } 2052849c42SDave May 2177048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldStringFindInList(const char name[],const PetscInt N,const DMSwarmDataField gfield[],PetscInt *index) 2252849c42SDave May { 235c18a9d6SDave May PetscInt i; 2452849c42SDave May 25521f74f9SMatthew G. Knepley PetscFunctionBegin; 2652849c42SDave May *index = -1; 27521f74f9SMatthew G. Knepley for (i = 0; i < N; ++i) { 28521f74f9SMatthew G. Knepley PetscBool flg; 299566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(name, gfield[i]->name, &flg)); 30521f74f9SMatthew G. Knepley if (flg) { 3152849c42SDave May *index = i; 322eac95f8SDave May PetscFunctionReturn(0); 3352849c42SDave May } 3452849c42SDave May } 352eac95f8SDave May PetscFunctionReturn(0); 3652849c42SDave May } 3752849c42SDave May 38ee71fbaeSPatrick Sanan PetscErrorCode DMSwarmDataFieldCreate(const char registration_function[],const char name[],const size_t size,const PetscInt L,DMSwarmDataField *DF) 3952849c42SDave May { 4077048351SPatrick Sanan DMSwarmDataField df; 4152849c42SDave May 42521f74f9SMatthew G. Knepley PetscFunctionBegin; 439566063dSJacob Faibussowitsch PetscCall(PetscNew(&df)); 449566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(registration_function, &df->registration_function)); 459566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &df->name)); 4652849c42SDave May df->atomic_size = size; 4752849c42SDave May df->L = L; 4852c3ed93SDave May df->bs = 1; 49521f74f9SMatthew G. Knepley /* allocate something so we don't have to reallocate */ 509566063dSJacob Faibussowitsch PetscCall(PetscMalloc(size * L, &df->data)); 519566063dSJacob Faibussowitsch PetscCall(PetscMemzero(df->data, size * L)); 5252849c42SDave May *DF = df; 532eac95f8SDave May PetscFunctionReturn(0); 5452849c42SDave May } 5552849c42SDave May 5677048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldDestroy(DMSwarmDataField *DF) 5752849c42SDave May { 5877048351SPatrick Sanan DMSwarmDataField df = *DF; 5952849c42SDave May 60521f74f9SMatthew G. Knepley PetscFunctionBegin; 619566063dSJacob Faibussowitsch PetscCall(PetscFree(df->registration_function)); 629566063dSJacob Faibussowitsch PetscCall(PetscFree(df->name)); 639566063dSJacob Faibussowitsch PetscCall(PetscFree(df->data)); 649566063dSJacob Faibussowitsch PetscCall(PetscFree(df)); 6552849c42SDave May *DF = NULL; 662eac95f8SDave May PetscFunctionReturn(0); 6752849c42SDave May } 6852849c42SDave May 6952849c42SDave May /* data bucket */ 7077048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketCreate(DMSwarmDataBucket *DB) 7152849c42SDave May { 7277048351SPatrick Sanan DMSwarmDataBucket db; 7352849c42SDave May 74521f74f9SMatthew G. Knepley PetscFunctionBegin; 759566063dSJacob Faibussowitsch PetscCall(PetscNew(&db)); 7652849c42SDave May 7752849c42SDave May db->finalised = PETSC_FALSE; 7852849c42SDave May /* create empty spaces for fields */ 793454631fSDave May db->L = -1; 8052849c42SDave May db->buffer = 1; 8152849c42SDave May db->allocated = 1; 8252849c42SDave May db->nfields = 0; 839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1, &db->field)); 8452849c42SDave May *DB = db; 852eac95f8SDave May PetscFunctionReturn(0); 8652849c42SDave May } 8752849c42SDave May 8877048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketDestroy(DMSwarmDataBucket *DB) 8952849c42SDave May { 9077048351SPatrick Sanan DMSwarmDataBucket db = *DB; 915c18a9d6SDave May PetscInt f; 9252849c42SDave May 93521f74f9SMatthew G. Knepley PetscFunctionBegin; 9452849c42SDave May /* release fields */ 95521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 969566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldDestroy(&db->field[f])); 9752849c42SDave May } 9852849c42SDave May /* this will catch the initially allocated objects in the event that no fields are registered */ 9952849c42SDave May if (db->field != NULL) { 1009566063dSJacob Faibussowitsch PetscCall(PetscFree(db->field)); 10152849c42SDave May } 1029566063dSJacob Faibussowitsch PetscCall(PetscFree(db)); 10352849c42SDave May *DB = NULL; 1042eac95f8SDave May PetscFunctionReturn(0); 10552849c42SDave May } 10652849c42SDave May 10777048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketQueryForActiveFields(DMSwarmDataBucket db,PetscBool *any_active_fields) 1080cb17e37SDave May { 1090cb17e37SDave May PetscInt f; 110521f74f9SMatthew G. Knepley 111521f74f9SMatthew G. Knepley PetscFunctionBegin; 1120cb17e37SDave May *any_active_fields = PETSC_FALSE; 113521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 1140cb17e37SDave May if (db->field[f]->active) { 1150cb17e37SDave May *any_active_fields = PETSC_TRUE; 1160cb17e37SDave May PetscFunctionReturn(0); 1170cb17e37SDave May } 1180cb17e37SDave May } 1190cb17e37SDave May PetscFunctionReturn(0); 1200cb17e37SDave May } 1210cb17e37SDave May 1225627991aSBarry Smith PetscErrorCode DMSwarmDataBucketRegisterField(DMSwarmDataBucket db,const char registration_function[],const char field_name[],size_t atomic_size, DMSwarmDataField *_gfield) 12352849c42SDave May { 12452849c42SDave May PetscBool val; 12577048351SPatrick Sanan DMSwarmDataField fp; 12652849c42SDave May 127521f74f9SMatthew G. Knepley PetscFunctionBegin; 12852849c42SDave May /* check we haven't finalised the registration of fields */ 12952849c42SDave May /* 13052849c42SDave May if (db->finalised==PETSC_TRUE) { 13177048351SPatrick Sanan printf("ERROR: DMSwarmDataBucketFinalize() has been called. Cannot register more fields\n"); 13252849c42SDave May ERROR(); 13352849c42SDave May } 13452849c42SDave May */ 13552849c42SDave May /* check for repeated name */ 1369566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldStringInList(field_name, db->nfields, (const DMSwarmDataField*) db->field, &val)); 13708401ef6SPierre Jolivet PetscCheck(val != PETSC_TRUE,PETSC_COMM_SELF,PETSC_ERR_USER,"Field %s already exists. Cannot add same field twice",field_name); 13852849c42SDave May /* create new space for data */ 1399566063dSJacob Faibussowitsch PetscCall(PetscRealloc(sizeof(DMSwarmDataField)*(db->nfields+1), &db->field)); 14052849c42SDave May /* add field */ 1419566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldCreate(registration_function, field_name, atomic_size, db->allocated, &fp)); 14252849c42SDave May db->field[db->nfields] = fp; 14352849c42SDave May db->nfields++; 14452849c42SDave May if (_gfield != NULL) { 14552849c42SDave May *_gfield = fp; 14652849c42SDave May } 1472eac95f8SDave May PetscFunctionReturn(0); 14852849c42SDave May } 14952849c42SDave May 15052849c42SDave May /* 15177048351SPatrick Sanan #define DMSwarmDataBucketRegisterField(db,name,size,k) {\ 15252849c42SDave May char *location;\ 15352849c42SDave May asprintf(&location,"Registered by %s() at line %d within file %s", __FUNCTION__, __LINE__, __FILE__);\ 15477048351SPatrick Sanan _DMSwarmDataBucketRegisterField( (db), location, (name), (size), (k));\ 155521f74f9SMatthew G. Knepley ierr = PetscFree(location);\ 15652849c42SDave May } 15752849c42SDave May */ 15852849c42SDave May 159*2e956fe4SStefano Zampini PetscErrorCode DMSwarmDataBucketGetDMSwarmDataFieldIdByName(DMSwarmDataBucket db,const char name[],PetscInt *idx) 160*2e956fe4SStefano Zampini { 161*2e956fe4SStefano Zampini PetscBool found; 162*2e956fe4SStefano Zampini 163*2e956fe4SStefano Zampini PetscFunctionBegin; 164*2e956fe4SStefano Zampini *idx = -1; 165*2e956fe4SStefano Zampini PetscCall(DMSwarmDataFieldStringInList(name,db->nfields,(const DMSwarmDataField*)db->field,&found)); 166*2e956fe4SStefano Zampini PetscCheck(found,PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot find DMSwarmDataField with name %s",name); 167*2e956fe4SStefano Zampini PetscCall(DMSwarmDataFieldStringFindInList(name,db->nfields,(const DMSwarmDataField*)db->field,idx)); 168*2e956fe4SStefano Zampini PetscFunctionReturn(0); 169*2e956fe4SStefano Zampini } 170*2e956fe4SStefano Zampini 17177048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketGetDMSwarmDataFieldByName(DMSwarmDataBucket db,const char name[],DMSwarmDataField *gfield) 17252849c42SDave May { 1735c18a9d6SDave May PetscInt idx; 17452849c42SDave May PetscBool found; 17552849c42SDave May 176521f74f9SMatthew G. Knepley PetscFunctionBegin; 1779566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldStringInList(name,db->nfields,(const DMSwarmDataField*)db->field,&found)); 17828b400f6SJacob Faibussowitsch PetscCheck(found,PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot find DMSwarmDataField with name %s",name); 1799566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldStringFindInList(name,db->nfields,(const DMSwarmDataField*)db->field,&idx)); 18052849c42SDave May *gfield = db->field[idx]; 1812eac95f8SDave May PetscFunctionReturn(0); 18252849c42SDave May } 18352849c42SDave May 18477048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketQueryDMSwarmDataFieldByName(DMSwarmDataBucket db,const char name[],PetscBool *found) 18552849c42SDave May { 186521f74f9SMatthew G. Knepley PetscFunctionBegin; 18752849c42SDave May *found = PETSC_FALSE; 1889566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldStringInList(name,db->nfields,(const DMSwarmDataField*)db->field,found)); 1897fbf63aeSDave May PetscFunctionReturn(0); 19052849c42SDave May } 19152849c42SDave May 19277048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketFinalize(DMSwarmDataBucket db) 19352849c42SDave May { 194521f74f9SMatthew G. Knepley PetscFunctionBegin; 19552849c42SDave May db->finalised = PETSC_TRUE; 1967fbf63aeSDave May PetscFunctionReturn(0); 19752849c42SDave May } 19852849c42SDave May 19977048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldGetNumEntries(DMSwarmDataField df,PetscInt *sum) 20052849c42SDave May { 201521f74f9SMatthew G. Knepley PetscFunctionBegin; 20252849c42SDave May *sum = df->L; 2037fbf63aeSDave May PetscFunctionReturn(0); 20452849c42SDave May } 20552849c42SDave May 20677048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldSetBlockSize(DMSwarmDataField df,PetscInt blocksize) 20752c3ed93SDave May { 208521f74f9SMatthew G. Knepley PetscFunctionBegin; 20952c3ed93SDave May df->bs = blocksize; 21052c3ed93SDave May PetscFunctionReturn(0); 21152c3ed93SDave May } 21252c3ed93SDave May 21377048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldSetSize(DMSwarmDataField df,const PetscInt new_L) 21452849c42SDave May { 215521f74f9SMatthew G. Knepley PetscFunctionBegin; 21608401ef6SPierre Jolivet PetscCheck(new_L >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot set size of DMSwarmDataField to be < 0"); 2177fbf63aeSDave May if (new_L == df->L) PetscFunctionReturn(0); 21852849c42SDave May if (new_L > df->L) { 2199566063dSJacob Faibussowitsch PetscCall(PetscRealloc(df->atomic_size * (new_L), &df->data)); 22052849c42SDave May /* init new contents */ 2219566063dSJacob Faibussowitsch PetscCall(PetscMemzero(( ((char*)df->data)+df->L*df->atomic_size), (new_L-df->L)*df->atomic_size)); 2227fbf63aeSDave May } else { 22352849c42SDave May /* reallocate pointer list, add +1 in case new_L = 0 */ 2249566063dSJacob Faibussowitsch PetscCall(PetscRealloc(df->atomic_size * (new_L+1), &df->data)); 22552849c42SDave May } 22652849c42SDave May df->L = new_L; 2277fbf63aeSDave May PetscFunctionReturn(0); 22852849c42SDave May } 22952849c42SDave May 23077048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldZeroBlock(DMSwarmDataField df,const PetscInt start,const PetscInt end) 23152849c42SDave May { 232521f74f9SMatthew G. Knepley PetscFunctionBegin; 23363a3b9bcSJacob Faibussowitsch PetscCheck(start <= end,PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if start(%" PetscInt_FMT ") > end(%" PetscInt_FMT ")",start,end); 23463a3b9bcSJacob Faibussowitsch PetscCheck(start >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if start(%" PetscInt_FMT ") < 0",start); 23563a3b9bcSJacob Faibussowitsch PetscCheck(end <= df->L,PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if end(%" PetscInt_FMT ") >= array size(%" PetscInt_FMT ")",end,df->L); 2369566063dSJacob Faibussowitsch PetscCall(PetscMemzero((((char*)df->data)+start*df->atomic_size), (end-start)*df->atomic_size)); 2377fbf63aeSDave May PetscFunctionReturn(0); 23852849c42SDave May } 23952849c42SDave May 24052849c42SDave May /* 24152849c42SDave May A negative buffer value will simply be ignored and the old buffer value will be used. 24252849c42SDave May */ 24377048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketSetSizes(DMSwarmDataBucket db,const PetscInt L,const PetscInt buffer) 24452849c42SDave May { 2455c18a9d6SDave May PetscInt current_allocated,new_used,new_unused,new_buffer,new_allocated,f; 2460cb17e37SDave May PetscBool any_active_fields; 24752849c42SDave May 248521f74f9SMatthew G. Knepley PetscFunctionBegin; 24908401ef6SPierre Jolivet PetscCheck(db->finalised != PETSC_FALSE,PETSC_COMM_SELF,PETSC_ERR_USER,"You must call DMSwarmDataBucketFinalize() before DMSwarmDataBucketSetSizes()"); 2509566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketQueryForActiveFields(db,&any_active_fields)); 25128b400f6SJacob Faibussowitsch PetscCheck(!any_active_fields,PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot safely re-size as at least one DMSwarmDataField is currently being accessed"); 2520cb17e37SDave May 25352849c42SDave May current_allocated = db->allocated; 25452849c42SDave May new_used = L; 25552849c42SDave May new_unused = current_allocated - new_used; 25652849c42SDave May new_buffer = db->buffer; 25752849c42SDave May if (buffer >= 0) { /* update the buffer value */ 25852849c42SDave May new_buffer = buffer; 25952849c42SDave May } 26052849c42SDave May new_allocated = new_used + new_buffer; 26152849c42SDave May /* action */ 26252849c42SDave May if (new_allocated > current_allocated) { 26352849c42SDave May /* increase size to new_used + new_buffer */ 26452849c42SDave May for (f=0; f<db->nfields; f++) { 2659566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldSetSize(db->field[f], new_allocated)); 26652849c42SDave May } 26752849c42SDave May db->L = new_used; 26852849c42SDave May db->buffer = new_buffer; 26952849c42SDave May db->allocated = new_used + new_buffer; 270521f74f9SMatthew G. Knepley } else { 27152849c42SDave May if (new_unused > 2 * new_buffer) { 27252849c42SDave May /* shrink array to new_used + new_buffer */ 273521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 2749566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldSetSize(db->field[f], new_allocated)); 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 db->L = new_used; 28152849c42SDave May db->buffer = new_buffer; 28252849c42SDave May } 28352849c42SDave May } 28452849c42SDave May /* zero all entries from db->L to db->allocated */ 285521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 28677048351SPatrick Sanan DMSwarmDataField field = db->field[f]; 2879566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldZeroBlock(field, db->L,db->allocated)); 28852849c42SDave May } 2897fbf63aeSDave May PetscFunctionReturn(0); 29052849c42SDave May } 29152849c42SDave May 29277048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketSetInitialSizes(DMSwarmDataBucket db,const PetscInt L,const PetscInt buffer) 29352849c42SDave May { 2945c18a9d6SDave May PetscInt f; 295dbe06d34SDave May 296521f74f9SMatthew G. Knepley PetscFunctionBegin; 2979566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketSetSizes(db,L,buffer)); 298521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 29977048351SPatrick Sanan DMSwarmDataField field = db->field[f]; 3009566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldZeroBlock(field,0,db->allocated)); 30152849c42SDave May } 3027fbf63aeSDave May PetscFunctionReturn(0); 30352849c42SDave May } 30452849c42SDave May 30577048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketGetSizes(DMSwarmDataBucket db,PetscInt *L,PetscInt *buffer,PetscInt *allocated) 30652849c42SDave May { 307521f74f9SMatthew G. Knepley PetscFunctionBegin; 30852849c42SDave May if (L) {*L = db->L;} 30952849c42SDave May if (buffer) {*buffer = db->buffer;} 31052849c42SDave May if (allocated) {*allocated = db->allocated;} 3117fbf63aeSDave May PetscFunctionReturn(0); 31252849c42SDave May } 31352849c42SDave May 31477048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketGetGlobalSizes(MPI_Comm comm,DMSwarmDataBucket db,PetscInt *L,PetscInt *buffer,PetscInt *allocated) 31552849c42SDave May { 316521f74f9SMatthew G. Knepley PetscFunctionBegin; 3179566063dSJacob Faibussowitsch if (L) PetscCallMPI(MPI_Allreduce(&db->L,L,1,MPIU_INT,MPI_SUM,comm)); 3189566063dSJacob Faibussowitsch if (buffer) PetscCallMPI(MPI_Allreduce(&db->buffer,buffer,1,MPIU_INT,MPI_SUM,comm)); 3199566063dSJacob Faibussowitsch if (allocated) PetscCallMPI(MPI_Allreduce(&db->allocated,allocated,1,MPIU_INT,MPI_SUM,comm)); 3207fbf63aeSDave May PetscFunctionReturn(0); 32152849c42SDave May } 32252849c42SDave May 32377048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketGetDMSwarmDataFields(DMSwarmDataBucket db,PetscInt *L,DMSwarmDataField *fields[]) 32452849c42SDave May { 325521f74f9SMatthew G. Knepley PetscFunctionBegin; 32652849c42SDave May if (L) {*L = db->nfields;} 32752849c42SDave May if (fields) {*fields = db->field;} 3287fbf63aeSDave May PetscFunctionReturn(0); 32952849c42SDave May } 33052849c42SDave May 33177048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldGetAccess(const DMSwarmDataField gfield) 33252849c42SDave May { 333521f74f9SMatthew G. Knepley PetscFunctionBegin; 33428b400f6SJacob Faibussowitsch PetscCheck(!gfield->active,PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is already active. You must call DMSwarmDataFieldRestoreAccess()",gfield->name); 33552849c42SDave May gfield->active = PETSC_TRUE; 3367fbf63aeSDave May PetscFunctionReturn(0); 33752849c42SDave May } 33852849c42SDave May 33977048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldAccessPoint(const DMSwarmDataField gfield,const PetscInt pid,void **ctx_p) 34052849c42SDave May { 341521f74f9SMatthew G. Knepley PetscFunctionBegin; 34272505a4dSJed Brown *ctx_p = NULL; 343497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD) 34452849c42SDave May /* debug mode */ 34584bcda08SDave May /* check point is valid */ 34608401ef6SPierre Jolivet PetscCheck(pid >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0"); 34763a3b9bcSJacob Faibussowitsch PetscCheck(pid < gfield->L,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %" PetscInt_FMT,gfield->L); 34808401ef6SPierre Jolivet PetscCheck(gfield->active != PETSC_FALSE,PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is not active. You must call DMSwarmDataFieldGetAccess() before point data can be retrivied",gfield->name); 34952849c42SDave May #endif 35077048351SPatrick Sanan *ctx_p = DMSWARM_DATAFIELD_point_access(gfield->data,pid,gfield->atomic_size); 3517fbf63aeSDave May PetscFunctionReturn(0); 35252849c42SDave May } 35352849c42SDave May 35477048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldAccessPointOffset(const DMSwarmDataField gfield,const size_t offset,const PetscInt pid,void **ctx_p) 35552849c42SDave May { 356521f74f9SMatthew G. Knepley PetscFunctionBegin; 357497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD) 35852849c42SDave May /* debug mode */ 35984bcda08SDave May /* check point is valid */ 36008401ef6SPierre Jolivet /* PetscCheck(offset >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be >= 0");*/ 361521f74f9SMatthew G. Knepley /* Note compiler realizes this can never happen with an unsigned PetscInt */ 36208401ef6SPierre Jolivet PetscCheck(offset < gfield->atomic_size,PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be < %zu",gfield->atomic_size); 36384bcda08SDave May /* check point is valid */ 36408401ef6SPierre Jolivet PetscCheck(pid >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0"); 36563a3b9bcSJacob Faibussowitsch PetscCheck(pid < gfield->L,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %" PetscInt_FMT,gfield->L); 36608401ef6SPierre Jolivet PetscCheck(gfield->active != PETSC_FALSE,PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is not active. You must call DMSwarmDataFieldGetAccess() before point data can be retrivied",gfield->name); 36752849c42SDave May #endif 36877048351SPatrick Sanan *ctx_p = DMSWARM_DATAFIELD_point_access_offset(gfield->data,pid,gfield->atomic_size,offset); 3697fbf63aeSDave May PetscFunctionReturn(0); 37052849c42SDave May } 37152849c42SDave May 37277048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldRestoreAccess(DMSwarmDataField gfield) 37352849c42SDave May { 374521f74f9SMatthew G. Knepley PetscFunctionBegin; 37508401ef6SPierre Jolivet PetscCheck(gfield->active != PETSC_FALSE,PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is not active. You must call DMSwarmDataFieldGetAccess()", gfield->name); 37652849c42SDave May gfield->active = PETSC_FALSE; 3777fbf63aeSDave May PetscFunctionReturn(0); 37852849c42SDave May } 37952849c42SDave May 38077048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldVerifyAccess(const DMSwarmDataField gfield,const size_t size) 38152849c42SDave May { 382521f74f9SMatthew G. Knepley PetscFunctionBegin; 383497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD) 38408401ef6SPierre Jolivet PetscCheck(gfield->atomic_size == size,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); 38552849c42SDave May #endif 3867fbf63aeSDave May PetscFunctionReturn(0); 38752849c42SDave May } 38852849c42SDave May 38977048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldGetAtomicSize(const DMSwarmDataField gfield,size_t *size) 39052849c42SDave May { 391521f74f9SMatthew G. Knepley PetscFunctionBegin; 39252849c42SDave May if (size) {*size = gfield->atomic_size;} 3937fbf63aeSDave May PetscFunctionReturn(0); 39452849c42SDave May } 39552849c42SDave May 39677048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldGetEntries(const DMSwarmDataField gfield,void **data) 39752849c42SDave May { 398521f74f9SMatthew G. Knepley PetscFunctionBegin; 399521f74f9SMatthew G. Knepley if (data) {*data = gfield->data;} 4007fbf63aeSDave May PetscFunctionReturn(0); 40152849c42SDave May } 40252849c42SDave May 40377048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldRestoreEntries(const DMSwarmDataField gfield,void **data) 40452849c42SDave May { 405521f74f9SMatthew G. Knepley PetscFunctionBegin; 406521f74f9SMatthew G. Knepley if (data) {*data = NULL;} 4077fbf63aeSDave May PetscFunctionReturn(0); 40852849c42SDave May } 40952849c42SDave May 41052849c42SDave May /* y = x */ 4115627991aSBarry Smith PetscErrorCode DMSwarmDataBucketCopyPoint(const DMSwarmDataBucket xb,const PetscInt pid_x,const DMSwarmDataBucket yb,const PetscInt pid_y) 41252849c42SDave May { 4135c18a9d6SDave May PetscInt f; 414dbe06d34SDave May 415521f74f9SMatthew G. Knepley PetscFunctionBegin; 416521f74f9SMatthew G. Knepley for (f = 0; f < xb->nfields; ++f) { 41752849c42SDave May void *dest; 41852849c42SDave May void *src; 41952849c42SDave May 4209566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldGetAccess(xb->field[f])); 4219566063dSJacob Faibussowitsch if (xb != yb) PetscCall(DMSwarmDataFieldGetAccess( yb->field[f])); 4229566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldAccessPoint(xb->field[f],pid_x, &src)); 4239566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldAccessPoint(yb->field[f],pid_y, &dest)); 4249566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(dest, src, xb->field[f]->atomic_size)); 4259566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldRestoreAccess(xb->field[f])); 4269566063dSJacob Faibussowitsch if (xb != yb) PetscCall(DMSwarmDataFieldRestoreAccess(yb->field[f])); 42752849c42SDave May } 4287fbf63aeSDave May PetscFunctionReturn(0); 42952849c42SDave May } 43052849c42SDave May 43177048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketCreateFromSubset(DMSwarmDataBucket DBIn,const PetscInt N,const PetscInt list[],DMSwarmDataBucket *DB) 43252849c42SDave May { 4335c18a9d6SDave May PetscInt nfields; 43477048351SPatrick Sanan DMSwarmDataField *fields; 4355c18a9d6SDave May PetscInt f,L,buffer,allocated,p; 43652849c42SDave May 437521f74f9SMatthew G. Knepley PetscFunctionBegin; 4389566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketCreate(DB)); 43952849c42SDave May /* copy contents of DBIn */ 4409566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetDMSwarmDataFields(DBIn,&nfields,&fields)); 4419566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(DBIn,&L,&buffer,&allocated)); 442521f74f9SMatthew G. Knepley for (f = 0; f < nfields; ++f) { 4439566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketRegisterField(*DB,"DMSwarmDataBucketCreateFromSubset",fields[f]->name,fields[f]->atomic_size,NULL)); 44452849c42SDave May } 4459566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketFinalize(*DB)); 4469566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketSetSizes(*DB,L,buffer)); 44752849c42SDave May /* now copy the desired guys from DBIn => DB */ 448521f74f9SMatthew G. Knepley for (p = 0; p < N; ++p) { 4499566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketCopyPoint(DBIn,list[p], *DB,list[p])); 45052849c42SDave May } 4517fbf63aeSDave May PetscFunctionReturn(0); 45252849c42SDave May } 45352849c42SDave May 454521f74f9SMatthew G. Knepley /* insert into an exisitng location */ 45577048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldInsertPoint(const DMSwarmDataField field,const PetscInt index,const void *ctx) 45652849c42SDave May { 457521f74f9SMatthew G. Knepley PetscFunctionBegin; 458497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD) 45984bcda08SDave May /* check point is valid */ 46008401ef6SPierre Jolivet PetscCheck(index >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0"); 46163a3b9bcSJacob Faibussowitsch PetscCheck(index < field->L,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %" PetscInt_FMT,field->L); 46252849c42SDave May #endif 4639566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(DMSWARM_DATAFIELD_point_access(field->data,index,field->atomic_size), ctx, field->atomic_size)); 4647fbf63aeSDave May PetscFunctionReturn(0); 46552849c42SDave May } 46652849c42SDave May 467521f74f9SMatthew G. Knepley /* remove data at index - replace with last point */ 46877048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketRemovePointAtIndex(const DMSwarmDataBucket db,const PetscInt index) 46952849c42SDave May { 4705c18a9d6SDave May PetscInt f; 4710cb17e37SDave May PetscBool any_active_fields; 47252849c42SDave May 473521f74f9SMatthew G. Knepley PetscFunctionBegin; 474497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD) 47584bcda08SDave May /* check point is valid */ 47608401ef6SPierre Jolivet PetscCheck(index >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0"); 47763a3b9bcSJacob Faibussowitsch PetscCheck(index < db->allocated,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %" PetscInt_FMT,db->L+db->buffer); 47852849c42SDave May #endif 4799566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketQueryForActiveFields(db,&any_active_fields)); 48028b400f6SJacob Faibussowitsch PetscCheck(!any_active_fields,PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot safely remove point as at least one DMSwarmDataField is currently being accessed"); 481a233d522SDave May if (index >= db->L) { /* this point is not in the list - no need to error, but I will anyway */ 48263a3b9bcSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"You should not be trying to remove point at index=%" PetscInt_FMT " since it's < db->L = %" PetscInt_FMT, index, db->L); 48352849c42SDave May } 484a233d522SDave May if (index != db->L-1) { /* not last point in list */ 485521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 48677048351SPatrick Sanan DMSwarmDataField field = db->field[f]; 48752849c42SDave May 48852849c42SDave May /* copy then remove */ 4899566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldCopyPoint(db->L-1, field, index, field)); 49077048351SPatrick Sanan /* DMSwarmDataFieldZeroPoint(field,index); */ 49152849c42SDave May } 49252849c42SDave May } 49352849c42SDave May /* decrement size */ 49452849c42SDave May /* this will zero out an crap at the end of the list */ 4959566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketRemovePoint(db)); 4967fbf63aeSDave May PetscFunctionReturn(0); 49752849c42SDave May } 49852849c42SDave May 49952849c42SDave May /* copy x into y */ 5005627991aSBarry Smith PetscErrorCode DMSwarmDataFieldCopyPoint(const PetscInt pid_x,const DMSwarmDataField field_x,const PetscInt pid_y,const DMSwarmDataField field_y) 50152849c42SDave May { 502521f74f9SMatthew G. Knepley PetscFunctionBegin; 503497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD) 50484bcda08SDave May /* check point is valid */ 50508401ef6SPierre Jolivet PetscCheck(pid_x >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"(IN) index must be >= 0"); 50663a3b9bcSJacob Faibussowitsch PetscCheck(pid_x < field_x->L,PETSC_COMM_SELF,PETSC_ERR_USER,"(IN) index must be < %" PetscInt_FMT,field_x->L); 50708401ef6SPierre Jolivet PetscCheck(pid_y >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"(OUT) index must be >= 0"); 50863a3b9bcSJacob Faibussowitsch PetscCheck(pid_y < field_y->L,PETSC_COMM_SELF,PETSC_ERR_USER,"(OUT) index must be < %" PetscInt_FMT,field_y->L); 50908401ef6SPierre Jolivet PetscCheck(field_y->atomic_size == field_x->atomic_size,PETSC_COMM_SELF,PETSC_ERR_USER,"atomic size must match"); 51052849c42SDave May #endif 5119566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(DMSWARM_DATAFIELD_point_access(field_y->data,pid_y,field_y->atomic_size),DMSWARM_DATAFIELD_point_access(field_x->data,pid_x,field_x->atomic_size),field_y->atomic_size)); 5127fbf63aeSDave May PetscFunctionReturn(0); 51352849c42SDave May } 51452849c42SDave May 515521f74f9SMatthew G. Knepley /* zero only the datafield at this point */ 51677048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldZeroPoint(const DMSwarmDataField field,const PetscInt index) 51752849c42SDave May { 518521f74f9SMatthew G. Knepley PetscFunctionBegin; 519497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD) 52084bcda08SDave May /* check point is valid */ 52108401ef6SPierre Jolivet PetscCheck(index >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0"); 52263a3b9bcSJacob Faibussowitsch PetscCheck(index < field->L,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %" PetscInt_FMT,field->L); 52352849c42SDave May #endif 5249566063dSJacob Faibussowitsch PetscCall(PetscMemzero(DMSWARM_DATAFIELD_point_access(field->data,index,field->atomic_size), field->atomic_size)); 5257fbf63aeSDave May PetscFunctionReturn(0); 52652849c42SDave May } 52752849c42SDave May 528521f74f9SMatthew G. Knepley /* zero ALL data for this point */ 52977048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketZeroPoint(const DMSwarmDataBucket db,const PetscInt index) 53052849c42SDave May { 5315c18a9d6SDave May PetscInt f; 53252849c42SDave May 533521f74f9SMatthew G. Knepley PetscFunctionBegin; 53484bcda08SDave May /* check point is valid */ 53508401ef6SPierre Jolivet PetscCheck(index >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0"); 53663a3b9bcSJacob Faibussowitsch PetscCheck(index < db->allocated,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %" PetscInt_FMT,db->allocated); 537521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 53877048351SPatrick Sanan DMSwarmDataField field = db->field[f]; 5399566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldZeroPoint(field,index)); 54052849c42SDave May } 5417fbf63aeSDave May PetscFunctionReturn(0); 54252849c42SDave May } 54352849c42SDave May 54452849c42SDave May /* increment */ 54577048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketAddPoint(DMSwarmDataBucket db) 54652849c42SDave May { 547521f74f9SMatthew G. Knepley PetscFunctionBegin; 5489566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketSetSizes(db,db->L+1,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT)); 5497fbf63aeSDave May PetscFunctionReturn(0); 55052849c42SDave May } 55152849c42SDave May 5527fbf63aeSDave May /* decrement */ 55377048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketRemovePoint(DMSwarmDataBucket db) 5547fbf63aeSDave May { 555521f74f9SMatthew G. Knepley PetscFunctionBegin; 5569566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketSetSizes(db,db->L-1,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT)); 5577fbf63aeSDave May PetscFunctionReturn(0); 5587fbf63aeSDave May } 5597fbf63aeSDave May 5605627991aSBarry Smith /* Should be redone to user PetscViewer */ 56177048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketView_stdout(MPI_Comm comm,DMSwarmDataBucket db) 562ab2be9a3SDave May { 563ab2be9a3SDave May PetscInt f; 564ab2be9a3SDave May double memory_usage_total,memory_usage_total_local = 0.0; 565ab2be9a3SDave May 566ab2be9a3SDave May PetscFunctionBegin; 5679566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"DMSwarmDataBucketView: \n")); 56863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm," L = %" PetscInt_FMT " \n", db->L)); 56963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm," buffer = %" PetscInt_FMT " \n", db->buffer)); 57063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm," allocated = %" PetscInt_FMT " \n", db->allocated)); 57163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm," nfields registered = %" PetscInt_FMT " \n", db->nfields)); 572ab2be9a3SDave May 573ab2be9a3SDave May for (f = 0; f < db->nfields; ++f) { 574ab2be9a3SDave May double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6; 575ab2be9a3SDave May memory_usage_total_local += memory_usage_f; 576ab2be9a3SDave May } 5779566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&memory_usage_total_local,&memory_usage_total,1,MPI_DOUBLE,MPI_SUM,comm)); 578ab2be9a3SDave May 579ab2be9a3SDave May for (f = 0; f < db->nfields; ++f) { 580ab2be9a3SDave May double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6; 58163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm," [%3" PetscInt_FMT "] %15s : Mem. usage = %1.2e (MB) [rank0]\n", f, db->field[f]->name, memory_usage_f)); 58263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm," blocksize = %" PetscInt_FMT " \n", db->field[f]->bs)); 583ab2be9a3SDave May if (db->field[f]->bs != 1) { 58463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm," atomic size = %zu [full block, bs=%" PetscInt_FMT "]\n", db->field[f]->atomic_size,db->field[f]->bs)); 58563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm," atomic size/item = %zu \n", (size_t)(db->field[f]->atomic_size/db->field[f]->bs))); 586ab2be9a3SDave May } else { 5879566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm," atomic size = %zu \n", db->field[f]->atomic_size)); 588ab2be9a3SDave May } 589ab2be9a3SDave May } 5909566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm," Total mem. usage = %1.2e (MB) (collective)\n", memory_usage_total)); 591ab2be9a3SDave May PetscFunctionReturn(0); 592ab2be9a3SDave May } 593ab2be9a3SDave May 5945627991aSBarry Smith PetscErrorCode DMSwarmDataBucketView_Seq(MPI_Comm comm,DMSwarmDataBucket db,const char filename[],DMSwarmDataBucketViewType type) 59552849c42SDave May { 596521f74f9SMatthew G. Knepley PetscFunctionBegin; 59752849c42SDave May switch (type) { 59852849c42SDave May case DATABUCKET_VIEW_STDOUT: 5999566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketView_stdout(PETSC_COMM_SELF,db)); 60052849c42SDave May break; 60152849c42SDave May case DATABUCKET_VIEW_ASCII: 602071ce369SDave May SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for ascii output"); 60352849c42SDave May case DATABUCKET_VIEW_BINARY: 604071ce369SDave May SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for binary output"); 60552849c42SDave May case DATABUCKET_VIEW_HDF5: 606071ce369SDave May SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for HDF5 output"); 607521f74f9SMatthew G. Knepley default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer method requested"); 60852849c42SDave May } 6097fbf63aeSDave May PetscFunctionReturn(0); 61052849c42SDave May } 61152849c42SDave May 61277048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketView_MPI(MPI_Comm comm,DMSwarmDataBucket db,const char filename[],DMSwarmDataBucketViewType type) 61352849c42SDave May { 614521f74f9SMatthew G. Knepley PetscFunctionBegin; 61552849c42SDave May switch (type) { 61652849c42SDave May case DATABUCKET_VIEW_STDOUT: 6179566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketView_stdout(comm,db)); 61852849c42SDave May break; 61952849c42SDave May case DATABUCKET_VIEW_ASCII: 620071ce369SDave May SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for ascii output"); 62152849c42SDave May case DATABUCKET_VIEW_BINARY: 622071ce369SDave May SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for binary output"); 62352849c42SDave May case DATABUCKET_VIEW_HDF5: 624071ce369SDave May SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for HDF5 output"); 625521f74f9SMatthew G. Knepley default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer method requested"); 62652849c42SDave May } 6277fbf63aeSDave May PetscFunctionReturn(0); 62852849c42SDave May } 62952849c42SDave May 63077048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketView(MPI_Comm comm,DMSwarmDataBucket db,const char filename[],DMSwarmDataBucketViewType type) 63152849c42SDave May { 632d7d19db6SBarry Smith PetscMPIInt size; 63352849c42SDave May 634521f74f9SMatthew G. Knepley PetscFunctionBegin; 6359566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 636d7d19db6SBarry Smith if (size == 1) { 6379566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketView_Seq(comm,db,filename,type)); 63852849c42SDave May } else { 6399566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketView_MPI(comm,db,filename,type)); 64052849c42SDave May } 6417fbf63aeSDave May PetscFunctionReturn(0); 64252849c42SDave May } 64352849c42SDave May 64477048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketDuplicateFields(DMSwarmDataBucket dbA,DMSwarmDataBucket *dbB) 64552849c42SDave May { 64677048351SPatrick Sanan DMSwarmDataBucket db2; 6475c18a9d6SDave May PetscInt f; 64852849c42SDave May 649521f74f9SMatthew G. Knepley PetscFunctionBegin; 6509566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketCreate(&db2)); 65152849c42SDave May /* copy contents from dbA into db2 */ 652521f74f9SMatthew G. Knepley for (f = 0; f < dbA->nfields; ++f) { 65377048351SPatrick Sanan DMSwarmDataField field; 65452849c42SDave May size_t atomic_size; 65552849c42SDave May char *name; 65652849c42SDave May 65752849c42SDave May field = dbA->field[f]; 65852849c42SDave May atomic_size = field->atomic_size; 65952849c42SDave May name = field->name; 6609566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketRegisterField(db2,"DMSwarmDataBucketDuplicateFields",name,atomic_size,NULL)); 66152849c42SDave May } 6629566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketFinalize(db2)); 6639566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketSetInitialSizes(db2,0,1000)); 66452849c42SDave May *dbB = db2; 6657fbf63aeSDave May PetscFunctionReturn(0); 66652849c42SDave May } 66752849c42SDave May 66852849c42SDave May /* 66952849c42SDave May Insert points from db2 into db1 67052849c42SDave May db1 <<== db2 67152849c42SDave May */ 67277048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketInsertValues(DMSwarmDataBucket db1,DMSwarmDataBucket db2) 67352849c42SDave May { 6745c18a9d6SDave May PetscInt n_mp_points1,n_mp_points2; 6755c18a9d6SDave May PetscInt n_mp_points1_new,p; 67652849c42SDave May 677521f74f9SMatthew G. Knepley PetscFunctionBegin; 6789566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(db1,&n_mp_points1,NULL,NULL)); 6799566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(db2,&n_mp_points2,NULL,NULL)); 68052849c42SDave May n_mp_points1_new = n_mp_points1 + n_mp_points2; 6819566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketSetSizes(db1,n_mp_points1_new,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT)); 682521f74f9SMatthew G. Knepley for (p = 0; p < n_mp_points2; ++p) { 683521f74f9SMatthew G. Knepley /* db1 <<== db2 */ 6849566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketCopyPoint(db2,p, db1,(n_mp_points1 + p))); 68552849c42SDave May } 6867fbf63aeSDave May PetscFunctionReturn(0); 68752849c42SDave May } 68852849c42SDave May 68952849c42SDave May /* helpers for parallel send/recv */ 69077048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketCreatePackedArray(DMSwarmDataBucket db,size_t *bytes,void **buf) 69152849c42SDave May { 6925c18a9d6SDave May PetscInt f; 69352849c42SDave May size_t sizeof_marker_contents; 69452849c42SDave May void *buffer; 69552849c42SDave May 696521f74f9SMatthew G. Knepley PetscFunctionBegin; 69752849c42SDave May sizeof_marker_contents = 0; 698521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 69977048351SPatrick Sanan DMSwarmDataField df = db->field[f]; 70052849c42SDave May sizeof_marker_contents += df->atomic_size; 70152849c42SDave May } 7029566063dSJacob Faibussowitsch PetscCall(PetscMalloc(sizeof_marker_contents, &buffer)); 7039566063dSJacob Faibussowitsch PetscCall(PetscMemzero(buffer, sizeof_marker_contents)); 70452849c42SDave May if (bytes) {*bytes = sizeof_marker_contents;} 70552849c42SDave May if (buf) {*buf = buffer;} 7067fbf63aeSDave May PetscFunctionReturn(0); 70752849c42SDave May } 70852849c42SDave May 70977048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketDestroyPackedArray(DMSwarmDataBucket db,void **buf) 71052849c42SDave May { 711521f74f9SMatthew G. Knepley PetscFunctionBegin; 71252849c42SDave May if (buf) { 7139566063dSJacob Faibussowitsch PetscCall(PetscFree(*buf)); 71452849c42SDave May *buf = NULL; 71552849c42SDave May } 7167fbf63aeSDave May PetscFunctionReturn(0); 71752849c42SDave May } 71852849c42SDave May 71977048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketFillPackedArray(DMSwarmDataBucket db,const PetscInt index,void *buf) 72052849c42SDave May { 7215c18a9d6SDave May PetscInt f; 72252849c42SDave May void *data, *data_p; 72352849c42SDave May size_t asize, offset; 72452849c42SDave May 725521f74f9SMatthew G. Knepley PetscFunctionBegin; 72652849c42SDave May offset = 0; 727521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 72877048351SPatrick Sanan DMSwarmDataField df = db->field[f]; 72952849c42SDave May 73052849c42SDave May asize = df->atomic_size; 73152849c42SDave May data = (void*)( df->data); 73252849c42SDave May data_p = (void*)( (char*)data + index*asize); 7339566063dSJacob Faibussowitsch PetscCall(PetscMemcpy((void*)((char*)buf + offset), data_p, asize)); 73452849c42SDave May offset = offset + asize; 73552849c42SDave May } 7367fbf63aeSDave May PetscFunctionReturn(0); 73752849c42SDave May } 73852849c42SDave May 73977048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketInsertPackedArray(DMSwarmDataBucket db,const PetscInt idx,void *data) 74052849c42SDave May { 7415c18a9d6SDave May PetscInt f; 74252849c42SDave May void *data_p; 74352849c42SDave May size_t offset; 74452849c42SDave May 745521f74f9SMatthew G. Knepley PetscFunctionBegin; 74652849c42SDave May offset = 0; 747521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 74877048351SPatrick Sanan DMSwarmDataField df = db->field[f]; 74952849c42SDave May 75052849c42SDave May data_p = (void*)( (char*)data + offset); 7519566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldInsertPoint(df, idx, (void*)data_p)); 75252849c42SDave May offset = offset + df->atomic_size; 75352849c42SDave May } 7547fbf63aeSDave May PetscFunctionReturn(0); 75552849c42SDave May } 756