1279f676cSBarry Smith #include "../src/dm/impls/swarm/data_bucket.h" 252849c42SDave May 352849c42SDave May /* string helpers */ 4d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldStringInList(const char name[], const PetscInt N, const DMSwarmDataField gfield[], PetscBool *val) 5d71ae5a4SJacob Faibussowitsch { 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; 153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1652849c42SDave May } 1752849c42SDave May } 183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1952849c42SDave May } 2052849c42SDave May 21d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldStringFindInList(const char name[], const PetscInt N, const DMSwarmDataField gfield[], PetscInt *index) 22d71ae5a4SJacob Faibussowitsch { 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; 323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3352849c42SDave May } 3452849c42SDave May } 353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3652849c42SDave May } 3752849c42SDave May 38d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldCreate(const char registration_function[], const char name[], const size_t size, const PetscInt L, DMSwarmDataField *DF) 39d71ae5a4SJacob Faibussowitsch { 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; 533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5452849c42SDave May } 5552849c42SDave May 56d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldDestroy(DMSwarmDataField *DF) 57d71ae5a4SJacob Faibussowitsch { 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; 663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6752849c42SDave May } 6852849c42SDave May 6952849c42SDave May /* data bucket */ 70d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketCreate(DMSwarmDataBucket *DB) 71d71ae5a4SJacob Faibussowitsch { 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; 853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8652849c42SDave May } 8752849c42SDave May 88d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketDestroy(DMSwarmDataBucket *DB) 89d71ae5a4SJacob Faibussowitsch { 9077048351SPatrick Sanan DMSwarmDataBucket db = *DB; 915c18a9d6SDave May PetscInt f; 9252849c42SDave May 93521f74f9SMatthew G. Knepley PetscFunctionBegin; 9452849c42SDave May /* release fields */ 9548a46eb9SPierre Jolivet for (f = 0; f < db->nfields; ++f) PetscCall(DMSwarmDataFieldDestroy(&db->field[f])); 9652849c42SDave May /* this will catch the initially allocated objects in the event that no fields are registered */ 9748a46eb9SPierre Jolivet if (db->field != NULL) PetscCall(PetscFree(db->field)); 989566063dSJacob Faibussowitsch PetscCall(PetscFree(db)); 9952849c42SDave May *DB = NULL; 1003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10152849c42SDave May } 10252849c42SDave May 103d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketQueryForActiveFields(DMSwarmDataBucket db, PetscBool *any_active_fields) 104d71ae5a4SJacob Faibussowitsch { 1050cb17e37SDave May PetscInt f; 106521f74f9SMatthew G. Knepley 107521f74f9SMatthew G. Knepley PetscFunctionBegin; 1080cb17e37SDave May *any_active_fields = PETSC_FALSE; 109521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 1100cb17e37SDave May if (db->field[f]->active) { 1110cb17e37SDave May *any_active_fields = PETSC_TRUE; 1123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1130cb17e37SDave May } 1140cb17e37SDave May } 1153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1160cb17e37SDave May } 1170cb17e37SDave May 118d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketRegisterField(DMSwarmDataBucket db, const char registration_function[], const char field_name[], size_t atomic_size, DMSwarmDataField *_gfield) 119d71ae5a4SJacob Faibussowitsch { 12052849c42SDave May PetscBool val; 12177048351SPatrick Sanan DMSwarmDataField fp; 12252849c42SDave May 123521f74f9SMatthew G. Knepley PetscFunctionBegin; 12452849c42SDave May /* check we haven't finalised the registration of fields */ 12552849c42SDave May /* 12652849c42SDave May if (db->finalised==PETSC_TRUE) { 12777048351SPatrick Sanan printf("ERROR: DMSwarmDataBucketFinalize() has been called. Cannot register more fields\n"); 12852849c42SDave May ERROR(); 12952849c42SDave May } 13052849c42SDave May */ 13152849c42SDave May /* check for repeated name */ 1329566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldStringInList(field_name, db->nfields, (const DMSwarmDataField *)db->field, &val)); 13308401ef6SPierre Jolivet PetscCheck(val != PETSC_TRUE, PETSC_COMM_SELF, PETSC_ERR_USER, "Field %s already exists. Cannot add same field twice", field_name); 13452849c42SDave May /* create new space for data */ 1359566063dSJacob Faibussowitsch PetscCall(PetscRealloc(sizeof(DMSwarmDataField) * (db->nfields + 1), &db->field)); 13652849c42SDave May /* add field */ 1379566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldCreate(registration_function, field_name, atomic_size, db->allocated, &fp)); 13852849c42SDave May db->field[db->nfields] = fp; 13952849c42SDave May db->nfields++; 140ad540459SPierre Jolivet if (_gfield != NULL) *_gfield = fp; 1413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14252849c42SDave May } 14352849c42SDave May 14452849c42SDave May /* 14577048351SPatrick Sanan #define DMSwarmDataBucketRegisterField(db,name,size,k) {\ 14652849c42SDave May char *location;\ 14752849c42SDave May asprintf(&location,"Registered by %s() at line %d within file %s", __FUNCTION__, __LINE__, __FILE__);\ 14877048351SPatrick Sanan _DMSwarmDataBucketRegisterField( (db), location, (name), (size), (k));\ 149521f74f9SMatthew G. Knepley ierr = PetscFree(location);\ 15052849c42SDave May } 15152849c42SDave May */ 15252849c42SDave May 153d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketGetDMSwarmDataFieldIdByName(DMSwarmDataBucket db, const char name[], PetscInt *idx) 154d71ae5a4SJacob Faibussowitsch { 1552e956fe4SStefano Zampini PetscBool found; 1562e956fe4SStefano Zampini 1572e956fe4SStefano Zampini PetscFunctionBegin; 1582e956fe4SStefano Zampini *idx = -1; 1592e956fe4SStefano Zampini PetscCall(DMSwarmDataFieldStringInList(name, db->nfields, (const DMSwarmDataField *)db->field, &found)); 1602e956fe4SStefano Zampini PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot find DMSwarmDataField with name %s", name); 1612e956fe4SStefano Zampini PetscCall(DMSwarmDataFieldStringFindInList(name, db->nfields, (const DMSwarmDataField *)db->field, idx)); 1623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1632e956fe4SStefano Zampini } 1642e956fe4SStefano Zampini 165d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketGetDMSwarmDataFieldByName(DMSwarmDataBucket db, const char name[], DMSwarmDataField *gfield) 166d71ae5a4SJacob Faibussowitsch { 1675c18a9d6SDave May PetscInt idx; 16852849c42SDave May PetscBool found; 16952849c42SDave May 170521f74f9SMatthew G. Knepley PetscFunctionBegin; 1719566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldStringInList(name, db->nfields, (const DMSwarmDataField *)db->field, &found)); 17228b400f6SJacob Faibussowitsch PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot find DMSwarmDataField with name %s", name); 1739566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldStringFindInList(name, db->nfields, (const DMSwarmDataField *)db->field, &idx)); 17452849c42SDave May *gfield = db->field[idx]; 1753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 17652849c42SDave May } 17752849c42SDave May 178d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketQueryDMSwarmDataFieldByName(DMSwarmDataBucket db, const char name[], PetscBool *found) 179d71ae5a4SJacob Faibussowitsch { 180521f74f9SMatthew G. Knepley PetscFunctionBegin; 18152849c42SDave May *found = PETSC_FALSE; 1829566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldStringInList(name, db->nfields, (const DMSwarmDataField *)db->field, found)); 1833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 18452849c42SDave May } 18552849c42SDave May 186d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketFinalize(DMSwarmDataBucket db) 187d71ae5a4SJacob Faibussowitsch { 188521f74f9SMatthew G. Knepley PetscFunctionBegin; 18952849c42SDave May db->finalised = PETSC_TRUE; 1903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 19152849c42SDave May } 19252849c42SDave May 193d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldGetNumEntries(DMSwarmDataField df, PetscInt *sum) 194d71ae5a4SJacob Faibussowitsch { 195521f74f9SMatthew G. Knepley PetscFunctionBegin; 19652849c42SDave May *sum = df->L; 1973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 19852849c42SDave May } 19952849c42SDave May 200d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldSetBlockSize(DMSwarmDataField df, PetscInt blocksize) 201d71ae5a4SJacob Faibussowitsch { 202521f74f9SMatthew G. Knepley PetscFunctionBegin; 20352c3ed93SDave May df->bs = blocksize; 2043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 20552c3ed93SDave May } 20652c3ed93SDave May 207d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldSetSize(DMSwarmDataField df, const PetscInt new_L) 208d71ae5a4SJacob Faibussowitsch { 209521f74f9SMatthew G. Knepley PetscFunctionBegin; 21008401ef6SPierre Jolivet PetscCheck(new_L >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot set size of DMSwarmDataField to be < 0"); 2113ba16761SJacob Faibussowitsch if (new_L == df->L) PetscFunctionReturn(PETSC_SUCCESS); 21252849c42SDave May if (new_L > df->L) { 2139566063dSJacob Faibussowitsch PetscCall(PetscRealloc(df->atomic_size * (new_L), &df->data)); 21452849c42SDave May /* init new contents */ 21557508eceSPierre Jolivet PetscCall(PetscMemzero(((char *)df->data) + df->L * df->atomic_size, (new_L - df->L) * df->atomic_size)); 2167fbf63aeSDave May } else { 21752849c42SDave May /* reallocate pointer list, add +1 in case new_L = 0 */ 2189566063dSJacob Faibussowitsch PetscCall(PetscRealloc(df->atomic_size * (new_L + 1), &df->data)); 21952849c42SDave May } 22052849c42SDave May df->L = new_L; 2213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 22252849c42SDave May } 22352849c42SDave May 224d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldZeroBlock(DMSwarmDataField df, const PetscInt start, const PetscInt end) 225d71ae5a4SJacob Faibussowitsch { 226521f74f9SMatthew G. Knepley PetscFunctionBegin; 22763a3b9bcSJacob 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); 22863a3b9bcSJacob Faibussowitsch PetscCheck(start >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot zero a block of entries if start(%" PetscInt_FMT ") < 0", start); 22963a3b9bcSJacob 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); 23057508eceSPierre Jolivet PetscCall(PetscMemzero(((char *)df->data) + start * df->atomic_size, (end - start) * df->atomic_size)); 2313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 23252849c42SDave May } 23352849c42SDave May 23452849c42SDave May /* 23552849c42SDave May A negative buffer value will simply be ignored and the old buffer value will be used. 23652849c42SDave May */ 237d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketSetSizes(DMSwarmDataBucket db, const PetscInt L, const PetscInt buffer) 238d71ae5a4SJacob Faibussowitsch { 239ccc1fad5SZach Atkins PetscInt current_allocated, current_used, new_used, new_unused, new_buffer, new_allocated, f, end; 2400cb17e37SDave May PetscBool any_active_fields; 24152849c42SDave May 242521f74f9SMatthew G. Knepley PetscFunctionBegin; 24308401ef6SPierre Jolivet PetscCheck(db->finalised != PETSC_FALSE, PETSC_COMM_SELF, PETSC_ERR_USER, "You must call DMSwarmDataBucketFinalize() before DMSwarmDataBucketSetSizes()"); 2449566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketQueryForActiveFields(db, &any_active_fields)); 24528b400f6SJacob Faibussowitsch PetscCheck(!any_active_fields, PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot safely re-size as at least one DMSwarmDataField is currently being accessed"); 2460cb17e37SDave May 24752849c42SDave May current_allocated = db->allocated; 248ccc1fad5SZach Atkins current_used = PetscMax(db->L, 0); 24952849c42SDave May new_used = L; 25052849c42SDave May new_unused = current_allocated - new_used; 25152849c42SDave May new_buffer = db->buffer; 25252849c42SDave May if (buffer >= 0) { /* update the buffer value */ 25352849c42SDave May new_buffer = buffer; 25452849c42SDave May } 25552849c42SDave May new_allocated = new_used + new_buffer; 25652849c42SDave May /* action */ 25752849c42SDave May if (new_allocated > current_allocated) { 258ccc1fad5SZach Atkins /* increase size to new_used + new_buffer and zero new space */ 259ccc1fad5SZach Atkins for (f = 0; f < db->nfields; f++) { 260ccc1fad5SZach Atkins PetscCall(DMSwarmDataFieldSetSize(db->field[f], new_allocated)); 261ccc1fad5SZach Atkins PetscCall(DMSwarmDataFieldZeroBlock(db->field[f], current_allocated, new_allocated)); 262ccc1fad5SZach Atkins } 26352849c42SDave May db->L = new_used; 26452849c42SDave May db->buffer = new_buffer; 26552849c42SDave May db->allocated = new_used + new_buffer; 266521f74f9SMatthew G. Knepley } else { 26752849c42SDave May if (new_unused > 2 * new_buffer) { 26852849c42SDave May /* shrink array to new_used + new_buffer */ 26948a46eb9SPierre Jolivet for (f = 0; f < db->nfields; ++f) PetscCall(DMSwarmDataFieldSetSize(db->field[f], new_allocated)); 27052849c42SDave May db->L = new_used; 27152849c42SDave May db->buffer = new_buffer; 27252849c42SDave May db->allocated = new_used + new_buffer; 273521f74f9SMatthew G. Knepley } else { 27452849c42SDave May db->L = new_used; 27552849c42SDave May db->buffer = new_buffer; 27652849c42SDave May } 27752849c42SDave May } 278ccc1fad5SZach Atkins /* if we shrunk, zero old entries from new_used to current_used or end of array */ 279ccc1fad5SZach Atkins end = PetscMin(current_used, new_allocated); 280ccc1fad5SZach Atkins if (end > new_used) { 281521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 28277048351SPatrick Sanan DMSwarmDataField field = db->field[f]; 283ccc1fad5SZach Atkins PetscCall(DMSwarmDataFieldZeroBlock(field, new_used, end)); 284ccc1fad5SZach Atkins } 28552849c42SDave May } 2863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 28752849c42SDave May } 28852849c42SDave May 289d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketSetInitialSizes(DMSwarmDataBucket db, const PetscInt L, const PetscInt buffer) 290d71ae5a4SJacob Faibussowitsch { 2915c18a9d6SDave May PetscInt f; 292dbe06d34SDave May 293521f74f9SMatthew G. Knepley PetscFunctionBegin; 2949566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketSetSizes(db, L, buffer)); 295521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 29677048351SPatrick Sanan DMSwarmDataField field = db->field[f]; 2979566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldZeroBlock(field, 0, db->allocated)); 29852849c42SDave May } 2993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 30052849c42SDave May } 30152849c42SDave May 302d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketGetSizes(DMSwarmDataBucket db, PetscInt *L, PetscInt *buffer, PetscInt *allocated) 303d71ae5a4SJacob Faibussowitsch { 304521f74f9SMatthew G. Knepley PetscFunctionBegin; 305ad540459SPierre Jolivet if (L) *L = db->L; 306ad540459SPierre Jolivet if (buffer) *buffer = db->buffer; 307ad540459SPierre Jolivet if (allocated) *allocated = db->allocated; 3083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 30952849c42SDave May } 31052849c42SDave May 311d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketGetGlobalSizes(MPI_Comm comm, DMSwarmDataBucket db, PetscInt *L, PetscInt *buffer, PetscInt *allocated) 312d71ae5a4SJacob Faibussowitsch { 313521f74f9SMatthew G. Knepley PetscFunctionBegin; 314462c564dSBarry Smith if (L) PetscCallMPI(MPIU_Allreduce(&db->L, L, 1, MPIU_INT, MPI_SUM, comm)); 315462c564dSBarry Smith if (buffer) PetscCallMPI(MPIU_Allreduce(&db->buffer, buffer, 1, MPIU_INT, MPI_SUM, comm)); 316462c564dSBarry Smith if (allocated) PetscCallMPI(MPIU_Allreduce(&db->allocated, allocated, 1, MPIU_INT, MPI_SUM, comm)); 3173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 31852849c42SDave May } 31952849c42SDave May 320d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketGetDMSwarmDataFields(DMSwarmDataBucket db, PetscInt *L, DMSwarmDataField *fields[]) 321d71ae5a4SJacob Faibussowitsch { 322521f74f9SMatthew G. Knepley PetscFunctionBegin; 323ad540459SPierre Jolivet if (L) *L = db->nfields; 324ad540459SPierre Jolivet if (fields) *fields = db->field; 3253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 32652849c42SDave May } 32752849c42SDave May 328d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldGetAccess(const DMSwarmDataField gfield) 329d71ae5a4SJacob Faibussowitsch { 330521f74f9SMatthew G. Knepley PetscFunctionBegin; 33128b400f6SJacob Faibussowitsch PetscCheck(!gfield->active, PETSC_COMM_SELF, PETSC_ERR_USER, "Field \"%s\" is already active. You must call DMSwarmDataFieldRestoreAccess()", gfield->name); 33252849c42SDave May gfield->active = PETSC_TRUE; 3333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 33452849c42SDave May } 33552849c42SDave May 336d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldAccessPoint(const DMSwarmDataField gfield, const PetscInt pid, void **ctx_p) 337d71ae5a4SJacob Faibussowitsch { 338521f74f9SMatthew G. Knepley PetscFunctionBegin; 33972505a4dSJed Brown *ctx_p = NULL; 340497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD) 34152849c42SDave May /* debug mode */ 34284bcda08SDave May /* check point is valid */ 34308401ef6SPierre Jolivet PetscCheck(pid >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be >= 0"); 34463a3b9bcSJacob Faibussowitsch PetscCheck(pid < gfield->L, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be < %" PetscInt_FMT, gfield->L); 34508401ef6SPierre 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); 34652849c42SDave May #endif 34777048351SPatrick Sanan *ctx_p = DMSWARM_DATAFIELD_point_access(gfield->data, pid, gfield->atomic_size); 3483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 34952849c42SDave May } 35052849c42SDave May 351d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldAccessPointOffset(const DMSwarmDataField gfield, const size_t offset, const PetscInt pid, void **ctx_p) 352d71ae5a4SJacob Faibussowitsch { 353521f74f9SMatthew G. Knepley PetscFunctionBegin; 354497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD) 35552849c42SDave May /* debug mode */ 35684bcda08SDave May /* check point is valid */ 35708401ef6SPierre Jolivet /* PetscCheck(offset >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be >= 0");*/ 358521f74f9SMatthew G. Knepley /* Note compiler realizes this can never happen with an unsigned PetscInt */ 35908401ef6SPierre Jolivet PetscCheck(offset < gfield->atomic_size, PETSC_COMM_SELF, PETSC_ERR_USER, "offset must be < %zu", gfield->atomic_size); 36084bcda08SDave May /* check point is valid */ 36108401ef6SPierre Jolivet PetscCheck(pid >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be >= 0"); 36263a3b9bcSJacob Faibussowitsch PetscCheck(pid < gfield->L, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be < %" PetscInt_FMT, gfield->L); 36308401ef6SPierre 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); 36452849c42SDave May #endif 36577048351SPatrick Sanan *ctx_p = DMSWARM_DATAFIELD_point_access_offset(gfield->data, pid, gfield->atomic_size, offset); 3663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 36752849c42SDave May } 36852849c42SDave May 369d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldRestoreAccess(DMSwarmDataField gfield) 370d71ae5a4SJacob Faibussowitsch { 371521f74f9SMatthew G. Knepley PetscFunctionBegin; 37208401ef6SPierre Jolivet PetscCheck(gfield->active != PETSC_FALSE, PETSC_COMM_SELF, PETSC_ERR_USER, "Field \"%s\" is not active. You must call DMSwarmDataFieldGetAccess()", gfield->name); 37352849c42SDave May gfield->active = PETSC_FALSE; 3743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 37552849c42SDave May } 37652849c42SDave May 377d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldVerifyAccess(const DMSwarmDataField gfield, const size_t size) 378d71ae5a4SJacob Faibussowitsch { 379521f74f9SMatthew G. Knepley PetscFunctionBegin; 380497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD) 38108401ef6SPierre 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); 38252849c42SDave May #endif 3833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 38452849c42SDave May } 38552849c42SDave May 386d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldGetAtomicSize(const DMSwarmDataField gfield, size_t *size) 387d71ae5a4SJacob Faibussowitsch { 388521f74f9SMatthew G. Knepley PetscFunctionBegin; 389ad540459SPierre Jolivet if (size) *size = gfield->atomic_size; 3903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 39152849c42SDave May } 39252849c42SDave May 393d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldGetEntries(const DMSwarmDataField gfield, void **data) 394d71ae5a4SJacob Faibussowitsch { 395521f74f9SMatthew G. Knepley PetscFunctionBegin; 396ad540459SPierre Jolivet if (data) *data = gfield->data; 3973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 39852849c42SDave May } 39952849c42SDave May 400d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldRestoreEntries(const DMSwarmDataField gfield, void **data) 401d71ae5a4SJacob Faibussowitsch { 402521f74f9SMatthew G. Knepley PetscFunctionBegin; 403ad540459SPierre Jolivet if (data) *data = NULL; 4043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 40552849c42SDave May } 40652849c42SDave May 40752849c42SDave May /* y = x */ 408d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketCopyPoint(const DMSwarmDataBucket xb, const PetscInt pid_x, const DMSwarmDataBucket yb, const PetscInt pid_y) 409d71ae5a4SJacob Faibussowitsch { 4105c18a9d6SDave May PetscInt f; 411dbe06d34SDave May 412521f74f9SMatthew G. Knepley PetscFunctionBegin; 413521f74f9SMatthew G. Knepley for (f = 0; f < xb->nfields; ++f) { 41452849c42SDave May void *dest; 41552849c42SDave May void *src; 41652849c42SDave May 4179566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldGetAccess(xb->field[f])); 4189566063dSJacob Faibussowitsch if (xb != yb) PetscCall(DMSwarmDataFieldGetAccess(yb->field[f])); 4199566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldAccessPoint(xb->field[f], pid_x, &src)); 4209566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldAccessPoint(yb->field[f], pid_y, &dest)); 4219566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(dest, src, xb->field[f]->atomic_size)); 4229566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldRestoreAccess(xb->field[f])); 4239566063dSJacob Faibussowitsch if (xb != yb) PetscCall(DMSwarmDataFieldRestoreAccess(yb->field[f])); 42452849c42SDave May } 4253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 42652849c42SDave May } 42752849c42SDave May 428d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketCreateFromSubset(DMSwarmDataBucket DBIn, const PetscInt N, const PetscInt list[], DMSwarmDataBucket *DB) 429d71ae5a4SJacob Faibussowitsch { 4305c18a9d6SDave May PetscInt nfields; 43177048351SPatrick Sanan DMSwarmDataField *fields; 4325c18a9d6SDave May PetscInt f, L, buffer, allocated, p; 43352849c42SDave May 434521f74f9SMatthew G. Knepley PetscFunctionBegin; 4359566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketCreate(DB)); 43652849c42SDave May /* copy contents of DBIn */ 4379566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetDMSwarmDataFields(DBIn, &nfields, &fields)); 4389566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(DBIn, &L, &buffer, &allocated)); 43948a46eb9SPierre Jolivet for (f = 0; f < nfields; ++f) PetscCall(DMSwarmDataBucketRegisterField(*DB, "DMSwarmDataBucketCreateFromSubset", fields[f]->name, fields[f]->atomic_size, NULL)); 4409566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketFinalize(*DB)); 4419566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketSetSizes(*DB, L, buffer)); 44219307e5cSMatthew G. Knepley for (f = 0; f < nfields; ++f) { 44319307e5cSMatthew G. Knepley DMSwarmDataField gfield; 44419307e5cSMatthew G. Knepley 44519307e5cSMatthew G. Knepley PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(*DB, fields[f]->name, &gfield)); 44619307e5cSMatthew G. Knepley PetscCall(DMSwarmDataFieldSetBlockSize(gfield, fields[f]->bs)); 44719307e5cSMatthew G. Knepley gfield->petsc_type = fields[f]->petsc_type; 44819307e5cSMatthew G. Knepley } 44952849c42SDave May /* now copy the desired guys from DBIn => DB */ 45048a46eb9SPierre Jolivet for (p = 0; p < N; ++p) PetscCall(DMSwarmDataBucketCopyPoint(DBIn, list[p], *DB, list[p])); 4513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 45252849c42SDave May } 45352849c42SDave May 454d5b43468SJose E. Roman /* insert into an existing location */ 455d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldInsertPoint(const DMSwarmDataField field, const PetscInt index, const void *ctx) 456d71ae5a4SJacob Faibussowitsch { 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)); 4643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 46552849c42SDave May } 46652849c42SDave May 467521f74f9SMatthew G. Knepley /* remove data at index - replace with last point */ 468d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketRemovePointAtIndex(const DMSwarmDataBucket db, const PetscInt index) 469d71ae5a4SJacob Faibussowitsch { 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)); 4963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 49752849c42SDave May } 49852849c42SDave May 49952849c42SDave May /* copy x into y */ 500d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldCopyPoint(const PetscInt pid_x, const DMSwarmDataField field_x, const PetscInt pid_y, const DMSwarmDataField field_y) 501d71ae5a4SJacob Faibussowitsch { 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)); 5123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 51352849c42SDave May } 51452849c42SDave May 515521f74f9SMatthew G. Knepley /* zero only the datafield at this point */ 516d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldZeroPoint(const DMSwarmDataField field, const PetscInt index) 517d71ae5a4SJacob Faibussowitsch { 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)); 5253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 52652849c42SDave May } 52752849c42SDave May 528521f74f9SMatthew G. Knepley /* zero ALL data for this point */ 529d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketZeroPoint(const DMSwarmDataBucket db, const PetscInt index) 530d71ae5a4SJacob Faibussowitsch { 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 } 5413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 54252849c42SDave May } 54352849c42SDave May 54452849c42SDave May /* increment */ 545d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketAddPoint(DMSwarmDataBucket db) 546d71ae5a4SJacob Faibussowitsch { 547521f74f9SMatthew G. Knepley PetscFunctionBegin; 548670292f5SMatthew Knepley PetscCall(DMSwarmDataBucketSetSizes(db, PetscMax(db->L, 0) + 1, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT)); 5493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 55052849c42SDave May } 55152849c42SDave May 5527fbf63aeSDave May /* decrement */ 553d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketRemovePoint(DMSwarmDataBucket db) 554d71ae5a4SJacob Faibussowitsch { 555521f74f9SMatthew G. Knepley PetscFunctionBegin; 556670292f5SMatthew Knepley PetscCheck(db->L > 0, PetscObjectComm((PetscObject)db), PETSC_ERR_ARG_WRONG, "Swarm has no points to be removed"); 5579566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketSetSizes(db, db->L - 1, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT)); 5583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5597fbf63aeSDave May } 5607fbf63aeSDave May 5615627991aSBarry Smith /* Should be redone to user PetscViewer */ 56266976f2fSJacob Faibussowitsch static PetscErrorCode DMSwarmDataBucketView_stdout(MPI_Comm comm, DMSwarmDataBucket db) 563d71ae5a4SJacob Faibussowitsch { 564ab2be9a3SDave May PetscInt f; 565*f1957bc3SPierre Jolivet double memory_usage_total = 0.0; 566ab2be9a3SDave May 567ab2be9a3SDave May PetscFunctionBegin; 5689566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "DMSwarmDataBucketView: \n")); 56963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, " L = %" PetscInt_FMT " \n", db->L)); 57063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, " buffer = %" PetscInt_FMT " \n", db->buffer)); 57163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, " allocated = %" PetscInt_FMT " \n", db->allocated)); 57263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, " nfields registered = %" PetscInt_FMT " \n", db->nfields)); 573ab2be9a3SDave May 574ab2be9a3SDave May for (f = 0; f < db->nfields; ++f) { 575ab2be9a3SDave May double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6; 576*f1957bc3SPierre Jolivet memory_usage_total += memory_usage_f; 577ab2be9a3SDave May } 578*f1957bc3SPierre Jolivet PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &memory_usage_total, 1, MPI_DOUBLE, MPI_SUM, comm)); 579ab2be9a3SDave May 580ab2be9a3SDave May for (f = 0; f < db->nfields; ++f) { 581ab2be9a3SDave May double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6; 58263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, " [%3" PetscInt_FMT "] %15s : Mem. usage = %1.2e (MB) [rank0]\n", f, db->field[f]->name, memory_usage_f)); 58363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, " blocksize = %" PetscInt_FMT " \n", db->field[f]->bs)); 584ab2be9a3SDave May if (db->field[f]->bs != 1) { 58563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, " atomic size = %zu [full block, bs=%" PetscInt_FMT "]\n", db->field[f]->atomic_size, db->field[f]->bs)); 58663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, " atomic size/item = %zu \n", (size_t)(db->field[f]->atomic_size / db->field[f]->bs))); 587ab2be9a3SDave May } else { 5889566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, " atomic size = %zu \n", db->field[f]->atomic_size)); 589ab2be9a3SDave May } 590ab2be9a3SDave May } 5919566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, " Total mem. usage = %1.2e (MB) (collective)\n", memory_usage_total)); 5923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 593ab2be9a3SDave May } 594ab2be9a3SDave May 59566976f2fSJacob Faibussowitsch static PetscErrorCode DMSwarmDataBucketView_Seq(MPI_Comm comm, DMSwarmDataBucket db, const char filename[], DMSwarmDataBucketViewType type) 596d71ae5a4SJacob Faibussowitsch { 597521f74f9SMatthew G. Knepley PetscFunctionBegin; 59852849c42SDave May switch (type) { 599d71ae5a4SJacob Faibussowitsch case DATABUCKET_VIEW_STDOUT: 600d71ae5a4SJacob Faibussowitsch PetscCall(DMSwarmDataBucketView_stdout(PETSC_COMM_SELF, db)); 601d71ae5a4SJacob Faibussowitsch break; 602d71ae5a4SJacob Faibussowitsch case DATABUCKET_VIEW_ASCII: 603d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for ascii output"); 604d71ae5a4SJacob Faibussowitsch case DATABUCKET_VIEW_BINARY: 605d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for binary output"); 606d71ae5a4SJacob Faibussowitsch case DATABUCKET_VIEW_HDF5: 607d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for HDF5 output"); 608d71ae5a4SJacob Faibussowitsch default: 609d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown viewer method requested"); 61052849c42SDave May } 6113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 61252849c42SDave May } 61352849c42SDave May 61466976f2fSJacob Faibussowitsch static PetscErrorCode DMSwarmDataBucketView_MPI(MPI_Comm comm, DMSwarmDataBucket db, const char filename[], DMSwarmDataBucketViewType type) 615d71ae5a4SJacob Faibussowitsch { 616521f74f9SMatthew G. Knepley PetscFunctionBegin; 61752849c42SDave May switch (type) { 618d71ae5a4SJacob Faibussowitsch case DATABUCKET_VIEW_STDOUT: 619d71ae5a4SJacob Faibussowitsch PetscCall(DMSwarmDataBucketView_stdout(comm, db)); 620d71ae5a4SJacob Faibussowitsch break; 621d71ae5a4SJacob Faibussowitsch case DATABUCKET_VIEW_ASCII: 622d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for ascii output"); 623d71ae5a4SJacob Faibussowitsch case DATABUCKET_VIEW_BINARY: 624d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for binary output"); 625d71ae5a4SJacob Faibussowitsch case DATABUCKET_VIEW_HDF5: 626d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for HDF5 output"); 627d71ae5a4SJacob Faibussowitsch default: 628d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown viewer method requested"); 62952849c42SDave May } 6303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 63152849c42SDave May } 63252849c42SDave May 633d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketView(MPI_Comm comm, DMSwarmDataBucket db, const char filename[], DMSwarmDataBucketViewType type) 634d71ae5a4SJacob Faibussowitsch { 635d7d19db6SBarry Smith PetscMPIInt size; 63652849c42SDave May 637521f74f9SMatthew G. Knepley PetscFunctionBegin; 6389566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 639d7d19db6SBarry Smith if (size == 1) { 6409566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketView_Seq(comm, db, filename, type)); 64152849c42SDave May } else { 6429566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketView_MPI(comm, db, filename, type)); 64352849c42SDave May } 6443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 64552849c42SDave May } 64652849c42SDave May 647d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketDuplicateFields(DMSwarmDataBucket dbA, DMSwarmDataBucket *dbB) 648d71ae5a4SJacob Faibussowitsch { 64977048351SPatrick Sanan DMSwarmDataBucket db2; 6505c18a9d6SDave May PetscInt f; 65152849c42SDave May 652521f74f9SMatthew G. Knepley PetscFunctionBegin; 6539566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketCreate(&db2)); 65452849c42SDave May /* copy contents from dbA into db2 */ 655521f74f9SMatthew G. Knepley for (f = 0; f < dbA->nfields; ++f) { 65677048351SPatrick Sanan DMSwarmDataField field; 65752849c42SDave May size_t atomic_size; 65852849c42SDave May char *name; 65952849c42SDave May 66052849c42SDave May field = dbA->field[f]; 66152849c42SDave May atomic_size = field->atomic_size; 66252849c42SDave May name = field->name; 6639566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketRegisterField(db2, "DMSwarmDataBucketDuplicateFields", name, atomic_size, NULL)); 66452849c42SDave May } 6659566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketFinalize(db2)); 6669566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketSetInitialSizes(db2, 0, 1000)); 66752849c42SDave May *dbB = db2; 6683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 66952849c42SDave May } 67052849c42SDave May 67152849c42SDave May /* 67252849c42SDave May Insert points from db2 into db1 67352849c42SDave May db1 <<== db2 67452849c42SDave May */ 675d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketInsertValues(DMSwarmDataBucket db1, DMSwarmDataBucket db2) 676d71ae5a4SJacob Faibussowitsch { 6775c18a9d6SDave May PetscInt n_mp_points1, n_mp_points2; 6785c18a9d6SDave May PetscInt n_mp_points1_new, p; 67952849c42SDave May 680521f74f9SMatthew G. Knepley PetscFunctionBegin; 6819566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(db1, &n_mp_points1, NULL, NULL)); 6829566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(db2, &n_mp_points2, NULL, NULL)); 68352849c42SDave May n_mp_points1_new = n_mp_points1 + n_mp_points2; 6849566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketSetSizes(db1, n_mp_points1_new, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT)); 685521f74f9SMatthew G. Knepley for (p = 0; p < n_mp_points2; ++p) { 686521f74f9SMatthew G. Knepley /* db1 <<== db2 */ 68757508eceSPierre Jolivet PetscCall(DMSwarmDataBucketCopyPoint(db2, p, db1, n_mp_points1 + p)); 68852849c42SDave May } 6893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 69052849c42SDave May } 69152849c42SDave May 69252849c42SDave May /* helpers for parallel send/recv */ 693d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketCreatePackedArray(DMSwarmDataBucket db, size_t *bytes, void **buf) 694d71ae5a4SJacob Faibussowitsch { 6955c18a9d6SDave May PetscInt f; 69652849c42SDave May size_t sizeof_marker_contents; 69752849c42SDave May void *buffer; 69852849c42SDave May 699521f74f9SMatthew G. Knepley PetscFunctionBegin; 70052849c42SDave May sizeof_marker_contents = 0; 701521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 70277048351SPatrick Sanan DMSwarmDataField df = db->field[f]; 70352849c42SDave May sizeof_marker_contents += df->atomic_size; 70452849c42SDave May } 7059566063dSJacob Faibussowitsch PetscCall(PetscMalloc(sizeof_marker_contents, &buffer)); 7069566063dSJacob Faibussowitsch PetscCall(PetscMemzero(buffer, sizeof_marker_contents)); 707ad540459SPierre Jolivet if (bytes) *bytes = sizeof_marker_contents; 708ad540459SPierre Jolivet if (buf) *buf = buffer; 7093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 71052849c42SDave May } 71152849c42SDave May 712d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketDestroyPackedArray(DMSwarmDataBucket db, void **buf) 713d71ae5a4SJacob Faibussowitsch { 714521f74f9SMatthew G. Knepley PetscFunctionBegin; 71552849c42SDave May if (buf) { 7169566063dSJacob Faibussowitsch PetscCall(PetscFree(*buf)); 71752849c42SDave May *buf = NULL; 71852849c42SDave May } 7193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 72052849c42SDave May } 72152849c42SDave May 722d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketFillPackedArray(DMSwarmDataBucket db, const PetscInt index, void *buf) 723d71ae5a4SJacob Faibussowitsch { 7245c18a9d6SDave May PetscInt f; 72552849c42SDave May void *data, *data_p; 72652849c42SDave May size_t asize, offset; 72752849c42SDave May 728521f74f9SMatthew G. Knepley PetscFunctionBegin; 72952849c42SDave May offset = 0; 730521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 73177048351SPatrick Sanan DMSwarmDataField df = db->field[f]; 73252849c42SDave May 73352849c42SDave May asize = df->atomic_size; 734835f2295SStefano Zampini data = df->data; 73552849c42SDave May data_p = (void *)((char *)data + index * asize); 7369566063dSJacob Faibussowitsch PetscCall(PetscMemcpy((void *)((char *)buf + offset), data_p, asize)); 73752849c42SDave May offset = offset + asize; 73852849c42SDave May } 7393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 74052849c42SDave May } 74152849c42SDave May 742d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketInsertPackedArray(DMSwarmDataBucket db, const PetscInt idx, void *data) 743d71ae5a4SJacob Faibussowitsch { 7445c18a9d6SDave May PetscInt f; 74552849c42SDave May void *data_p; 74652849c42SDave May size_t offset; 74752849c42SDave May 748521f74f9SMatthew G. Knepley PetscFunctionBegin; 74952849c42SDave May offset = 0; 750521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 75177048351SPatrick Sanan DMSwarmDataField df = db->field[f]; 75252849c42SDave May 75352849c42SDave May data_p = (void *)((char *)data + offset); 754835f2295SStefano Zampini PetscCall(DMSwarmDataFieldInsertPoint(df, idx, data_p)); 75552849c42SDave May offset = offset + df->atomic_size; 75652849c42SDave May } 7573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 75852849c42SDave May } 759