1 #include "../src/dm/impls/swarm/data_bucket.h" 2 3 /* string helpers */ 4 PetscErrorCode DMSwarmDataFieldStringInList(const char name[], const PetscInt N, const DMSwarmDataField gfield[], PetscBool *val) 5 { 6 PetscInt i; 7 8 PetscFunctionBegin; 9 *val = PETSC_FALSE; 10 for (i = 0; i < N; ++i) { 11 PetscBool flg; 12 PetscCall(PetscStrcmp(name, gfield[i]->name, &flg)); 13 if (flg) { 14 *val = PETSC_TRUE; 15 PetscFunctionReturn(PETSC_SUCCESS); 16 } 17 } 18 PetscFunctionReturn(PETSC_SUCCESS); 19 } 20 21 PetscErrorCode DMSwarmDataFieldStringFindInList(const char name[], const PetscInt N, const DMSwarmDataField gfield[], PetscInt *index) 22 { 23 PetscInt i; 24 25 PetscFunctionBegin; 26 *index = -1; 27 for (i = 0; i < N; ++i) { 28 PetscBool flg; 29 PetscCall(PetscStrcmp(name, gfield[i]->name, &flg)); 30 if (flg) { 31 *index = i; 32 PetscFunctionReturn(PETSC_SUCCESS); 33 } 34 } 35 PetscFunctionReturn(PETSC_SUCCESS); 36 } 37 38 PetscErrorCode DMSwarmDataFieldCreate(const char registration_function[], const char name[], const size_t size, const PetscInt L, DMSwarmDataField *DF) 39 { 40 DMSwarmDataField df; 41 42 PetscFunctionBegin; 43 PetscCall(PetscNew(&df)); 44 PetscCall(PetscStrallocpy(registration_function, &df->registration_function)); 45 PetscCall(PetscStrallocpy(name, &df->name)); 46 df->atomic_size = size; 47 df->L = L; 48 df->bs = 1; 49 /* allocate something so we don't have to reallocate */ 50 PetscCall(PetscMalloc(size * L, &df->data)); 51 PetscCall(PetscMemzero(df->data, size * L)); 52 *DF = df; 53 PetscFunctionReturn(PETSC_SUCCESS); 54 } 55 56 PetscErrorCode DMSwarmDataFieldDestroy(DMSwarmDataField *DF) 57 { 58 DMSwarmDataField df = *DF; 59 60 PetscFunctionBegin; 61 PetscCall(PetscFree(df->registration_function)); 62 PetscCall(PetscFree(df->name)); 63 PetscCall(PetscFree(df->data)); 64 PetscCall(PetscFree(df)); 65 *DF = NULL; 66 PetscFunctionReturn(PETSC_SUCCESS); 67 } 68 69 /* data bucket */ 70 PetscErrorCode DMSwarmDataBucketCreate(DMSwarmDataBucket *DB) 71 { 72 DMSwarmDataBucket db; 73 74 PetscFunctionBegin; 75 PetscCall(PetscNew(&db)); 76 77 db->finalised = PETSC_FALSE; 78 /* create empty spaces for fields */ 79 db->L = -1; 80 db->buffer = 1; 81 db->allocated = 1; 82 db->nfields = 0; 83 PetscCall(PetscMalloc1(1, &db->field)); 84 *DB = db; 85 PetscFunctionReturn(PETSC_SUCCESS); 86 } 87 88 PetscErrorCode DMSwarmDataBucketDestroy(DMSwarmDataBucket *DB) 89 { 90 DMSwarmDataBucket db = *DB; 91 PetscInt f; 92 93 PetscFunctionBegin; 94 /* release fields */ 95 for (f = 0; f < db->nfields; ++f) PetscCall(DMSwarmDataFieldDestroy(&db->field[f])); 96 /* this will catch the initially allocated objects in the event that no fields are registered */ 97 if (db->field != NULL) PetscCall(PetscFree(db->field)); 98 PetscCall(PetscFree(db)); 99 *DB = NULL; 100 PetscFunctionReturn(PETSC_SUCCESS); 101 } 102 103 PetscErrorCode DMSwarmDataBucketQueryForActiveFields(DMSwarmDataBucket db, PetscBool *any_active_fields) 104 { 105 PetscInt f; 106 107 PetscFunctionBegin; 108 *any_active_fields = PETSC_FALSE; 109 for (f = 0; f < db->nfields; ++f) { 110 if (db->field[f]->active) { 111 *any_active_fields = PETSC_TRUE; 112 PetscFunctionReturn(PETSC_SUCCESS); 113 } 114 } 115 PetscFunctionReturn(PETSC_SUCCESS); 116 } 117 118 PetscErrorCode DMSwarmDataBucketRegisterField(DMSwarmDataBucket db, const char registration_function[], const char field_name[], size_t atomic_size, DMSwarmDataField *_gfield) 119 { 120 PetscBool val; 121 DMSwarmDataField fp; 122 123 PetscFunctionBegin; 124 /* check we haven't finalised the registration of fields */ 125 /* 126 if (db->finalised==PETSC_TRUE) { 127 printf("ERROR: DMSwarmDataBucketFinalize() has been called. Cannot register more fields\n"); 128 ERROR(); 129 } 130 */ 131 /* check for repeated name */ 132 PetscCall(DMSwarmDataFieldStringInList(field_name, db->nfields, (const DMSwarmDataField *)db->field, &val)); 133 PetscCheck(val != PETSC_TRUE, PETSC_COMM_SELF, PETSC_ERR_USER, "Field %s already exists. Cannot add same field twice", field_name); 134 /* create new space for data */ 135 PetscCall(PetscRealloc(sizeof(DMSwarmDataField) * (db->nfields + 1), &db->field)); 136 /* add field */ 137 PetscCall(DMSwarmDataFieldCreate(registration_function, field_name, atomic_size, db->allocated, &fp)); 138 db->field[db->nfields] = fp; 139 db->nfields++; 140 if (_gfield != NULL) *_gfield = fp; 141 PetscFunctionReturn(PETSC_SUCCESS); 142 } 143 144 /* 145 #define DMSwarmDataBucketRegisterField(db,name,size,k) {\ 146 char *location;\ 147 asprintf(&location,"Registered by %s() at line %d within file %s", __FUNCTION__, __LINE__, __FILE__);\ 148 _DMSwarmDataBucketRegisterField( (db), location, (name), (size), (k));\ 149 ierr = PetscFree(location);\ 150 } 151 */ 152 153 PetscErrorCode DMSwarmDataBucketGetDMSwarmDataFieldIdByName(DMSwarmDataBucket db, const char name[], PetscInt *idx) 154 { 155 PetscBool found; 156 157 PetscFunctionBegin; 158 *idx = -1; 159 PetscCall(DMSwarmDataFieldStringInList(name, db->nfields, (const DMSwarmDataField *)db->field, &found)); 160 PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot find DMSwarmDataField with name %s", name); 161 PetscCall(DMSwarmDataFieldStringFindInList(name, db->nfields, (const DMSwarmDataField *)db->field, idx)); 162 PetscFunctionReturn(PETSC_SUCCESS); 163 } 164 165 PetscErrorCode DMSwarmDataBucketGetDMSwarmDataFieldByName(DMSwarmDataBucket db, const char name[], DMSwarmDataField *gfield) 166 { 167 PetscInt idx; 168 PetscBool found; 169 170 PetscFunctionBegin; 171 PetscCall(DMSwarmDataFieldStringInList(name, db->nfields, (const DMSwarmDataField *)db->field, &found)); 172 PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot find DMSwarmDataField with name %s", name); 173 PetscCall(DMSwarmDataFieldStringFindInList(name, db->nfields, (const DMSwarmDataField *)db->field, &idx)); 174 *gfield = db->field[idx]; 175 PetscFunctionReturn(PETSC_SUCCESS); 176 } 177 178 PetscErrorCode DMSwarmDataBucketQueryDMSwarmDataFieldByName(DMSwarmDataBucket db, const char name[], PetscBool *found) 179 { 180 PetscFunctionBegin; 181 *found = PETSC_FALSE; 182 PetscCall(DMSwarmDataFieldStringInList(name, db->nfields, (const DMSwarmDataField *)db->field, found)); 183 PetscFunctionReturn(PETSC_SUCCESS); 184 } 185 186 PetscErrorCode DMSwarmDataBucketFinalize(DMSwarmDataBucket db) 187 { 188 PetscFunctionBegin; 189 db->finalised = PETSC_TRUE; 190 PetscFunctionReturn(PETSC_SUCCESS); 191 } 192 193 PetscErrorCode DMSwarmDataFieldGetNumEntries(DMSwarmDataField df, PetscInt *sum) 194 { 195 PetscFunctionBegin; 196 *sum = df->L; 197 PetscFunctionReturn(PETSC_SUCCESS); 198 } 199 200 PetscErrorCode DMSwarmDataFieldSetBlockSize(DMSwarmDataField df, PetscInt blocksize) 201 { 202 PetscFunctionBegin; 203 df->bs = blocksize; 204 PetscFunctionReturn(PETSC_SUCCESS); 205 } 206 207 PetscErrorCode DMSwarmDataFieldSetSize(DMSwarmDataField df, const PetscInt new_L) 208 { 209 PetscFunctionBegin; 210 PetscCheck(new_L >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot set size of DMSwarmDataField to be < 0"); 211 if (new_L == df->L) PetscFunctionReturn(PETSC_SUCCESS); 212 if (new_L > df->L) { 213 PetscCall(PetscRealloc(df->atomic_size * (new_L), &df->data)); 214 /* init new contents */ 215 PetscCall(PetscMemzero((((char *)df->data) + df->L * df->atomic_size), (new_L - df->L) * df->atomic_size)); 216 } else { 217 /* reallocate pointer list, add +1 in case new_L = 0 */ 218 PetscCall(PetscRealloc(df->atomic_size * (new_L + 1), &df->data)); 219 } 220 df->L = new_L; 221 PetscFunctionReturn(PETSC_SUCCESS); 222 } 223 224 PetscErrorCode DMSwarmDataFieldZeroBlock(DMSwarmDataField df, const PetscInt start, const PetscInt end) 225 { 226 PetscFunctionBegin; 227 PetscCheck(start <= end, PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot zero a block of entries if start(%" PetscInt_FMT ") > end(%" PetscInt_FMT ")", start, end); 228 PetscCheck(start >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot zero a block of entries if start(%" PetscInt_FMT ") < 0", start); 229 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); 230 PetscCall(PetscMemzero((((char *)df->data) + start * df->atomic_size), (end - start) * df->atomic_size)); 231 PetscFunctionReturn(PETSC_SUCCESS); 232 } 233 234 /* 235 A negative buffer value will simply be ignored and the old buffer value will be used. 236 */ 237 PetscErrorCode DMSwarmDataBucketSetSizes(DMSwarmDataBucket db, const PetscInt L, const PetscInt buffer) 238 { 239 PetscInt current_allocated, new_used, new_unused, new_buffer, new_allocated, f; 240 PetscBool any_active_fields; 241 242 PetscFunctionBegin; 243 PetscCheck(db->finalised != PETSC_FALSE, PETSC_COMM_SELF, PETSC_ERR_USER, "You must call DMSwarmDataBucketFinalize() before DMSwarmDataBucketSetSizes()"); 244 PetscCall(DMSwarmDataBucketQueryForActiveFields(db, &any_active_fields)); 245 PetscCheck(!any_active_fields, PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot safely re-size as at least one DMSwarmDataField is currently being accessed"); 246 247 current_allocated = db->allocated; 248 new_used = L; 249 new_unused = current_allocated - new_used; 250 new_buffer = db->buffer; 251 if (buffer >= 0) { /* update the buffer value */ 252 new_buffer = buffer; 253 } 254 new_allocated = new_used + new_buffer; 255 /* action */ 256 if (new_allocated > current_allocated) { 257 /* increase size to new_used + new_buffer */ 258 for (f = 0; f < db->nfields; f++) PetscCall(DMSwarmDataFieldSetSize(db->field[f], new_allocated)); 259 db->L = new_used; 260 db->buffer = new_buffer; 261 db->allocated = new_used + new_buffer; 262 } else { 263 if (new_unused > 2 * new_buffer) { 264 /* shrink array to new_used + new_buffer */ 265 for (f = 0; f < db->nfields; ++f) PetscCall(DMSwarmDataFieldSetSize(db->field[f], new_allocated)); 266 db->L = new_used; 267 db->buffer = new_buffer; 268 db->allocated = new_used + new_buffer; 269 } else { 270 db->L = new_used; 271 db->buffer = new_buffer; 272 } 273 } 274 /* zero all entries from db->L to db->allocated */ 275 for (f = 0; f < db->nfields; ++f) { 276 DMSwarmDataField field = db->field[f]; 277 PetscCall(DMSwarmDataFieldZeroBlock(field, db->L, db->allocated)); 278 } 279 PetscFunctionReturn(PETSC_SUCCESS); 280 } 281 282 PetscErrorCode DMSwarmDataBucketSetInitialSizes(DMSwarmDataBucket db, const PetscInt L, const PetscInt buffer) 283 { 284 PetscInt f; 285 286 PetscFunctionBegin; 287 PetscCall(DMSwarmDataBucketSetSizes(db, L, buffer)); 288 for (f = 0; f < db->nfields; ++f) { 289 DMSwarmDataField field = db->field[f]; 290 PetscCall(DMSwarmDataFieldZeroBlock(field, 0, db->allocated)); 291 } 292 PetscFunctionReturn(PETSC_SUCCESS); 293 } 294 295 PetscErrorCode DMSwarmDataBucketGetSizes(DMSwarmDataBucket db, PetscInt *L, PetscInt *buffer, PetscInt *allocated) 296 { 297 PetscFunctionBegin; 298 if (L) *L = db->L; 299 if (buffer) *buffer = db->buffer; 300 if (allocated) *allocated = db->allocated; 301 PetscFunctionReturn(PETSC_SUCCESS); 302 } 303 304 PetscErrorCode DMSwarmDataBucketGetGlobalSizes(MPI_Comm comm, DMSwarmDataBucket db, PetscInt *L, PetscInt *buffer, PetscInt *allocated) 305 { 306 PetscFunctionBegin; 307 if (L) PetscCallMPI(MPI_Allreduce(&db->L, L, 1, MPIU_INT, MPI_SUM, comm)); 308 if (buffer) PetscCallMPI(MPI_Allreduce(&db->buffer, buffer, 1, MPIU_INT, MPI_SUM, comm)); 309 if (allocated) PetscCallMPI(MPI_Allreduce(&db->allocated, allocated, 1, MPIU_INT, MPI_SUM, comm)); 310 PetscFunctionReturn(PETSC_SUCCESS); 311 } 312 313 PetscErrorCode DMSwarmDataBucketGetDMSwarmDataFields(DMSwarmDataBucket db, PetscInt *L, DMSwarmDataField *fields[]) 314 { 315 PetscFunctionBegin; 316 if (L) *L = db->nfields; 317 if (fields) *fields = db->field; 318 PetscFunctionReturn(PETSC_SUCCESS); 319 } 320 321 PetscErrorCode DMSwarmDataFieldGetAccess(const DMSwarmDataField gfield) 322 { 323 PetscFunctionBegin; 324 PetscCheck(!gfield->active, PETSC_COMM_SELF, PETSC_ERR_USER, "Field \"%s\" is already active. You must call DMSwarmDataFieldRestoreAccess()", gfield->name); 325 gfield->active = PETSC_TRUE; 326 PetscFunctionReturn(PETSC_SUCCESS); 327 } 328 329 PetscErrorCode DMSwarmDataFieldAccessPoint(const DMSwarmDataField gfield, const PetscInt pid, void **ctx_p) 330 { 331 PetscFunctionBegin; 332 *ctx_p = NULL; 333 #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD) 334 /* debug mode */ 335 /* check point is valid */ 336 PetscCheck(pid >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be >= 0"); 337 PetscCheck(pid < gfield->L, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be < %" PetscInt_FMT, gfield->L); 338 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); 339 #endif 340 *ctx_p = DMSWARM_DATAFIELD_point_access(gfield->data, pid, gfield->atomic_size); 341 PetscFunctionReturn(PETSC_SUCCESS); 342 } 343 344 PetscErrorCode DMSwarmDataFieldAccessPointOffset(const DMSwarmDataField gfield, const size_t offset, const PetscInt pid, void **ctx_p) 345 { 346 PetscFunctionBegin; 347 #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD) 348 /* debug mode */ 349 /* check point is valid */ 350 /* PetscCheck(offset >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be >= 0");*/ 351 /* Note compiler realizes this can never happen with an unsigned PetscInt */ 352 PetscCheck(offset < gfield->atomic_size, PETSC_COMM_SELF, PETSC_ERR_USER, "offset must be < %zu", gfield->atomic_size); 353 /* check point is valid */ 354 PetscCheck(pid >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be >= 0"); 355 PetscCheck(pid < gfield->L, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be < %" PetscInt_FMT, gfield->L); 356 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); 357 #endif 358 *ctx_p = DMSWARM_DATAFIELD_point_access_offset(gfield->data, pid, gfield->atomic_size, offset); 359 PetscFunctionReturn(PETSC_SUCCESS); 360 } 361 362 PetscErrorCode DMSwarmDataFieldRestoreAccess(DMSwarmDataField gfield) 363 { 364 PetscFunctionBegin; 365 PetscCheck(gfield->active != PETSC_FALSE, PETSC_COMM_SELF, PETSC_ERR_USER, "Field \"%s\" is not active. You must call DMSwarmDataFieldGetAccess()", gfield->name); 366 gfield->active = PETSC_FALSE; 367 PetscFunctionReturn(PETSC_SUCCESS); 368 } 369 370 PetscErrorCode DMSwarmDataFieldVerifyAccess(const DMSwarmDataField gfield, const size_t size) 371 { 372 PetscFunctionBegin; 373 #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD) 374 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); 375 #endif 376 PetscFunctionReturn(PETSC_SUCCESS); 377 } 378 379 PetscErrorCode DMSwarmDataFieldGetAtomicSize(const DMSwarmDataField gfield, size_t *size) 380 { 381 PetscFunctionBegin; 382 if (size) *size = gfield->atomic_size; 383 PetscFunctionReturn(PETSC_SUCCESS); 384 } 385 386 PetscErrorCode DMSwarmDataFieldGetEntries(const DMSwarmDataField gfield, void **data) 387 { 388 PetscFunctionBegin; 389 if (data) *data = gfield->data; 390 PetscFunctionReturn(PETSC_SUCCESS); 391 } 392 393 PetscErrorCode DMSwarmDataFieldRestoreEntries(const DMSwarmDataField gfield, void **data) 394 { 395 PetscFunctionBegin; 396 if (data) *data = NULL; 397 PetscFunctionReturn(PETSC_SUCCESS); 398 } 399 400 /* y = x */ 401 PetscErrorCode DMSwarmDataBucketCopyPoint(const DMSwarmDataBucket xb, const PetscInt pid_x, const DMSwarmDataBucket yb, const PetscInt pid_y) 402 { 403 PetscInt f; 404 405 PetscFunctionBegin; 406 for (f = 0; f < xb->nfields; ++f) { 407 void *dest; 408 void *src; 409 410 PetscCall(DMSwarmDataFieldGetAccess(xb->field[f])); 411 if (xb != yb) PetscCall(DMSwarmDataFieldGetAccess(yb->field[f])); 412 PetscCall(DMSwarmDataFieldAccessPoint(xb->field[f], pid_x, &src)); 413 PetscCall(DMSwarmDataFieldAccessPoint(yb->field[f], pid_y, &dest)); 414 PetscCall(PetscMemcpy(dest, src, xb->field[f]->atomic_size)); 415 PetscCall(DMSwarmDataFieldRestoreAccess(xb->field[f])); 416 if (xb != yb) PetscCall(DMSwarmDataFieldRestoreAccess(yb->field[f])); 417 } 418 PetscFunctionReturn(PETSC_SUCCESS); 419 } 420 421 PetscErrorCode DMSwarmDataBucketCreateFromSubset(DMSwarmDataBucket DBIn, const PetscInt N, const PetscInt list[], DMSwarmDataBucket *DB) 422 { 423 PetscInt nfields; 424 DMSwarmDataField *fields; 425 PetscInt f, L, buffer, allocated, p; 426 427 PetscFunctionBegin; 428 PetscCall(DMSwarmDataBucketCreate(DB)); 429 /* copy contents of DBIn */ 430 PetscCall(DMSwarmDataBucketGetDMSwarmDataFields(DBIn, &nfields, &fields)); 431 PetscCall(DMSwarmDataBucketGetSizes(DBIn, &L, &buffer, &allocated)); 432 for (f = 0; f < nfields; ++f) PetscCall(DMSwarmDataBucketRegisterField(*DB, "DMSwarmDataBucketCreateFromSubset", fields[f]->name, fields[f]->atomic_size, NULL)); 433 PetscCall(DMSwarmDataBucketFinalize(*DB)); 434 PetscCall(DMSwarmDataBucketSetSizes(*DB, L, buffer)); 435 /* now copy the desired guys from DBIn => DB */ 436 for (p = 0; p < N; ++p) PetscCall(DMSwarmDataBucketCopyPoint(DBIn, list[p], *DB, list[p])); 437 PetscFunctionReturn(PETSC_SUCCESS); 438 } 439 440 /* insert into an existing location */ 441 PetscErrorCode DMSwarmDataFieldInsertPoint(const DMSwarmDataField field, const PetscInt index, const void *ctx) 442 { 443 PetscFunctionBegin; 444 #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD) 445 /* check point is valid */ 446 PetscCheck(index >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be >= 0"); 447 PetscCheck(index < field->L, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be < %" PetscInt_FMT, field->L); 448 #endif 449 PetscCall(PetscMemcpy(DMSWARM_DATAFIELD_point_access(field->data, index, field->atomic_size), ctx, field->atomic_size)); 450 PetscFunctionReturn(PETSC_SUCCESS); 451 } 452 453 /* remove data at index - replace with last point */ 454 PetscErrorCode DMSwarmDataBucketRemovePointAtIndex(const DMSwarmDataBucket db, const PetscInt index) 455 { 456 PetscInt f; 457 PetscBool any_active_fields; 458 459 PetscFunctionBegin; 460 #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD) 461 /* check point is valid */ 462 PetscCheck(index >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be >= 0"); 463 PetscCheck(index < db->allocated, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be < %" PetscInt_FMT, db->L + db->buffer); 464 #endif 465 PetscCall(DMSwarmDataBucketQueryForActiveFields(db, &any_active_fields)); 466 PetscCheck(!any_active_fields, PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot safely remove point as at least one DMSwarmDataField is currently being accessed"); 467 if (index >= db->L) { /* this point is not in the list - no need to error, but I will anyway */ 468 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); 469 } 470 if (index != db->L - 1) { /* not last point in list */ 471 for (f = 0; f < db->nfields; ++f) { 472 DMSwarmDataField field = db->field[f]; 473 474 /* copy then remove */ 475 PetscCall(DMSwarmDataFieldCopyPoint(db->L - 1, field, index, field)); 476 /* DMSwarmDataFieldZeroPoint(field,index); */ 477 } 478 } 479 /* decrement size */ 480 /* this will zero out an crap at the end of the list */ 481 PetscCall(DMSwarmDataBucketRemovePoint(db)); 482 PetscFunctionReturn(PETSC_SUCCESS); 483 } 484 485 /* copy x into y */ 486 PetscErrorCode DMSwarmDataFieldCopyPoint(const PetscInt pid_x, const DMSwarmDataField field_x, const PetscInt pid_y, const DMSwarmDataField field_y) 487 { 488 PetscFunctionBegin; 489 #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD) 490 /* check point is valid */ 491 PetscCheck(pid_x >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "(IN) index must be >= 0"); 492 PetscCheck(pid_x < field_x->L, PETSC_COMM_SELF, PETSC_ERR_USER, "(IN) index must be < %" PetscInt_FMT, field_x->L); 493 PetscCheck(pid_y >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "(OUT) index must be >= 0"); 494 PetscCheck(pid_y < field_y->L, PETSC_COMM_SELF, PETSC_ERR_USER, "(OUT) index must be < %" PetscInt_FMT, field_y->L); 495 PetscCheck(field_y->atomic_size == field_x->atomic_size, PETSC_COMM_SELF, PETSC_ERR_USER, "atomic size must match"); 496 #endif 497 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)); 498 PetscFunctionReturn(PETSC_SUCCESS); 499 } 500 501 /* zero only the datafield at this point */ 502 PetscErrorCode DMSwarmDataFieldZeroPoint(const DMSwarmDataField field, const PetscInt index) 503 { 504 PetscFunctionBegin; 505 #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD) 506 /* check point is valid */ 507 PetscCheck(index >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be >= 0"); 508 PetscCheck(index < field->L, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be < %" PetscInt_FMT, field->L); 509 #endif 510 PetscCall(PetscMemzero(DMSWARM_DATAFIELD_point_access(field->data, index, field->atomic_size), field->atomic_size)); 511 PetscFunctionReturn(PETSC_SUCCESS); 512 } 513 514 /* zero ALL data for this point */ 515 PetscErrorCode DMSwarmDataBucketZeroPoint(const DMSwarmDataBucket db, const PetscInt index) 516 { 517 PetscInt f; 518 519 PetscFunctionBegin; 520 /* check point is valid */ 521 PetscCheck(index >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be >= 0"); 522 PetscCheck(index < db->allocated, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be < %" PetscInt_FMT, db->allocated); 523 for (f = 0; f < db->nfields; ++f) { 524 DMSwarmDataField field = db->field[f]; 525 PetscCall(DMSwarmDataFieldZeroPoint(field, index)); 526 } 527 PetscFunctionReturn(PETSC_SUCCESS); 528 } 529 530 /* increment */ 531 PetscErrorCode DMSwarmDataBucketAddPoint(DMSwarmDataBucket db) 532 { 533 PetscFunctionBegin; 534 PetscCall(DMSwarmDataBucketSetSizes(db, db->L + 1, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT)); 535 PetscFunctionReturn(PETSC_SUCCESS); 536 } 537 538 /* decrement */ 539 PetscErrorCode DMSwarmDataBucketRemovePoint(DMSwarmDataBucket db) 540 { 541 PetscFunctionBegin; 542 PetscCall(DMSwarmDataBucketSetSizes(db, db->L - 1, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT)); 543 PetscFunctionReturn(PETSC_SUCCESS); 544 } 545 546 /* Should be redone to user PetscViewer */ 547 PetscErrorCode DMSwarmDataBucketView_stdout(MPI_Comm comm, DMSwarmDataBucket db) 548 { 549 PetscInt f; 550 double memory_usage_total, memory_usage_total_local = 0.0; 551 552 PetscFunctionBegin; 553 PetscCall(PetscPrintf(comm, "DMSwarmDataBucketView: \n")); 554 PetscCall(PetscPrintf(comm, " L = %" PetscInt_FMT " \n", db->L)); 555 PetscCall(PetscPrintf(comm, " buffer = %" PetscInt_FMT " \n", db->buffer)); 556 PetscCall(PetscPrintf(comm, " allocated = %" PetscInt_FMT " \n", db->allocated)); 557 PetscCall(PetscPrintf(comm, " nfields registered = %" PetscInt_FMT " \n", db->nfields)); 558 559 for (f = 0; f < db->nfields; ++f) { 560 double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6; 561 memory_usage_total_local += memory_usage_f; 562 } 563 PetscCallMPI(MPI_Allreduce(&memory_usage_total_local, &memory_usage_total, 1, MPI_DOUBLE, MPI_SUM, comm)); 564 565 for (f = 0; f < db->nfields; ++f) { 566 double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6; 567 PetscCall(PetscPrintf(comm, " [%3" PetscInt_FMT "] %15s : Mem. usage = %1.2e (MB) [rank0]\n", f, db->field[f]->name, memory_usage_f)); 568 PetscCall(PetscPrintf(comm, " blocksize = %" PetscInt_FMT " \n", db->field[f]->bs)); 569 if (db->field[f]->bs != 1) { 570 PetscCall(PetscPrintf(comm, " atomic size = %zu [full block, bs=%" PetscInt_FMT "]\n", db->field[f]->atomic_size, db->field[f]->bs)); 571 PetscCall(PetscPrintf(comm, " atomic size/item = %zu \n", (size_t)(db->field[f]->atomic_size / db->field[f]->bs))); 572 } else { 573 PetscCall(PetscPrintf(comm, " atomic size = %zu \n", db->field[f]->atomic_size)); 574 } 575 } 576 PetscCall(PetscPrintf(comm, " Total mem. usage = %1.2e (MB) (collective)\n", memory_usage_total)); 577 PetscFunctionReturn(PETSC_SUCCESS); 578 } 579 580 PetscErrorCode DMSwarmDataBucketView_Seq(MPI_Comm comm, DMSwarmDataBucket db, const char filename[], DMSwarmDataBucketViewType type) 581 { 582 PetscFunctionBegin; 583 switch (type) { 584 case DATABUCKET_VIEW_STDOUT: 585 PetscCall(DMSwarmDataBucketView_stdout(PETSC_COMM_SELF, db)); 586 break; 587 case DATABUCKET_VIEW_ASCII: 588 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for ascii output"); 589 case DATABUCKET_VIEW_BINARY: 590 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for binary output"); 591 case DATABUCKET_VIEW_HDF5: 592 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for HDF5 output"); 593 default: 594 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown viewer method requested"); 595 } 596 PetscFunctionReturn(PETSC_SUCCESS); 597 } 598 599 PetscErrorCode DMSwarmDataBucketView_MPI(MPI_Comm comm, DMSwarmDataBucket db, const char filename[], DMSwarmDataBucketViewType type) 600 { 601 PetscFunctionBegin; 602 switch (type) { 603 case DATABUCKET_VIEW_STDOUT: 604 PetscCall(DMSwarmDataBucketView_stdout(comm, db)); 605 break; 606 case DATABUCKET_VIEW_ASCII: 607 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for ascii output"); 608 case DATABUCKET_VIEW_BINARY: 609 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for binary output"); 610 case DATABUCKET_VIEW_HDF5: 611 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for HDF5 output"); 612 default: 613 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown viewer method requested"); 614 } 615 PetscFunctionReturn(PETSC_SUCCESS); 616 } 617 618 PetscErrorCode DMSwarmDataBucketView(MPI_Comm comm, DMSwarmDataBucket db, const char filename[], DMSwarmDataBucketViewType type) 619 { 620 PetscMPIInt size; 621 622 PetscFunctionBegin; 623 PetscCallMPI(MPI_Comm_size(comm, &size)); 624 if (size == 1) { 625 PetscCall(DMSwarmDataBucketView_Seq(comm, db, filename, type)); 626 } else { 627 PetscCall(DMSwarmDataBucketView_MPI(comm, db, filename, type)); 628 } 629 PetscFunctionReturn(PETSC_SUCCESS); 630 } 631 632 PetscErrorCode DMSwarmDataBucketDuplicateFields(DMSwarmDataBucket dbA, DMSwarmDataBucket *dbB) 633 { 634 DMSwarmDataBucket db2; 635 PetscInt f; 636 637 PetscFunctionBegin; 638 PetscCall(DMSwarmDataBucketCreate(&db2)); 639 /* copy contents from dbA into db2 */ 640 for (f = 0; f < dbA->nfields; ++f) { 641 DMSwarmDataField field; 642 size_t atomic_size; 643 char *name; 644 645 field = dbA->field[f]; 646 atomic_size = field->atomic_size; 647 name = field->name; 648 PetscCall(DMSwarmDataBucketRegisterField(db2, "DMSwarmDataBucketDuplicateFields", name, atomic_size, NULL)); 649 } 650 PetscCall(DMSwarmDataBucketFinalize(db2)); 651 PetscCall(DMSwarmDataBucketSetInitialSizes(db2, 0, 1000)); 652 *dbB = db2; 653 PetscFunctionReturn(PETSC_SUCCESS); 654 } 655 656 /* 657 Insert points from db2 into db1 658 db1 <<== db2 659 */ 660 PetscErrorCode DMSwarmDataBucketInsertValues(DMSwarmDataBucket db1, DMSwarmDataBucket db2) 661 { 662 PetscInt n_mp_points1, n_mp_points2; 663 PetscInt n_mp_points1_new, p; 664 665 PetscFunctionBegin; 666 PetscCall(DMSwarmDataBucketGetSizes(db1, &n_mp_points1, NULL, NULL)); 667 PetscCall(DMSwarmDataBucketGetSizes(db2, &n_mp_points2, NULL, NULL)); 668 n_mp_points1_new = n_mp_points1 + n_mp_points2; 669 PetscCall(DMSwarmDataBucketSetSizes(db1, n_mp_points1_new, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT)); 670 for (p = 0; p < n_mp_points2; ++p) { 671 /* db1 <<== db2 */ 672 PetscCall(DMSwarmDataBucketCopyPoint(db2, p, db1, (n_mp_points1 + p))); 673 } 674 PetscFunctionReturn(PETSC_SUCCESS); 675 } 676 677 /* helpers for parallel send/recv */ 678 PetscErrorCode DMSwarmDataBucketCreatePackedArray(DMSwarmDataBucket db, size_t *bytes, void **buf) 679 { 680 PetscInt f; 681 size_t sizeof_marker_contents; 682 void *buffer; 683 684 PetscFunctionBegin; 685 sizeof_marker_contents = 0; 686 for (f = 0; f < db->nfields; ++f) { 687 DMSwarmDataField df = db->field[f]; 688 sizeof_marker_contents += df->atomic_size; 689 } 690 PetscCall(PetscMalloc(sizeof_marker_contents, &buffer)); 691 PetscCall(PetscMemzero(buffer, sizeof_marker_contents)); 692 if (bytes) *bytes = sizeof_marker_contents; 693 if (buf) *buf = buffer; 694 PetscFunctionReturn(PETSC_SUCCESS); 695 } 696 697 PetscErrorCode DMSwarmDataBucketDestroyPackedArray(DMSwarmDataBucket db, void **buf) 698 { 699 PetscFunctionBegin; 700 if (buf) { 701 PetscCall(PetscFree(*buf)); 702 *buf = NULL; 703 } 704 PetscFunctionReturn(PETSC_SUCCESS); 705 } 706 707 PetscErrorCode DMSwarmDataBucketFillPackedArray(DMSwarmDataBucket db, const PetscInt index, void *buf) 708 { 709 PetscInt f; 710 void *data, *data_p; 711 size_t asize, offset; 712 713 PetscFunctionBegin; 714 offset = 0; 715 for (f = 0; f < db->nfields; ++f) { 716 DMSwarmDataField df = db->field[f]; 717 718 asize = df->atomic_size; 719 data = (void *)(df->data); 720 data_p = (void *)((char *)data + index * asize); 721 PetscCall(PetscMemcpy((void *)((char *)buf + offset), data_p, asize)); 722 offset = offset + asize; 723 } 724 PetscFunctionReturn(PETSC_SUCCESS); 725 } 726 727 PetscErrorCode DMSwarmDataBucketInsertPackedArray(DMSwarmDataBucket db, const PetscInt idx, void *data) 728 { 729 PetscInt f; 730 void *data_p; 731 size_t offset; 732 733 PetscFunctionBegin; 734 offset = 0; 735 for (f = 0; f < db->nfields; ++f) { 736 DMSwarmDataField df = db->field[f]; 737 738 data_p = (void *)((char *)data + offset); 739 PetscCall(DMSwarmDataFieldInsertPoint(df, idx, (void *)data_p)); 740 offset = offset + df->atomic_size; 741 } 742 PetscFunctionReturn(PETSC_SUCCESS); 743 } 744