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