1279f676cSBarry Smith #include "../src/dm/impls/swarm/data_bucket.h" 252849c42SDave May 352849c42SDave May /* string helpers */ 49371c9d4SSatish Balay PetscErrorCode DMSwarmDataFieldStringInList(const char name[], const PetscInt N, const DMSwarmDataField gfield[], PetscBool *val) { 55c18a9d6SDave May PetscInt i; 652849c42SDave May 7521f74f9SMatthew G. Knepley PetscFunctionBegin; 852849c42SDave May *val = PETSC_FALSE; 9521f74f9SMatthew G. Knepley for (i = 0; i < N; ++i) { 10521f74f9SMatthew G. Knepley PetscBool flg; 119566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(name, gfield[i]->name, &flg)); 12521f74f9SMatthew G. Knepley if (flg) { 1352849c42SDave May *val = PETSC_TRUE; 142eac95f8SDave May PetscFunctionReturn(0); 1552849c42SDave May } 1652849c42SDave May } 172eac95f8SDave May PetscFunctionReturn(0); 1852849c42SDave May } 1952849c42SDave May 209371c9d4SSatish Balay PetscErrorCode DMSwarmDataFieldStringFindInList(const char name[], const PetscInt N, const DMSwarmDataField gfield[], PetscInt *index) { 215c18a9d6SDave May PetscInt i; 2252849c42SDave May 23521f74f9SMatthew G. Knepley PetscFunctionBegin; 2452849c42SDave May *index = -1; 25521f74f9SMatthew G. Knepley for (i = 0; i < N; ++i) { 26521f74f9SMatthew G. Knepley PetscBool flg; 279566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(name, gfield[i]->name, &flg)); 28521f74f9SMatthew G. Knepley if (flg) { 2952849c42SDave May *index = i; 302eac95f8SDave May PetscFunctionReturn(0); 3152849c42SDave May } 3252849c42SDave May } 332eac95f8SDave May PetscFunctionReturn(0); 3452849c42SDave May } 3552849c42SDave May 369371c9d4SSatish Balay PetscErrorCode DMSwarmDataFieldCreate(const char registration_function[], const char name[], const size_t size, const PetscInt L, DMSwarmDataField *DF) { 3777048351SPatrick Sanan DMSwarmDataField df; 3852849c42SDave May 39521f74f9SMatthew G. Knepley PetscFunctionBegin; 409566063dSJacob Faibussowitsch PetscCall(PetscNew(&df)); 419566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(registration_function, &df->registration_function)); 429566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &df->name)); 4352849c42SDave May df->atomic_size = size; 4452849c42SDave May df->L = L; 4552c3ed93SDave May df->bs = 1; 46521f74f9SMatthew G. Knepley /* allocate something so we don't have to reallocate */ 479566063dSJacob Faibussowitsch PetscCall(PetscMalloc(size * L, &df->data)); 489566063dSJacob Faibussowitsch PetscCall(PetscMemzero(df->data, size * L)); 4952849c42SDave May *DF = df; 502eac95f8SDave May PetscFunctionReturn(0); 5152849c42SDave May } 5252849c42SDave May 539371c9d4SSatish Balay PetscErrorCode DMSwarmDataFieldDestroy(DMSwarmDataField *DF) { 5477048351SPatrick Sanan DMSwarmDataField df = *DF; 5552849c42SDave May 56521f74f9SMatthew G. Knepley PetscFunctionBegin; 579566063dSJacob Faibussowitsch PetscCall(PetscFree(df->registration_function)); 589566063dSJacob Faibussowitsch PetscCall(PetscFree(df->name)); 599566063dSJacob Faibussowitsch PetscCall(PetscFree(df->data)); 609566063dSJacob Faibussowitsch PetscCall(PetscFree(df)); 6152849c42SDave May *DF = NULL; 622eac95f8SDave May PetscFunctionReturn(0); 6352849c42SDave May } 6452849c42SDave May 6552849c42SDave May /* data bucket */ 669371c9d4SSatish Balay PetscErrorCode DMSwarmDataBucketCreate(DMSwarmDataBucket *DB) { 6777048351SPatrick Sanan DMSwarmDataBucket db; 6852849c42SDave May 69521f74f9SMatthew G. Knepley PetscFunctionBegin; 709566063dSJacob Faibussowitsch PetscCall(PetscNew(&db)); 7152849c42SDave May 7252849c42SDave May db->finalised = PETSC_FALSE; 7352849c42SDave May /* create empty spaces for fields */ 743454631fSDave May db->L = -1; 7552849c42SDave May db->buffer = 1; 7652849c42SDave May db->allocated = 1; 7752849c42SDave May db->nfields = 0; 789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1, &db->field)); 7952849c42SDave May *DB = db; 802eac95f8SDave May PetscFunctionReturn(0); 8152849c42SDave May } 8252849c42SDave May 839371c9d4SSatish Balay PetscErrorCode DMSwarmDataBucketDestroy(DMSwarmDataBucket *DB) { 8477048351SPatrick Sanan DMSwarmDataBucket db = *DB; 855c18a9d6SDave May PetscInt f; 8652849c42SDave May 87521f74f9SMatthew G. Knepley PetscFunctionBegin; 8852849c42SDave May /* release fields */ 89*48a46eb9SPierre Jolivet for (f = 0; f < db->nfields; ++f) PetscCall(DMSwarmDataFieldDestroy(&db->field[f])); 9052849c42SDave May /* this will catch the initially allocated objects in the event that no fields are registered */ 91*48a46eb9SPierre Jolivet if (db->field != NULL) PetscCall(PetscFree(db->field)); 929566063dSJacob Faibussowitsch PetscCall(PetscFree(db)); 9352849c42SDave May *DB = NULL; 942eac95f8SDave May PetscFunctionReturn(0); 9552849c42SDave May } 9652849c42SDave May 979371c9d4SSatish Balay PetscErrorCode DMSwarmDataBucketQueryForActiveFields(DMSwarmDataBucket db, PetscBool *any_active_fields) { 980cb17e37SDave May PetscInt f; 99521f74f9SMatthew G. Knepley 100521f74f9SMatthew G. Knepley PetscFunctionBegin; 1010cb17e37SDave May *any_active_fields = PETSC_FALSE; 102521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 1030cb17e37SDave May if (db->field[f]->active) { 1040cb17e37SDave May *any_active_fields = PETSC_TRUE; 1050cb17e37SDave May PetscFunctionReturn(0); 1060cb17e37SDave May } 1070cb17e37SDave May } 1080cb17e37SDave May PetscFunctionReturn(0); 1090cb17e37SDave May } 1100cb17e37SDave May 1119371c9d4SSatish Balay PetscErrorCode DMSwarmDataBucketRegisterField(DMSwarmDataBucket db, const char registration_function[], const char field_name[], size_t atomic_size, DMSwarmDataField *_gfield) { 11252849c42SDave May PetscBool val; 11377048351SPatrick Sanan DMSwarmDataField fp; 11452849c42SDave May 115521f74f9SMatthew G. Knepley PetscFunctionBegin; 11652849c42SDave May /* check we haven't finalised the registration of fields */ 11752849c42SDave May /* 11852849c42SDave May if (db->finalised==PETSC_TRUE) { 11977048351SPatrick Sanan printf("ERROR: DMSwarmDataBucketFinalize() has been called. Cannot register more fields\n"); 12052849c42SDave May ERROR(); 12152849c42SDave May } 12252849c42SDave May */ 12352849c42SDave May /* check for repeated name */ 1249566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldStringInList(field_name, db->nfields, (const DMSwarmDataField *)db->field, &val)); 12508401ef6SPierre Jolivet PetscCheck(val != PETSC_TRUE, PETSC_COMM_SELF, PETSC_ERR_USER, "Field %s already exists. Cannot add same field twice", field_name); 12652849c42SDave May /* create new space for data */ 1279566063dSJacob Faibussowitsch PetscCall(PetscRealloc(sizeof(DMSwarmDataField) * (db->nfields + 1), &db->field)); 12852849c42SDave May /* add field */ 1299566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldCreate(registration_function, field_name, atomic_size, db->allocated, &fp)); 13052849c42SDave May db->field[db->nfields] = fp; 13152849c42SDave May db->nfields++; 1329371c9d4SSatish Balay if (_gfield != NULL) { *_gfield = fp; } 1332eac95f8SDave May PetscFunctionReturn(0); 13452849c42SDave May } 13552849c42SDave May 13652849c42SDave May /* 13777048351SPatrick Sanan #define DMSwarmDataBucketRegisterField(db,name,size,k) {\ 13852849c42SDave May char *location;\ 13952849c42SDave May asprintf(&location,"Registered by %s() at line %d within file %s", __FUNCTION__, __LINE__, __FILE__);\ 14077048351SPatrick Sanan _DMSwarmDataBucketRegisterField( (db), location, (name), (size), (k));\ 141521f74f9SMatthew G. Knepley ierr = PetscFree(location);\ 14252849c42SDave May } 14352849c42SDave May */ 14452849c42SDave May 1459371c9d4SSatish Balay PetscErrorCode DMSwarmDataBucketGetDMSwarmDataFieldIdByName(DMSwarmDataBucket db, const char name[], PetscInt *idx) { 1462e956fe4SStefano Zampini PetscBool found; 1472e956fe4SStefano Zampini 1482e956fe4SStefano Zampini PetscFunctionBegin; 1492e956fe4SStefano Zampini *idx = -1; 1502e956fe4SStefano Zampini PetscCall(DMSwarmDataFieldStringInList(name, db->nfields, (const DMSwarmDataField *)db->field, &found)); 1512e956fe4SStefano Zampini PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot find DMSwarmDataField with name %s", name); 1522e956fe4SStefano Zampini PetscCall(DMSwarmDataFieldStringFindInList(name, db->nfields, (const DMSwarmDataField *)db->field, idx)); 1532e956fe4SStefano Zampini PetscFunctionReturn(0); 1542e956fe4SStefano Zampini } 1552e956fe4SStefano Zampini 1569371c9d4SSatish Balay PetscErrorCode DMSwarmDataBucketGetDMSwarmDataFieldByName(DMSwarmDataBucket db, const char name[], DMSwarmDataField *gfield) { 1575c18a9d6SDave May PetscInt idx; 15852849c42SDave May PetscBool found; 15952849c42SDave May 160521f74f9SMatthew G. Knepley PetscFunctionBegin; 1619566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldStringInList(name, db->nfields, (const DMSwarmDataField *)db->field, &found)); 16228b400f6SJacob Faibussowitsch PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot find DMSwarmDataField with name %s", name); 1639566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldStringFindInList(name, db->nfields, (const DMSwarmDataField *)db->field, &idx)); 16452849c42SDave May *gfield = db->field[idx]; 1652eac95f8SDave May PetscFunctionReturn(0); 16652849c42SDave May } 16752849c42SDave May 1689371c9d4SSatish Balay PetscErrorCode DMSwarmDataBucketQueryDMSwarmDataFieldByName(DMSwarmDataBucket db, const char name[], PetscBool *found) { 169521f74f9SMatthew G. Knepley PetscFunctionBegin; 17052849c42SDave May *found = PETSC_FALSE; 1719566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldStringInList(name, db->nfields, (const DMSwarmDataField *)db->field, found)); 1727fbf63aeSDave May PetscFunctionReturn(0); 17352849c42SDave May } 17452849c42SDave May 1759371c9d4SSatish Balay PetscErrorCode DMSwarmDataBucketFinalize(DMSwarmDataBucket db) { 176521f74f9SMatthew G. Knepley PetscFunctionBegin; 17752849c42SDave May db->finalised = PETSC_TRUE; 1787fbf63aeSDave May PetscFunctionReturn(0); 17952849c42SDave May } 18052849c42SDave May 1819371c9d4SSatish Balay PetscErrorCode DMSwarmDataFieldGetNumEntries(DMSwarmDataField df, PetscInt *sum) { 182521f74f9SMatthew G. Knepley PetscFunctionBegin; 18352849c42SDave May *sum = df->L; 1847fbf63aeSDave May PetscFunctionReturn(0); 18552849c42SDave May } 18652849c42SDave May 1879371c9d4SSatish Balay PetscErrorCode DMSwarmDataFieldSetBlockSize(DMSwarmDataField df, PetscInt blocksize) { 188521f74f9SMatthew G. Knepley PetscFunctionBegin; 18952c3ed93SDave May df->bs = blocksize; 19052c3ed93SDave May PetscFunctionReturn(0); 19152c3ed93SDave May } 19252c3ed93SDave May 1939371c9d4SSatish Balay PetscErrorCode DMSwarmDataFieldSetSize(DMSwarmDataField df, const PetscInt new_L) { 194521f74f9SMatthew G. Knepley PetscFunctionBegin; 19508401ef6SPierre Jolivet PetscCheck(new_L >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot set size of DMSwarmDataField to be < 0"); 1967fbf63aeSDave May if (new_L == df->L) PetscFunctionReturn(0); 19752849c42SDave May if (new_L > df->L) { 1989566063dSJacob Faibussowitsch PetscCall(PetscRealloc(df->atomic_size * (new_L), &df->data)); 19952849c42SDave May /* init new contents */ 2009566063dSJacob Faibussowitsch PetscCall(PetscMemzero((((char *)df->data) + df->L * df->atomic_size), (new_L - df->L) * df->atomic_size)); 2017fbf63aeSDave May } else { 20252849c42SDave May /* reallocate pointer list, add +1 in case new_L = 0 */ 2039566063dSJacob Faibussowitsch PetscCall(PetscRealloc(df->atomic_size * (new_L + 1), &df->data)); 20452849c42SDave May } 20552849c42SDave May df->L = new_L; 2067fbf63aeSDave May PetscFunctionReturn(0); 20752849c42SDave May } 20852849c42SDave May 2099371c9d4SSatish Balay PetscErrorCode DMSwarmDataFieldZeroBlock(DMSwarmDataField df, const PetscInt start, const PetscInt end) { 210521f74f9SMatthew G. Knepley PetscFunctionBegin; 21163a3b9bcSJacob 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); 21263a3b9bcSJacob Faibussowitsch PetscCheck(start >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot zero a block of entries if start(%" PetscInt_FMT ") < 0", start); 21363a3b9bcSJacob 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); 2149566063dSJacob Faibussowitsch PetscCall(PetscMemzero((((char *)df->data) + start * df->atomic_size), (end - start) * df->atomic_size)); 2157fbf63aeSDave May PetscFunctionReturn(0); 21652849c42SDave May } 21752849c42SDave May 21852849c42SDave May /* 21952849c42SDave May A negative buffer value will simply be ignored and the old buffer value will be used. 22052849c42SDave May */ 2219371c9d4SSatish Balay PetscErrorCode DMSwarmDataBucketSetSizes(DMSwarmDataBucket db, const PetscInt L, const PetscInt buffer) { 2225c18a9d6SDave May PetscInt current_allocated, new_used, new_unused, new_buffer, new_allocated, f; 2230cb17e37SDave May PetscBool any_active_fields; 22452849c42SDave May 225521f74f9SMatthew G. Knepley PetscFunctionBegin; 22608401ef6SPierre Jolivet PetscCheck(db->finalised != PETSC_FALSE, PETSC_COMM_SELF, PETSC_ERR_USER, "You must call DMSwarmDataBucketFinalize() before DMSwarmDataBucketSetSizes()"); 2279566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketQueryForActiveFields(db, &any_active_fields)); 22828b400f6SJacob Faibussowitsch PetscCheck(!any_active_fields, PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot safely re-size as at least one DMSwarmDataField is currently being accessed"); 2290cb17e37SDave May 23052849c42SDave May current_allocated = db->allocated; 23152849c42SDave May new_used = L; 23252849c42SDave May new_unused = current_allocated - new_used; 23352849c42SDave May new_buffer = db->buffer; 23452849c42SDave May if (buffer >= 0) { /* update the buffer value */ 23552849c42SDave May new_buffer = buffer; 23652849c42SDave May } 23752849c42SDave May new_allocated = new_used + new_buffer; 23852849c42SDave May /* action */ 23952849c42SDave May if (new_allocated > current_allocated) { 24052849c42SDave May /* increase size to new_used + new_buffer */ 241*48a46eb9SPierre Jolivet for (f = 0; f < db->nfields; f++) PetscCall(DMSwarmDataFieldSetSize(db->field[f], new_allocated)); 24252849c42SDave May db->L = new_used; 24352849c42SDave May db->buffer = new_buffer; 24452849c42SDave May db->allocated = new_used + new_buffer; 245521f74f9SMatthew G. Knepley } else { 24652849c42SDave May if (new_unused > 2 * new_buffer) { 24752849c42SDave May /* shrink array to new_used + new_buffer */ 248*48a46eb9SPierre Jolivet for (f = 0; f < db->nfields; ++f) PetscCall(DMSwarmDataFieldSetSize(db->field[f], new_allocated)); 24952849c42SDave May db->L = new_used; 25052849c42SDave May db->buffer = new_buffer; 25152849c42SDave May db->allocated = new_used + new_buffer; 252521f74f9SMatthew G. Knepley } else { 25352849c42SDave May db->L = new_used; 25452849c42SDave May db->buffer = new_buffer; 25552849c42SDave May } 25652849c42SDave May } 25752849c42SDave May /* zero all entries from db->L to db->allocated */ 258521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 25977048351SPatrick Sanan DMSwarmDataField field = db->field[f]; 2609566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldZeroBlock(field, db->L, db->allocated)); 26152849c42SDave May } 2627fbf63aeSDave May PetscFunctionReturn(0); 26352849c42SDave May } 26452849c42SDave May 2659371c9d4SSatish Balay PetscErrorCode DMSwarmDataBucketSetInitialSizes(DMSwarmDataBucket db, const PetscInt L, const PetscInt buffer) { 2665c18a9d6SDave May PetscInt f; 267dbe06d34SDave May 268521f74f9SMatthew G. Knepley PetscFunctionBegin; 2699566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketSetSizes(db, L, buffer)); 270521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 27177048351SPatrick Sanan DMSwarmDataField field = db->field[f]; 2729566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldZeroBlock(field, 0, db->allocated)); 27352849c42SDave May } 2747fbf63aeSDave May PetscFunctionReturn(0); 27552849c42SDave May } 27652849c42SDave May 2779371c9d4SSatish Balay PetscErrorCode DMSwarmDataBucketGetSizes(DMSwarmDataBucket db, PetscInt *L, PetscInt *buffer, PetscInt *allocated) { 278521f74f9SMatthew G. Knepley PetscFunctionBegin; 27952849c42SDave May if (L) { *L = db->L; } 28052849c42SDave May if (buffer) { *buffer = db->buffer; } 28152849c42SDave May if (allocated) { *allocated = db->allocated; } 2827fbf63aeSDave May PetscFunctionReturn(0); 28352849c42SDave May } 28452849c42SDave May 2859371c9d4SSatish Balay PetscErrorCode DMSwarmDataBucketGetGlobalSizes(MPI_Comm comm, DMSwarmDataBucket db, PetscInt *L, PetscInt *buffer, PetscInt *allocated) { 286521f74f9SMatthew G. Knepley PetscFunctionBegin; 2879566063dSJacob Faibussowitsch if (L) PetscCallMPI(MPI_Allreduce(&db->L, L, 1, MPIU_INT, MPI_SUM, comm)); 2889566063dSJacob Faibussowitsch if (buffer) PetscCallMPI(MPI_Allreduce(&db->buffer, buffer, 1, MPIU_INT, MPI_SUM, comm)); 2899566063dSJacob Faibussowitsch if (allocated) PetscCallMPI(MPI_Allreduce(&db->allocated, allocated, 1, MPIU_INT, MPI_SUM, comm)); 2907fbf63aeSDave May PetscFunctionReturn(0); 29152849c42SDave May } 29252849c42SDave May 2939371c9d4SSatish Balay PetscErrorCode DMSwarmDataBucketGetDMSwarmDataFields(DMSwarmDataBucket db, PetscInt *L, DMSwarmDataField *fields[]) { 294521f74f9SMatthew G. Knepley PetscFunctionBegin; 29552849c42SDave May if (L) { *L = db->nfields; } 29652849c42SDave May if (fields) { *fields = db->field; } 2977fbf63aeSDave May PetscFunctionReturn(0); 29852849c42SDave May } 29952849c42SDave May 3009371c9d4SSatish Balay PetscErrorCode DMSwarmDataFieldGetAccess(const DMSwarmDataField gfield) { 301521f74f9SMatthew G. Knepley PetscFunctionBegin; 30228b400f6SJacob Faibussowitsch PetscCheck(!gfield->active, PETSC_COMM_SELF, PETSC_ERR_USER, "Field \"%s\" is already active. You must call DMSwarmDataFieldRestoreAccess()", gfield->name); 30352849c42SDave May gfield->active = PETSC_TRUE; 3047fbf63aeSDave May PetscFunctionReturn(0); 30552849c42SDave May } 30652849c42SDave May 3079371c9d4SSatish Balay PetscErrorCode DMSwarmDataFieldAccessPoint(const DMSwarmDataField gfield, const PetscInt pid, void **ctx_p) { 308521f74f9SMatthew G. Knepley PetscFunctionBegin; 30972505a4dSJed Brown *ctx_p = NULL; 310497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD) 31152849c42SDave May /* debug mode */ 31284bcda08SDave May /* check point is valid */ 31308401ef6SPierre Jolivet PetscCheck(pid >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be >= 0"); 31463a3b9bcSJacob Faibussowitsch PetscCheck(pid < gfield->L, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be < %" PetscInt_FMT, gfield->L); 31508401ef6SPierre 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); 31652849c42SDave May #endif 31777048351SPatrick Sanan *ctx_p = DMSWARM_DATAFIELD_point_access(gfield->data, pid, gfield->atomic_size); 3187fbf63aeSDave May PetscFunctionReturn(0); 31952849c42SDave May } 32052849c42SDave May 3219371c9d4SSatish Balay PetscErrorCode DMSwarmDataFieldAccessPointOffset(const DMSwarmDataField gfield, const size_t offset, const PetscInt pid, void **ctx_p) { 322521f74f9SMatthew G. Knepley PetscFunctionBegin; 323497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD) 32452849c42SDave May /* debug mode */ 32584bcda08SDave May /* check point is valid */ 32608401ef6SPierre Jolivet /* PetscCheck(offset >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be >= 0");*/ 327521f74f9SMatthew G. Knepley /* Note compiler realizes this can never happen with an unsigned PetscInt */ 32808401ef6SPierre Jolivet PetscCheck(offset < gfield->atomic_size, PETSC_COMM_SELF, PETSC_ERR_USER, "offset must be < %zu", gfield->atomic_size); 32984bcda08SDave May /* check point is valid */ 33008401ef6SPierre Jolivet PetscCheck(pid >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be >= 0"); 33163a3b9bcSJacob Faibussowitsch PetscCheck(pid < gfield->L, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be < %" PetscInt_FMT, gfield->L); 33208401ef6SPierre 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); 33352849c42SDave May #endif 33477048351SPatrick Sanan *ctx_p = DMSWARM_DATAFIELD_point_access_offset(gfield->data, pid, gfield->atomic_size, offset); 3357fbf63aeSDave May PetscFunctionReturn(0); 33652849c42SDave May } 33752849c42SDave May 3389371c9d4SSatish Balay PetscErrorCode DMSwarmDataFieldRestoreAccess(DMSwarmDataField gfield) { 339521f74f9SMatthew G. Knepley PetscFunctionBegin; 34008401ef6SPierre Jolivet PetscCheck(gfield->active != PETSC_FALSE, PETSC_COMM_SELF, PETSC_ERR_USER, "Field \"%s\" is not active. You must call DMSwarmDataFieldGetAccess()", gfield->name); 34152849c42SDave May gfield->active = PETSC_FALSE; 3427fbf63aeSDave May PetscFunctionReturn(0); 34352849c42SDave May } 34452849c42SDave May 3459371c9d4SSatish Balay PetscErrorCode DMSwarmDataFieldVerifyAccess(const DMSwarmDataField gfield, const size_t size) { 346521f74f9SMatthew G. Knepley PetscFunctionBegin; 347497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD) 34808401ef6SPierre 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); 34952849c42SDave May #endif 3507fbf63aeSDave May PetscFunctionReturn(0); 35152849c42SDave May } 35252849c42SDave May 3539371c9d4SSatish Balay PetscErrorCode DMSwarmDataFieldGetAtomicSize(const DMSwarmDataField gfield, size_t *size) { 354521f74f9SMatthew G. Knepley PetscFunctionBegin; 35552849c42SDave May if (size) { *size = gfield->atomic_size; } 3567fbf63aeSDave May PetscFunctionReturn(0); 35752849c42SDave May } 35852849c42SDave May 3599371c9d4SSatish Balay PetscErrorCode DMSwarmDataFieldGetEntries(const DMSwarmDataField gfield, void **data) { 360521f74f9SMatthew G. Knepley PetscFunctionBegin; 361521f74f9SMatthew G. Knepley if (data) { *data = gfield->data; } 3627fbf63aeSDave May PetscFunctionReturn(0); 36352849c42SDave May } 36452849c42SDave May 3659371c9d4SSatish Balay PetscErrorCode DMSwarmDataFieldRestoreEntries(const DMSwarmDataField gfield, void **data) { 366521f74f9SMatthew G. Knepley PetscFunctionBegin; 367521f74f9SMatthew G. Knepley if (data) { *data = NULL; } 3687fbf63aeSDave May PetscFunctionReturn(0); 36952849c42SDave May } 37052849c42SDave May 37152849c42SDave May /* y = x */ 3729371c9d4SSatish Balay PetscErrorCode DMSwarmDataBucketCopyPoint(const DMSwarmDataBucket xb, const PetscInt pid_x, const DMSwarmDataBucket yb, const PetscInt pid_y) { 3735c18a9d6SDave May PetscInt f; 374dbe06d34SDave May 375521f74f9SMatthew G. Knepley PetscFunctionBegin; 376521f74f9SMatthew G. Knepley for (f = 0; f < xb->nfields; ++f) { 37752849c42SDave May void *dest; 37852849c42SDave May void *src; 37952849c42SDave May 3809566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldGetAccess(xb->field[f])); 3819566063dSJacob Faibussowitsch if (xb != yb) PetscCall(DMSwarmDataFieldGetAccess(yb->field[f])); 3829566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldAccessPoint(xb->field[f], pid_x, &src)); 3839566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldAccessPoint(yb->field[f], pid_y, &dest)); 3849566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(dest, src, xb->field[f]->atomic_size)); 3859566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldRestoreAccess(xb->field[f])); 3869566063dSJacob Faibussowitsch if (xb != yb) PetscCall(DMSwarmDataFieldRestoreAccess(yb->field[f])); 38752849c42SDave May } 3887fbf63aeSDave May PetscFunctionReturn(0); 38952849c42SDave May } 39052849c42SDave May 3919371c9d4SSatish Balay PetscErrorCode DMSwarmDataBucketCreateFromSubset(DMSwarmDataBucket DBIn, const PetscInt N, const PetscInt list[], DMSwarmDataBucket *DB) { 3925c18a9d6SDave May PetscInt nfields; 39377048351SPatrick Sanan DMSwarmDataField *fields; 3945c18a9d6SDave May PetscInt f, L, buffer, allocated, p; 39552849c42SDave May 396521f74f9SMatthew G. Knepley PetscFunctionBegin; 3979566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketCreate(DB)); 39852849c42SDave May /* copy contents of DBIn */ 3999566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetDMSwarmDataFields(DBIn, &nfields, &fields)); 4009566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(DBIn, &L, &buffer, &allocated)); 401*48a46eb9SPierre Jolivet for (f = 0; f < nfields; ++f) PetscCall(DMSwarmDataBucketRegisterField(*DB, "DMSwarmDataBucketCreateFromSubset", fields[f]->name, fields[f]->atomic_size, NULL)); 4029566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketFinalize(*DB)); 4039566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketSetSizes(*DB, L, buffer)); 40452849c42SDave May /* now copy the desired guys from DBIn => DB */ 405*48a46eb9SPierre Jolivet for (p = 0; p < N; ++p) PetscCall(DMSwarmDataBucketCopyPoint(DBIn, list[p], *DB, list[p])); 4067fbf63aeSDave May PetscFunctionReturn(0); 40752849c42SDave May } 40852849c42SDave May 409521f74f9SMatthew G. Knepley /* insert into an exisitng location */ 4109371c9d4SSatish Balay PetscErrorCode DMSwarmDataFieldInsertPoint(const DMSwarmDataField field, const PetscInt index, const void *ctx) { 411521f74f9SMatthew G. Knepley PetscFunctionBegin; 412497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD) 41384bcda08SDave May /* check point is valid */ 41408401ef6SPierre Jolivet PetscCheck(index >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be >= 0"); 41563a3b9bcSJacob Faibussowitsch PetscCheck(index < field->L, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be < %" PetscInt_FMT, field->L); 41652849c42SDave May #endif 4179566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(DMSWARM_DATAFIELD_point_access(field->data, index, field->atomic_size), ctx, field->atomic_size)); 4187fbf63aeSDave May PetscFunctionReturn(0); 41952849c42SDave May } 42052849c42SDave May 421521f74f9SMatthew G. Knepley /* remove data at index - replace with last point */ 4229371c9d4SSatish Balay PetscErrorCode DMSwarmDataBucketRemovePointAtIndex(const DMSwarmDataBucket db, const PetscInt index) { 4235c18a9d6SDave May PetscInt f; 4240cb17e37SDave May PetscBool any_active_fields; 42552849c42SDave May 426521f74f9SMatthew G. Knepley PetscFunctionBegin; 427497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD) 42884bcda08SDave May /* check point is valid */ 42908401ef6SPierre Jolivet PetscCheck(index >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be >= 0"); 43063a3b9bcSJacob Faibussowitsch PetscCheck(index < db->allocated, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be < %" PetscInt_FMT, db->L + db->buffer); 43152849c42SDave May #endif 4329566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketQueryForActiveFields(db, &any_active_fields)); 43328b400f6SJacob Faibussowitsch PetscCheck(!any_active_fields, PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot safely remove point as at least one DMSwarmDataField is currently being accessed"); 434a233d522SDave May if (index >= db->L) { /* this point is not in the list - no need to error, but I will anyway */ 43563a3b9bcSJacob 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); 43652849c42SDave May } 437a233d522SDave May if (index != db->L - 1) { /* not last point in list */ 438521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 43977048351SPatrick Sanan DMSwarmDataField field = db->field[f]; 44052849c42SDave May 44152849c42SDave May /* copy then remove */ 4429566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldCopyPoint(db->L - 1, field, index, field)); 44377048351SPatrick Sanan /* DMSwarmDataFieldZeroPoint(field,index); */ 44452849c42SDave May } 44552849c42SDave May } 44652849c42SDave May /* decrement size */ 44752849c42SDave May /* this will zero out an crap at the end of the list */ 4489566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketRemovePoint(db)); 4497fbf63aeSDave May PetscFunctionReturn(0); 45052849c42SDave May } 45152849c42SDave May 45252849c42SDave May /* copy x into y */ 4539371c9d4SSatish Balay PetscErrorCode DMSwarmDataFieldCopyPoint(const PetscInt pid_x, const DMSwarmDataField field_x, const PetscInt pid_y, const DMSwarmDataField field_y) { 454521f74f9SMatthew G. Knepley PetscFunctionBegin; 455497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD) 45684bcda08SDave May /* check point is valid */ 45708401ef6SPierre Jolivet PetscCheck(pid_x >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "(IN) index must be >= 0"); 45863a3b9bcSJacob Faibussowitsch PetscCheck(pid_x < field_x->L, PETSC_COMM_SELF, PETSC_ERR_USER, "(IN) index must be < %" PetscInt_FMT, field_x->L); 45908401ef6SPierre Jolivet PetscCheck(pid_y >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "(OUT) index must be >= 0"); 46063a3b9bcSJacob Faibussowitsch PetscCheck(pid_y < field_y->L, PETSC_COMM_SELF, PETSC_ERR_USER, "(OUT) index must be < %" PetscInt_FMT, field_y->L); 46108401ef6SPierre Jolivet PetscCheck(field_y->atomic_size == field_x->atomic_size, PETSC_COMM_SELF, PETSC_ERR_USER, "atomic size must match"); 46252849c42SDave May #endif 4639566063dSJacob 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)); 4647fbf63aeSDave May PetscFunctionReturn(0); 46552849c42SDave May } 46652849c42SDave May 467521f74f9SMatthew G. Knepley /* zero only the datafield at this point */ 4689371c9d4SSatish Balay PetscErrorCode DMSwarmDataFieldZeroPoint(const DMSwarmDataField field, const PetscInt index) { 469521f74f9SMatthew G. Knepley PetscFunctionBegin; 470497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD) 47184bcda08SDave May /* check point is valid */ 47208401ef6SPierre Jolivet PetscCheck(index >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be >= 0"); 47363a3b9bcSJacob Faibussowitsch PetscCheck(index < field->L, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be < %" PetscInt_FMT, field->L); 47452849c42SDave May #endif 4759566063dSJacob Faibussowitsch PetscCall(PetscMemzero(DMSWARM_DATAFIELD_point_access(field->data, index, field->atomic_size), field->atomic_size)); 4767fbf63aeSDave May PetscFunctionReturn(0); 47752849c42SDave May } 47852849c42SDave May 479521f74f9SMatthew G. Knepley /* zero ALL data for this point */ 4809371c9d4SSatish Balay PetscErrorCode DMSwarmDataBucketZeroPoint(const DMSwarmDataBucket db, const PetscInt index) { 4815c18a9d6SDave May PetscInt f; 48252849c42SDave May 483521f74f9SMatthew G. Knepley PetscFunctionBegin; 48484bcda08SDave May /* check point is valid */ 48508401ef6SPierre Jolivet PetscCheck(index >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be >= 0"); 48663a3b9bcSJacob Faibussowitsch PetscCheck(index < db->allocated, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be < %" PetscInt_FMT, db->allocated); 487521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 48877048351SPatrick Sanan DMSwarmDataField field = db->field[f]; 4899566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldZeroPoint(field, index)); 49052849c42SDave May } 4917fbf63aeSDave May PetscFunctionReturn(0); 49252849c42SDave May } 49352849c42SDave May 49452849c42SDave May /* increment */ 4959371c9d4SSatish Balay PetscErrorCode DMSwarmDataBucketAddPoint(DMSwarmDataBucket db) { 496521f74f9SMatthew G. Knepley PetscFunctionBegin; 4979566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketSetSizes(db, db->L + 1, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT)); 4987fbf63aeSDave May PetscFunctionReturn(0); 49952849c42SDave May } 50052849c42SDave May 5017fbf63aeSDave May /* decrement */ 5029371c9d4SSatish Balay PetscErrorCode DMSwarmDataBucketRemovePoint(DMSwarmDataBucket db) { 503521f74f9SMatthew G. Knepley PetscFunctionBegin; 5049566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketSetSizes(db, db->L - 1, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT)); 5057fbf63aeSDave May PetscFunctionReturn(0); 5067fbf63aeSDave May } 5077fbf63aeSDave May 5085627991aSBarry Smith /* Should be redone to user PetscViewer */ 5099371c9d4SSatish Balay PetscErrorCode DMSwarmDataBucketView_stdout(MPI_Comm comm, DMSwarmDataBucket db) { 510ab2be9a3SDave May PetscInt f; 511ab2be9a3SDave May double memory_usage_total, memory_usage_total_local = 0.0; 512ab2be9a3SDave May 513ab2be9a3SDave May PetscFunctionBegin; 5149566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "DMSwarmDataBucketView: \n")); 51563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, " L = %" PetscInt_FMT " \n", db->L)); 51663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, " buffer = %" PetscInt_FMT " \n", db->buffer)); 51763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, " allocated = %" PetscInt_FMT " \n", db->allocated)); 51863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, " nfields registered = %" PetscInt_FMT " \n", db->nfields)); 519ab2be9a3SDave May 520ab2be9a3SDave May for (f = 0; f < db->nfields; ++f) { 521ab2be9a3SDave May double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6; 522ab2be9a3SDave May memory_usage_total_local += memory_usage_f; 523ab2be9a3SDave May } 5249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&memory_usage_total_local, &memory_usage_total, 1, MPI_DOUBLE, MPI_SUM, comm)); 525ab2be9a3SDave May 526ab2be9a3SDave May for (f = 0; f < db->nfields; ++f) { 527ab2be9a3SDave May double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6; 52863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, " [%3" PetscInt_FMT "] %15s : Mem. usage = %1.2e (MB) [rank0]\n", f, db->field[f]->name, memory_usage_f)); 52963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, " blocksize = %" PetscInt_FMT " \n", db->field[f]->bs)); 530ab2be9a3SDave May if (db->field[f]->bs != 1) { 53163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, " atomic size = %zu [full block, bs=%" PetscInt_FMT "]\n", db->field[f]->atomic_size, db->field[f]->bs)); 53263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, " atomic size/item = %zu \n", (size_t)(db->field[f]->atomic_size / db->field[f]->bs))); 533ab2be9a3SDave May } else { 5349566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, " atomic size = %zu \n", db->field[f]->atomic_size)); 535ab2be9a3SDave May } 536ab2be9a3SDave May } 5379566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, " Total mem. usage = %1.2e (MB) (collective)\n", memory_usage_total)); 538ab2be9a3SDave May PetscFunctionReturn(0); 539ab2be9a3SDave May } 540ab2be9a3SDave May 5419371c9d4SSatish Balay PetscErrorCode DMSwarmDataBucketView_Seq(MPI_Comm comm, DMSwarmDataBucket db, const char filename[], DMSwarmDataBucketViewType type) { 542521f74f9SMatthew G. Knepley PetscFunctionBegin; 54352849c42SDave May switch (type) { 5449371c9d4SSatish Balay case DATABUCKET_VIEW_STDOUT: PetscCall(DMSwarmDataBucketView_stdout(PETSC_COMM_SELF, db)); break; 5459371c9d4SSatish Balay case DATABUCKET_VIEW_ASCII: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for ascii output"); 5469371c9d4SSatish Balay case DATABUCKET_VIEW_BINARY: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for binary output"); 5479371c9d4SSatish Balay case DATABUCKET_VIEW_HDF5: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for HDF5 output"); 548521f74f9SMatthew G. Knepley default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown viewer method requested"); 54952849c42SDave May } 5507fbf63aeSDave May PetscFunctionReturn(0); 55152849c42SDave May } 55252849c42SDave May 5539371c9d4SSatish Balay PetscErrorCode DMSwarmDataBucketView_MPI(MPI_Comm comm, DMSwarmDataBucket db, const char filename[], DMSwarmDataBucketViewType type) { 554521f74f9SMatthew G. Knepley PetscFunctionBegin; 55552849c42SDave May switch (type) { 5569371c9d4SSatish Balay case DATABUCKET_VIEW_STDOUT: PetscCall(DMSwarmDataBucketView_stdout(comm, db)); break; 5579371c9d4SSatish Balay case DATABUCKET_VIEW_ASCII: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for ascii output"); 5589371c9d4SSatish Balay case DATABUCKET_VIEW_BINARY: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for binary output"); 5599371c9d4SSatish Balay case DATABUCKET_VIEW_HDF5: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for HDF5 output"); 560521f74f9SMatthew G. Knepley default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown viewer method requested"); 56152849c42SDave May } 5627fbf63aeSDave May PetscFunctionReturn(0); 56352849c42SDave May } 56452849c42SDave May 5659371c9d4SSatish Balay PetscErrorCode DMSwarmDataBucketView(MPI_Comm comm, DMSwarmDataBucket db, const char filename[], DMSwarmDataBucketViewType type) { 566d7d19db6SBarry Smith PetscMPIInt size; 56752849c42SDave May 568521f74f9SMatthew G. Knepley PetscFunctionBegin; 5699566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 570d7d19db6SBarry Smith if (size == 1) { 5719566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketView_Seq(comm, db, filename, type)); 57252849c42SDave May } else { 5739566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketView_MPI(comm, db, filename, type)); 57452849c42SDave May } 5757fbf63aeSDave May PetscFunctionReturn(0); 57652849c42SDave May } 57752849c42SDave May 5789371c9d4SSatish Balay PetscErrorCode DMSwarmDataBucketDuplicateFields(DMSwarmDataBucket dbA, DMSwarmDataBucket *dbB) { 57977048351SPatrick Sanan DMSwarmDataBucket db2; 5805c18a9d6SDave May PetscInt f; 58152849c42SDave May 582521f74f9SMatthew G. Knepley PetscFunctionBegin; 5839566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketCreate(&db2)); 58452849c42SDave May /* copy contents from dbA into db2 */ 585521f74f9SMatthew G. Knepley for (f = 0; f < dbA->nfields; ++f) { 58677048351SPatrick Sanan DMSwarmDataField field; 58752849c42SDave May size_t atomic_size; 58852849c42SDave May char *name; 58952849c42SDave May 59052849c42SDave May field = dbA->field[f]; 59152849c42SDave May atomic_size = field->atomic_size; 59252849c42SDave May name = field->name; 5939566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketRegisterField(db2, "DMSwarmDataBucketDuplicateFields", name, atomic_size, NULL)); 59452849c42SDave May } 5959566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketFinalize(db2)); 5969566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketSetInitialSizes(db2, 0, 1000)); 59752849c42SDave May *dbB = db2; 5987fbf63aeSDave May PetscFunctionReturn(0); 59952849c42SDave May } 60052849c42SDave May 60152849c42SDave May /* 60252849c42SDave May Insert points from db2 into db1 60352849c42SDave May db1 <<== db2 60452849c42SDave May */ 6059371c9d4SSatish Balay PetscErrorCode DMSwarmDataBucketInsertValues(DMSwarmDataBucket db1, DMSwarmDataBucket db2) { 6065c18a9d6SDave May PetscInt n_mp_points1, n_mp_points2; 6075c18a9d6SDave May PetscInt n_mp_points1_new, p; 60852849c42SDave May 609521f74f9SMatthew G. Knepley PetscFunctionBegin; 6109566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(db1, &n_mp_points1, NULL, NULL)); 6119566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetSizes(db2, &n_mp_points2, NULL, NULL)); 61252849c42SDave May n_mp_points1_new = n_mp_points1 + n_mp_points2; 6139566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketSetSizes(db1, n_mp_points1_new, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT)); 614521f74f9SMatthew G. Knepley for (p = 0; p < n_mp_points2; ++p) { 615521f74f9SMatthew G. Knepley /* db1 <<== db2 */ 6169566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketCopyPoint(db2, p, db1, (n_mp_points1 + p))); 61752849c42SDave May } 6187fbf63aeSDave May PetscFunctionReturn(0); 61952849c42SDave May } 62052849c42SDave May 62152849c42SDave May /* helpers for parallel send/recv */ 6229371c9d4SSatish Balay PetscErrorCode DMSwarmDataBucketCreatePackedArray(DMSwarmDataBucket db, size_t *bytes, void **buf) { 6235c18a9d6SDave May PetscInt f; 62452849c42SDave May size_t sizeof_marker_contents; 62552849c42SDave May void *buffer; 62652849c42SDave May 627521f74f9SMatthew G. Knepley PetscFunctionBegin; 62852849c42SDave May sizeof_marker_contents = 0; 629521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 63077048351SPatrick Sanan DMSwarmDataField df = db->field[f]; 63152849c42SDave May sizeof_marker_contents += df->atomic_size; 63252849c42SDave May } 6339566063dSJacob Faibussowitsch PetscCall(PetscMalloc(sizeof_marker_contents, &buffer)); 6349566063dSJacob Faibussowitsch PetscCall(PetscMemzero(buffer, sizeof_marker_contents)); 63552849c42SDave May if (bytes) { *bytes = sizeof_marker_contents; } 63652849c42SDave May if (buf) { *buf = buffer; } 6377fbf63aeSDave May PetscFunctionReturn(0); 63852849c42SDave May } 63952849c42SDave May 6409371c9d4SSatish Balay PetscErrorCode DMSwarmDataBucketDestroyPackedArray(DMSwarmDataBucket db, void **buf) { 641521f74f9SMatthew G. Knepley PetscFunctionBegin; 64252849c42SDave May if (buf) { 6439566063dSJacob Faibussowitsch PetscCall(PetscFree(*buf)); 64452849c42SDave May *buf = NULL; 64552849c42SDave May } 6467fbf63aeSDave May PetscFunctionReturn(0); 64752849c42SDave May } 64852849c42SDave May 6499371c9d4SSatish Balay PetscErrorCode DMSwarmDataBucketFillPackedArray(DMSwarmDataBucket db, const PetscInt index, void *buf) { 6505c18a9d6SDave May PetscInt f; 65152849c42SDave May void *data, *data_p; 65252849c42SDave May size_t asize, offset; 65352849c42SDave May 654521f74f9SMatthew G. Knepley PetscFunctionBegin; 65552849c42SDave May offset = 0; 656521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 65777048351SPatrick Sanan DMSwarmDataField df = db->field[f]; 65852849c42SDave May 65952849c42SDave May asize = df->atomic_size; 66052849c42SDave May data = (void *)(df->data); 66152849c42SDave May data_p = (void *)((char *)data + index * asize); 6629566063dSJacob Faibussowitsch PetscCall(PetscMemcpy((void *)((char *)buf + offset), data_p, asize)); 66352849c42SDave May offset = offset + asize; 66452849c42SDave May } 6657fbf63aeSDave May PetscFunctionReturn(0); 66652849c42SDave May } 66752849c42SDave May 6689371c9d4SSatish Balay PetscErrorCode DMSwarmDataBucketInsertPackedArray(DMSwarmDataBucket db, const PetscInt idx, void *data) { 6695c18a9d6SDave May PetscInt f; 67052849c42SDave May void *data_p; 67152849c42SDave May size_t offset; 67252849c42SDave May 673521f74f9SMatthew G. Knepley PetscFunctionBegin; 67452849c42SDave May offset = 0; 675521f74f9SMatthew G. Knepley for (f = 0; f < db->nfields; ++f) { 67677048351SPatrick Sanan DMSwarmDataField df = db->field[f]; 67752849c42SDave May 67852849c42SDave May data_p = (void *)((char *)data + offset); 6799566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldInsertPoint(df, idx, (void *)data_p)); 68052849c42SDave May offset = offset + df->atomic_size; 68152849c42SDave May } 6827fbf63aeSDave May PetscFunctionReturn(0); 68352849c42SDave May } 684