1 2 #include "data_bucket.h" 3 4 /* string helpers */ 5 void StringInList( const char name[], const PetscInt N, const DataField gfield[], PetscBool *val ) 6 { 7 PetscInt i; 8 9 *val = PETSC_FALSE; 10 for( i=0; i<N; i++ ) { 11 if( strcmp( name, gfield[i]->name ) == 0 ) { 12 *val = PETSC_TRUE; 13 return; 14 } 15 } 16 } 17 18 void StringFindInList( const char name[], const PetscInt N, const DataField gfield[], PetscInt *index ) 19 { 20 PetscInt i; 21 22 *index = -1; 23 for( i=0; i<N; i++ ) { 24 if( strcmp( name, gfield[i]->name ) == 0 ) { 25 *index = i; 26 return; 27 } 28 } 29 } 30 31 void DataFieldCreate( const char registeration_function[], const char name[], const size_t size, const PetscInt L, DataField *DF ) 32 { 33 DataField df; 34 35 df = malloc( sizeof(struct _p_DataField) ); 36 memset( df, 0, sizeof(struct _p_DataField) ); 37 38 39 asprintf( &df->registeration_function, "%s", registeration_function ); 40 asprintf( &df->name, "%s", name ); 41 df->atomic_size = size; 42 df->L = L; 43 44 df->data = malloc( size * L ); /* allocate something so we don't have to reallocate */ 45 memset( df->data, 0, size * L ); 46 47 *DF = df; 48 } 49 50 void DataFieldDestroy( DataField *DF ) 51 { 52 DataField df = *DF; 53 54 free( df->registeration_function ); 55 free( df->name ); 56 free( df->data ); 57 free(df); 58 59 *DF = NULL; 60 } 61 62 /* data bucket */ 63 void DataBucketCreate( DataBucket *DB ) 64 { 65 DataBucket db; 66 67 68 db = malloc( sizeof(struct _p_DataBucket) ); 69 memset( db, 0, sizeof(struct _p_DataBucket) ); 70 71 db->finalised = PETSC_FALSE; 72 73 /* create empty spaces for fields */ 74 db->L = 0; 75 db->buffer = 1; 76 db->allocated = 1; 77 78 db->nfields = 0; 79 db->field = malloc(sizeof(DataField)); 80 81 *DB = db; 82 } 83 84 void DataBucketDestroy( DataBucket *DB ) 85 { 86 DataBucket db = *DB; 87 PetscInt f; 88 89 /* release fields */ 90 for( f=0; f<db->nfields; f++ ) { 91 DataFieldDestroy(&db->field[f]); 92 } 93 94 /* this will catch the initially allocated objects in the event that no fields are registered */ 95 if(db->field!=NULL) { 96 free(db->field); 97 } 98 99 free(db); 100 101 *DB = NULL; 102 } 103 104 void _DataBucketRegisterField( 105 DataBucket db, 106 const char registeration_function[], 107 const char field_name[], 108 size_t atomic_size, DataField *_gfield ) 109 { 110 PetscBool val; 111 DataField *field,fp; 112 113 /* check we haven't finalised the registration of fields */ 114 /* 115 if(db->finalised==PETSC_TRUE) { 116 printf("ERROR: DataBucketFinalize() has been called. Cannot register more fields\n"); 117 ERROR(); 118 } 119 */ 120 121 /* check for repeated name */ 122 StringInList( field_name, db->nfields, (const DataField*)db->field, &val ); 123 if(val == PETSC_TRUE ) { 124 printf("ERROR: Cannot add same field twice\n"); 125 ERROR(); 126 } 127 128 /* create new space for data */ 129 field = realloc( db->field, sizeof(DataField)*(db->nfields+1)); 130 db->field = field; 131 132 /* add field */ 133 DataFieldCreate( registeration_function, field_name, atomic_size, db->allocated, &fp ); 134 db->field[ db->nfields ] = fp; 135 136 db->nfields++; 137 138 if(_gfield!=NULL){ 139 *_gfield = fp; 140 } 141 } 142 143 /* 144 #define DataBucketRegisterField(db,name,size,k) {\ 145 char *location;\ 146 asprintf(&location,"Registered by %s() at line %d within file %s", __FUNCTION__, __LINE__, __FILE__);\ 147 _DataBucketRegisterField( (db), location, (name), (size), (k) );\ 148 free(location);\ 149 } 150 */ 151 152 void DataBucketGetDataFieldByName(DataBucket db,const char name[],DataField *gfield) 153 { 154 PetscInt idx; 155 PetscBool found; 156 157 StringInList(name,db->nfields,(const DataField*)db->field,&found); 158 if(found==PETSC_FALSE) { 159 printf("ERROR: Cannot find DataField with name %s \n", name ); 160 ERROR(); 161 } 162 StringFindInList(name,db->nfields,(const DataField*)db->field,&idx); 163 164 *gfield = db->field[idx]; 165 } 166 167 void DataBucketQueryDataFieldByName(DataBucket db,const char name[],PetscBool *found) 168 { 169 *found = PETSC_FALSE; 170 StringInList(name,db->nfields,(const DataField*)db->field,found); 171 } 172 173 void DataBucketFinalize(DataBucket db) 174 { 175 db->finalised = PETSC_TRUE; 176 } 177 178 void DataFieldGetNumEntries(DataField df, PetscInt *sum) 179 { 180 *sum = df->L; 181 } 182 183 void DataFieldSetSize( DataField df, const PetscInt new_L ) 184 { 185 void *tmp_data; 186 187 if( new_L <= 0 ) { 188 printf("ERROR: Cannot set size of DataField to be <= 0 \n"); 189 ERROR(); 190 } 191 if( new_L == df->L ) return; 192 193 if( new_L > df->L ) { 194 195 tmp_data = realloc( df->data, df->atomic_size * (new_L) ); 196 df->data = tmp_data; 197 198 /* init new contents */ 199 memset( ( ((char*)df->data)+df->L*df->atomic_size), 0, (new_L-df->L)*df->atomic_size ); 200 201 } 202 else { 203 /* reallocate pointer list, add +1 in case new_L = 0 */ 204 tmp_data = realloc( df->data, df->atomic_size * (new_L+1) ); 205 df->data = tmp_data; 206 } 207 208 df->L = new_L; 209 } 210 211 void DataFieldZeroBlock( DataField df, const PetscInt start, const PetscInt end ) 212 { 213 if( start > end ) { 214 printf("ERROR: Cannot zero a block of entries if start(%d) > end(%d) \n",start,end); 215 ERROR(); 216 } 217 if( start < 0 ) { 218 printf("ERROR: Cannot zero a block of entries if start(%d) < 0 \n",start); 219 ERROR(); 220 } 221 if( end > df->L ) { 222 printf("ERROR: Cannot zero a block of entries if end(%d) >= array size(%d) \n",end,df->L); 223 ERROR(); 224 } 225 226 memset( ( ((char*)df->data)+start*df->atomic_size), 0, (end-start)*df->atomic_size ); 227 } 228 229 /* 230 A negative buffer value will simply be ignored and the old buffer value will be used. 231 */ 232 void DataBucketSetSizes( DataBucket db, const PetscInt L, const PetscInt buffer ) 233 { 234 PetscInt current_allocated,new_used,new_unused,new_buffer,new_allocated,f; 235 236 237 if( db->finalised == PETSC_FALSE ) { 238 printf("ERROR: You must call DataBucketFinalize() before DataBucketSetSizes() \n"); 239 ERROR(); 240 } 241 242 current_allocated = db->allocated; 243 244 new_used = L; 245 new_unused = current_allocated - new_used; 246 new_buffer = db->buffer; 247 if( buffer >= 0 ) { /* update the buffer value */ 248 new_buffer = buffer; 249 } 250 new_allocated = new_used + new_buffer; 251 252 /* action */ 253 if ( new_allocated > current_allocated ) { 254 /* increase size to new_used + new_buffer */ 255 for( f=0; f<db->nfields; f++ ) { 256 DataFieldSetSize( db->field[f], new_allocated ); 257 } 258 259 db->L = new_used; 260 db->buffer = new_buffer; 261 db->allocated = new_used + new_buffer; 262 } 263 else { 264 if( new_unused > 2 * new_buffer ) { 265 266 /* shrink array to new_used + new_buffer */ 267 for( f=0; f<db->nfields; f++ ) { 268 DataFieldSetSize( db->field[f], new_allocated ); 269 } 270 271 db->L = new_used; 272 db->buffer = new_buffer; 273 db->allocated = new_used + new_buffer; 274 } 275 else { 276 db->L = new_used; 277 db->buffer = new_buffer; 278 } 279 } 280 281 /* zero all entries from db->L to db->allocated */ 282 for( f=0; f<db->nfields; f++ ) { 283 DataField field = db->field[f]; 284 DataFieldZeroBlock(field, db->L,db->allocated); 285 } 286 } 287 288 void DataBucketSetInitialSizes( DataBucket db, const PetscInt L, const PetscInt buffer ) 289 { 290 PetscInt f; 291 DataBucketSetSizes(db,L,buffer); 292 293 for( f=0; f<db->nfields; f++ ) { 294 DataField field = db->field[f]; 295 DataFieldZeroBlock(field,0,db->allocated); 296 } 297 } 298 299 void DataBucketGetSizes( DataBucket db, PetscInt *L, PetscInt *buffer, PetscInt *allocated ) 300 { 301 if (L) { *L = db->L; } 302 if (buffer) { *buffer = db->buffer; } 303 if (allocated) { *allocated = db->allocated; } 304 } 305 306 void DataBucketGetGlobalSizes(MPI_Comm comm, DataBucket db, PetscInt *L, PetscInt *buffer, PetscInt *allocated ) 307 { 308 PetscInt _L,_buffer,_allocated; 309 PetscInt ierr; 310 311 _L = db->L; 312 _buffer = db->buffer; 313 _allocated = db->allocated; 314 315 if (L) { ierr = MPI_Allreduce(&_L,L,1,MPIU_INT,MPI_SUM,comm); } 316 if (buffer) { ierr = MPI_Allreduce(&_buffer,buffer,1,MPIU_INT,MPI_SUM,comm); } 317 if (allocated) { ierr = MPI_Allreduce(&_allocated,allocated,1,MPIU_INT,MPI_SUM,comm); } 318 } 319 320 void DataBucketGetDataFields( DataBucket db, PetscInt *L, DataField *fields[] ) 321 { 322 if(L){ *L = db->nfields; } 323 if(fields){ *fields = db->field; } 324 } 325 326 void DataFieldGetAccess( const DataField gfield ) 327 { 328 if(gfield->active==PETSC_TRUE) { 329 printf("ERROR: Field \"%s\" is already active. You must call DataFieldRestoreAccess()\n", gfield->name ); 330 ERROR(); 331 } 332 gfield->active = PETSC_TRUE; 333 } 334 335 void DataFieldAccessPoint( const DataField gfield, const PetscInt pid, void **ctx_p ) 336 { 337 #ifdef DATAFIELD_POINT_ACCESS_GUARD 338 /* debug mode */ 339 /* check poPetscInt is valid */ 340 if( pid < 0 ){ printf("ERROR: index must be >= 0\n"); ERROR(); } 341 if( pid >= gfield->L ){ printf("ERROR: index must be < %d\n",gfield->L); ERROR(); } 342 343 if(gfield->active==PETSC_FALSE) { 344 printf("ERROR: Field \"%s\" is not active. You must call DataFieldGetAccess() before poPetscInt data can be retrivied\n",gfield->name); 345 ERROR(); 346 } 347 #endif 348 349 //*ctx_p = (void*)( ((char*)gfield->data) + pid * gfield->atomic_size); 350 *ctx_p = __DATATFIELD_point_access(gfield->data,pid,gfield->atomic_size); 351 } 352 353 void DataFieldAccessPointOffset( const DataField gfield, const size_t offset, const PetscInt pid, void **ctx_p ) 354 { 355 #ifdef DATAFIELD_POINT_ACCESS_GUARD 356 /* debug mode */ 357 358 /* check poPetscInt is valid */ 359 /* if( offset < 0 ){ printf("ERROR: offset must be >= 0\n"); ERROR(); } *//* Note compiler realizes this can never happen with an unsigned PetscInt */ 360 if( offset >= gfield->atomic_size ){ printf("ERROR: offset must be < %zu\n",gfield->atomic_size); ERROR(); } 361 362 /* check poPetscInt is valid */ 363 if( pid < 0 ){ printf("ERROR: index must be >= 0\n"); ERROR(); } 364 if( pid >= gfield->L ){ printf("ERROR: index must be < %d\n",gfield->L); ERROR(); } 365 366 if(gfield->active==PETSC_FALSE) { 367 printf("ERROR: Field \"%s\" is not active. You must call DataFieldGetAccess() before poPetscInt data can be retrivied\n",gfield->name); 368 ERROR(); 369 } 370 #endif 371 372 *ctx_p = __DATATFIELD_point_access_offset(gfield->data,pid,gfield->atomic_size,offset); 373 } 374 375 void DataFieldRestoreAccess( DataField gfield ) 376 { 377 if(gfield->active==PETSC_FALSE) { 378 printf("ERROR: Field \"%s\" is not active. You must call DataFieldGetAccess()\n", gfield->name ); 379 ERROR(); 380 } 381 gfield->active = PETSC_FALSE; 382 } 383 384 void DataFieldVerifyAccess( const DataField gfield, const size_t size) 385 { 386 #ifdef DATAFIELD_POINT_ACCESS_GUARD 387 if(gfield->atomic_size != size ) { 388 printf("ERROR: Field \"%s\" must be mapped to %zu bytes, your intended structure is %zu bytes in length.\n", 389 gfield->name, gfield->atomic_size, size ); 390 ERROR(); 391 } 392 #endif 393 } 394 395 void DataFieldGetAtomicSize(const DataField gfield,size_t *size) 396 { 397 if (size) { *size = gfield->atomic_size; } 398 } 399 400 void DataFieldGetEntries(const DataField gfield,void **data) 401 { 402 if (data) { 403 *data = gfield->data; 404 } 405 } 406 407 void DataFieldRestoreEntries(const DataField gfield,void **data) 408 { 409 if (data) { 410 *data = NULL; 411 } 412 } 413 414 /* y = x */ 415 void DataBucketCopyPoint( const DataBucket xb, const PetscInt pid_x, 416 const DataBucket yb, const PetscInt pid_y ) 417 { 418 PetscInt f; 419 for( f=0; f<xb->nfields; f++ ) { 420 void *dest; 421 void *src; 422 423 DataFieldGetAccess( xb->field[f] ); 424 if (xb!=yb) { DataFieldGetAccess( yb->field[f] ); } 425 426 DataFieldAccessPoint( xb->field[f],pid_x, &src ); 427 DataFieldAccessPoint( yb->field[f],pid_y, &dest ); 428 429 memcpy( dest, src, xb->field[f]->atomic_size ); 430 431 DataFieldRestoreAccess( xb->field[f] ); 432 if (xb!=yb) { DataFieldRestoreAccess( yb->field[f] ); } 433 } 434 435 } 436 437 void DataBucketCreateFromSubset( DataBucket DBIn, const PetscInt N, const PetscInt list[], DataBucket *DB ) 438 { 439 PetscInt nfields; 440 DataField *fields; 441 DataBucketCreate(DB); 442 PetscInt f,L,buffer,allocated,p; 443 444 /* copy contents of DBIn */ 445 DataBucketGetDataFields(DBIn,&nfields,&fields); 446 DataBucketGetSizes(DBIn,&L,&buffer,&allocated); 447 448 for(f=0;f<nfields;f++) { 449 DataBucketRegisterField(*DB,fields[f]->name,fields[f]->atomic_size,NULL); 450 } 451 DataBucketFinalize(*DB); 452 453 DataBucketSetSizes(*DB,L,buffer); 454 455 /* now copy the desired guys from DBIn => DB */ 456 for( p=0; p<N; p++ ) { 457 DataBucketCopyPoint(DBIn,list[p], *DB,p); 458 } 459 460 } 461 462 // insert into an exisitng location 463 void DataFieldInsertPoint( const DataField field, const PetscInt index, const void *ctx ) 464 { 465 466 #ifdef DATAFIELD_POINT_ACCESS_GUARD 467 /* check poPetscInt is valid */ 468 if( index < 0 ){ printf("ERROR: index must be >= 0\n"); ERROR(); } 469 if( index >= field->L ){ printf("ERROR: index must be < %d\n",field->L); ERROR(); } 470 #endif 471 472 // memcpy( (void*)((char*)field->data + index*field->atomic_size), ctx, field->atomic_size ); 473 memcpy( __DATATFIELD_point_access(field->data,index,field->atomic_size), ctx, field->atomic_size ); 474 } 475 476 // remove data at index - replace with last point 477 void DataBucketRemovePointAtIndex( const DataBucket db, const PetscInt index ) 478 { 479 PetscInt f; 480 481 #ifdef DATAFIELD_POINT_ACCESS_GUARD 482 /* check poPetscInt is valid */ 483 if( index < 0 ){ printf("ERROR: index must be >= 0\n"); ERROR(); } 484 if( index >= db->allocated ){ printf("ERROR: index must be < %d\n",db->L+db->buffer); ERROR(); } 485 #endif 486 487 if (index >= db->L) { /* this poPetscInt is not in the list - no need to error, but I will anyway */ 488 printf("ERROR: You should not be trying to remove poPetscInt at index=%d since it's < db->L = %d \n", index, db->L ); 489 ERROR(); 490 } 491 492 #if 0 493 if (index == db->L-1) { /* last poPetscInt in list */ 494 for( f=0; f<db->nfields; f++ ) { 495 DataField field = db->field[f]; 496 497 DataFieldZeroPoint(field,index); 498 } 499 } 500 else { 501 for( f=0; f<db->nfields; f++ ) { 502 DataField field = db->field[f]; 503 504 /* copy then remove */ 505 DataFieldCopyPoint( db->L-1,field, index,field ); 506 507 DataFieldZeroPoint(field,index); 508 } 509 } 510 #endif 511 512 if (index != db->L-1) { /* not last poPetscInt in list */ 513 for( f=0; f<db->nfields; f++ ) { 514 DataField field = db->field[f]; 515 516 /* copy then remove */ 517 DataFieldCopyPoint( db->L-1,field, index,field ); 518 519 //DataFieldZeroPoint(field,index); 520 } 521 } 522 523 /* decrement size */ 524 /* this will zero out an crap at the end of the list */ 525 DataBucketRemovePoint(db); 526 527 } 528 529 /* copy x into y */ 530 void DataFieldCopyPoint( const PetscInt pid_x, const DataField field_x, 531 const PetscInt pid_y, const DataField field_y ) 532 { 533 534 #ifdef DATAFIELD_POINT_ACCESS_GUARD 535 /* check poPetscInt is valid */ 536 if( pid_x < 0 ){ printf("ERROR: (IN) index must be >= 0\n"); ERROR(); } 537 if( pid_x >= field_x->L ){ printf("ERROR: (IN) index must be < %d\n",field_x->L); ERROR(); } 538 539 if( pid_y < 0 ){ printf("ERROR: (OUT) index must be >= 0\n"); ERROR(); } 540 if( pid_y >= field_y->L ){ printf("ERROR: (OUT) index must be < %d\n",field_y->L); ERROR(); } 541 542 if( field_y->atomic_size != field_x->atomic_size ) { 543 printf("ERROR: atomic size must match \n"); ERROR(); 544 } 545 #endif 546 /* 547 memcpy( (void*)((char*)field_y->data + pid_y*field_y->atomic_size), 548 (void*)((char*)field_x->data + pid_x*field_x->atomic_size), 549 field_x->atomic_size ); 550 */ 551 memcpy( __DATATFIELD_point_access(field_y->data,pid_y,field_y->atomic_size), 552 __DATATFIELD_point_access(field_x->data,pid_x,field_x->atomic_size), 553 field_y->atomic_size ); 554 555 } 556 557 558 // zero only the datafield at this point 559 void DataFieldZeroPoint( const DataField field, const PetscInt index ) 560 { 561 #ifdef DATAFIELD_POINT_ACCESS_GUARD 562 /* check poPetscInt is valid */ 563 if( index < 0 ){ printf("ERROR: index must be >= 0\n"); ERROR(); } 564 if( index >= field->L ){ printf("ERROR: index must be < %d\n",field->L); ERROR(); } 565 #endif 566 567 // memset( (void*)((char*)field->data + index*field->atomic_size), 0, field->atomic_size ); 568 memset( __DATATFIELD_point_access(field->data,index,field->atomic_size), 0, field->atomic_size ); 569 } 570 571 // zero ALL data for this point 572 void DataBucketZeroPoint( const DataBucket db, const PetscInt index ) 573 { 574 PetscInt f; 575 576 /* check poPetscInt is valid */ 577 if( index < 0 ){ printf("ERROR: index must be >= 0\n"); ERROR(); } 578 if( index >= db->allocated ){ printf("ERROR: index must be < %d\n",db->allocated); ERROR(); } 579 580 for(f=0;f<db->nfields;f++){ 581 DataField field = db->field[f]; 582 583 DataFieldZeroPoint(field,index); 584 } 585 } 586 587 /* increment */ 588 void DataBucketAddPoint( DataBucket db ) 589 { 590 DataBucketSetSizes( db, db->L+1, -1 ); 591 } 592 /* decrement */ 593 void DataBucketRemovePoint( DataBucket db ) 594 { 595 DataBucketSetSizes( db, db->L-1, -1 ); 596 } 597 598 void _DataFieldViewBinary(DataField field, FILE *fp ) 599 { 600 fprintf(fp,"<DataField>\n"); 601 fprintf(fp,"%d\n", field->L); 602 fprintf(fp,"%zu\n",field->atomic_size); 603 fprintf(fp,"%s\n", field->registeration_function); 604 fprintf(fp,"%s\n", field->name); 605 606 fwrite(field->data, field->atomic_size, field->L, fp); 607 /* 608 printf(" ** wrote %zu bytes for DataField \"%s\" \n", field->atomic_size * field->L, field->name ); 609 */ 610 fprintf(fp,"\n</DataField>\n"); 611 } 612 613 void _DataBucketRegisterFieldFromFile( FILE *fp, DataBucket db ) 614 { 615 PetscBool val; 616 DataField *field; 617 618 DataField gfield; 619 char dummy[100]; 620 char registeration_function[5000]; 621 char field_name[5000]; 622 PetscInt L; 623 size_t atomic_size,strL; 624 625 626 /* check we haven't finalised the registration of fields */ 627 /* 628 if(db->finalised==PETSC_TRUE) { 629 printf("ERROR: DataBucketFinalize() has been called. Cannot register more fields\n"); 630 ERROR(); 631 } 632 */ 633 634 635 /* read file contents */ 636 fgets(dummy,99,fp); //printf("read(header): %s", dummy ); 637 638 fscanf( fp, "%d\n",&L); //printf("read(L): %d\n", L); 639 640 fscanf( fp, "%zu\n",&atomic_size); //printf("read(size): %zu\n",atomic_size); 641 642 fgets(registeration_function,4999,fp); //printf("read(reg func): %s", registeration_function ); 643 strL = strlen(registeration_function); 644 if(strL>1){ 645 registeration_function[strL-1] = 0; 646 } 647 648 fgets(field_name,4999,fp); //printf("read(name): %s", field_name ); 649 strL = strlen(field_name); 650 if(strL>1){ 651 field_name[strL-1] = 0; 652 } 653 654 #ifdef DATA_BUCKET_LOG 655 printf(" ** read L=%d; atomic_size=%zu; reg_func=\"%s\"; name=\"%s\" \n", L,atomic_size,registeration_function,field_name); 656 #endif 657 658 659 /* check for repeated name */ 660 StringInList( field_name, db->nfields, (const DataField*)db->field, &val ); 661 if(val == PETSC_TRUE ) { 662 printf("ERROR: Cannot add same field twice\n"); 663 ERROR(); 664 } 665 666 /* create new space for data */ 667 field = realloc( db->field, sizeof(DataField)*(db->nfields+1)); 668 db->field = field; 669 670 /* add field */ 671 DataFieldCreate( registeration_function, field_name, atomic_size, L, &gfield ); 672 673 /* copy contents of file */ 674 fread(gfield->data, gfield->atomic_size, gfield->L, fp); 675 #ifdef DATA_BUCKET_LOG 676 printf(" ** read %zu bytes for DataField \"%s\" \n", gfield->atomic_size * gfield->L, field_name ); 677 #endif 678 /* finish reading meta data */ 679 fgets(dummy,99,fp); //printf("read(header): %s", dummy ); 680 fgets(dummy,99,fp); //printf("read(header): %s", dummy ); 681 682 db->field[ db->nfields ] = gfield; 683 684 db->nfields++; 685 686 } 687 688 void _DataBucketViewAscii_HeaderWrite_v00(FILE *fp) 689 { 690 fprintf(fp,"<DataBucketHeader>\n"); 691 fprintf(fp,"type=DataBucket\n"); 692 fprintf(fp,"format=ascii\n"); 693 fprintf(fp,"version=0.0\n"); 694 fprintf(fp,"options=\n"); 695 fprintf(fp,"</DataBucketHeader>\n"); 696 } 697 void _DataBucketViewAscii_HeaderRead_v00(FILE *fp) 698 { 699 char dummy[100]; 700 size_t strL; 701 702 // header open 703 fgets(dummy,99,fp); //printf("read(header): %s", dummy ); 704 705 // type 706 fgets(dummy,99,fp); //printf("read(header): %s", dummy ); 707 strL = strlen(dummy); 708 if(strL>1) { dummy[strL-1] = 0; } 709 if(strcmp(dummy,"type=DataBucket")!=0) { 710 printf("ERROR: Data file doesn't contain a DataBucket type\n"); 711 ERROR(); 712 } 713 714 // format 715 fgets(dummy,99,fp); //printf("read(header): %s", dummy ); 716 717 // version 718 fgets(dummy,99,fp); //printf("read(header): %s", dummy ); 719 strL = strlen(dummy); 720 if(strL>1) { dummy[strL-1] = 0; } 721 if(strcmp(dummy,"version=0.0")!=0) { 722 printf("ERROR: DataBucket file must be parsed with version=0.0 : You tried %s \n", dummy); 723 ERROR(); 724 } 725 726 // options 727 fgets(dummy,99,fp); //printf("read(header): %s", dummy ); 728 // header close 729 fgets(dummy,99,fp); //printf("read(header): %s", dummy ); 730 } 731 732 733 void _DataBucketLoadFromFileBinary_SEQ(const char filename[], DataBucket *_db) 734 { 735 DataBucket db; 736 FILE *fp; 737 PetscInt L,buffer,f,nfields; 738 739 740 #ifdef DATA_BUCKET_LOG 741 printf("** DataBucketLoadFromFile **\n"); 742 #endif 743 744 /* open file */ 745 fp = fopen(filename,"rb"); 746 if(fp==NULL){ 747 printf("ERROR: Cannot open file with name %s \n", filename); 748 ERROR(); 749 } 750 751 /* read header */ 752 _DataBucketViewAscii_HeaderRead_v00(fp); 753 754 fscanf(fp,"%d\n%d\n%d\n",&L,&buffer,&nfields); 755 756 DataBucketCreate(&db); 757 758 for( f=0; f<nfields; f++ ) { 759 _DataBucketRegisterFieldFromFile(fp,db); 760 } 761 fclose(fp); 762 763 DataBucketFinalize(db); 764 765 766 /* 767 DataBucketSetSizes(db,L,buffer); 768 */ 769 db->L = L; 770 db->buffer = buffer; 771 db->allocated = L + buffer; 772 773 *_db = db; 774 } 775 776 void DataBucketLoadFromFile(MPI_Comm comm,const char filename[], DataBucketViewType type, DataBucket *db) 777 { 778 PetscMPIInt nproc,rank; 779 780 MPI_Comm_size(comm,&nproc); 781 MPI_Comm_rank(comm,&rank); 782 783 #ifdef DATA_BUCKET_LOG 784 printf("** DataBucketLoadFromFile **\n"); 785 #endif 786 if(type==DATABUCKET_VIEW_STDOUT) { 787 788 } else if(type==DATABUCKET_VIEW_ASCII) { 789 printf("ERROR: Cannot be implemented as we don't know the underlying particle data structure\n"); 790 ERROR(); 791 } else if(type==DATABUCKET_VIEW_BINARY) { 792 if (nproc==1) { 793 _DataBucketLoadFromFileBinary_SEQ(filename,db); 794 } else { 795 char *name; 796 797 asprintf(&name,"%s_p%1.5d",filename, rank ); 798 _DataBucketLoadFromFileBinary_SEQ(name,db); 799 free(name); 800 } 801 } else { 802 printf("ERROR: Not implemented\n"); 803 ERROR(); 804 } 805 } 806 807 808 void _DataBucketViewBinary(DataBucket db,const char filename[]) 809 { 810 FILE *fp = NULL; 811 PetscInt f; 812 813 fp = fopen(filename,"wb"); 814 if(fp==NULL){ 815 printf("ERROR: Cannot open file with name %s \n", filename); 816 ERROR(); 817 } 818 819 /* db header */ 820 _DataBucketViewAscii_HeaderWrite_v00(fp); 821 822 /* meta-data */ 823 fprintf(fp,"%d\n%d\n%d\n", db->L,db->buffer,db->nfields); 824 825 for( f=0; f<db->nfields; f++ ) { 826 /* load datafields */ 827 _DataFieldViewBinary(db->field[f],fp); 828 } 829 830 fclose(fp); 831 } 832 833 void DataBucketView_SEQ(DataBucket db,const char filename[],DataBucketViewType type) 834 { 835 switch (type) { 836 case DATABUCKET_VIEW_STDOUT: 837 { 838 PetscInt f; 839 double memory_usage_total = 0.0; 840 841 printf("DataBucketView(SEQ): (\"%s\")\n",filename); 842 printf(" L = %d \n", db->L ); 843 printf(" buffer = %d \n", db->buffer ); 844 printf(" allocated = %d \n", db->allocated ); 845 846 printf(" nfields registered = %d \n", db->nfields ); 847 for( f=0; f<db->nfields; f++ ) { 848 double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6; 849 850 printf(" [%3d]: field name ==>> %30s : Mem. usage = %1.2e (MB) \n", f, db->field[f]->name, memory_usage_f ); 851 memory_usage_total += memory_usage_f; 852 } 853 printf(" Total mem. usage = %1.2e (MB) \n", memory_usage_total ); 854 } 855 break; 856 857 case DATABUCKET_VIEW_ASCII: 858 { 859 printf("ERROR: Cannot be implemented as we don't know the underlying particle data structure\n"); 860 ERROR(); 861 } 862 break; 863 864 case DATABUCKET_VIEW_BINARY: 865 { 866 _DataBucketViewBinary(db,filename); 867 } 868 break; 869 870 case DATABUCKET_VIEW_HDF5: 871 { 872 printf("ERROR: Has not been implemented \n"); 873 ERROR(); 874 } 875 break; 876 877 default: 878 printf("ERROR: Unknown method requested \n"); 879 ERROR(); 880 break; 881 } 882 } 883 884 void DataBucketView_MPI(MPI_Comm comm,DataBucket db,const char filename[],DataBucketViewType type) 885 { 886 switch (type) { 887 case DATABUCKET_VIEW_STDOUT: 888 { 889 PetscInt f; 890 PetscInt L,buffer,allocated; 891 double memory_usage_total,memory_usage_total_local = 0.0; 892 PetscMPIInt rank; 893 PetscInt ierr; 894 895 ierr = MPI_Comm_rank(comm,&rank); 896 897 DataBucketGetGlobalSizes(comm,db,&L,&buffer,&allocated); 898 899 for( f=0; f<db->nfields; f++ ) { 900 double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6; 901 902 memory_usage_total_local += memory_usage_f; 903 } 904 MPI_Allreduce(&memory_usage_total_local,&memory_usage_total,1,MPI_DOUBLE,MPI_SUM,comm); 905 906 if (rank==0) { 907 PetscPrintf(comm,"DataBucketView(MPI): (\"%s\")\n",filename); 908 PetscPrintf(comm," L = %D \n", L ); 909 PetscPrintf(comm," buffer (max) = %D \n", buffer ); 910 PetscPrintf(comm," allocated = %D \n", allocated ); 911 912 PetscPrintf(comm," nfields registered = %D \n", db->nfields ); 913 for( f=0; f<db->nfields; f++ ) { 914 double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6; 915 916 printf(" [%3d]: field name ==>> %30s : Mem. usage = %1.2e (MB) : rank0\n", f, db->field[f]->name, memory_usage_f ); 917 } 918 919 printf(" Total mem. usage = %1.2e (MB) : collective\n", memory_usage_total ); 920 } 921 922 } 923 break; 924 925 case DATABUCKET_VIEW_ASCII: 926 { 927 printf("ERROR: Cannot be implemented as we don't know the underlying particle data structure\n"); 928 ERROR(); 929 } 930 break; 931 932 case DATABUCKET_VIEW_BINARY: 933 { 934 char *name; 935 PetscMPIInt rank; 936 937 /* create correct extension */ 938 MPI_Comm_rank(comm,&rank); 939 asprintf(&name,"%s_p%1.5d",filename, rank ); 940 941 _DataBucketViewBinary(db,name); 942 943 free(name); 944 } 945 break; 946 947 case DATABUCKET_VIEW_HDF5: 948 { 949 printf("ERROR: Has not been implemented \n"); 950 ERROR(); 951 } 952 break; 953 954 default: 955 printf("ERROR: Unknown method requested \n"); 956 ERROR(); 957 break; 958 } 959 } 960 961 962 void DataBucketView(MPI_Comm comm,DataBucket db,const char filename[],DataBucketViewType type) 963 { 964 PetscMPIInt nproc; 965 966 MPI_Comm_size(comm,&nproc); 967 if (nproc==1) { 968 DataBucketView_SEQ(db,filename,type); 969 } else { 970 DataBucketView_MPI(comm,db,filename,type); 971 } 972 } 973 974 void DataBucketDuplicateFields(DataBucket dbA,DataBucket *dbB) 975 { 976 DataBucket db2; 977 PetscInt f; 978 979 DataBucketCreate(&db2); 980 981 /* copy contents from dbA into db2 */ 982 for (f=0; f<dbA->nfields; f++) { 983 DataField field; 984 size_t atomic_size; 985 char *name; 986 987 field = dbA->field[f]; 988 989 atomic_size = field->atomic_size; 990 name = field->name; 991 992 DataBucketRegisterField(db2,name,atomic_size,NULL); 993 } 994 DataBucketFinalize(db2); 995 DataBucketSetInitialSizes(db2,0,1000); 996 997 /* set pointer */ 998 *dbB = db2; 999 } 1000 1001 /* 1002 Insert points from db2 into db1 1003 db1 <<== db2 1004 */ 1005 void DataBucketInsertValues(DataBucket db1,DataBucket db2) 1006 { 1007 PetscInt n_mp_points1,n_mp_points2; 1008 PetscInt n_mp_points1_new,p; 1009 1010 DataBucketGetSizes(db1,&n_mp_points1,0,0); 1011 DataBucketGetSizes(db2,&n_mp_points2,0,0); 1012 1013 n_mp_points1_new = n_mp_points1 + n_mp_points2; 1014 DataBucketSetSizes(db1,n_mp_points1_new,-1); 1015 1016 for (p=0; p<n_mp_points2; p++) { 1017 // db1 <<== db2 // 1018 DataBucketCopyPoint( db2,p, db1,(n_mp_points1 + p) ); 1019 } 1020 } 1021 1022 /* helpers for parallel send/recv */ 1023 void DataBucketCreatePackedArray(DataBucket db,size_t *bytes,void **buf) 1024 { 1025 PetscInt f; 1026 size_t sizeof_marker_contents; 1027 void *buffer; 1028 1029 sizeof_marker_contents = 0; 1030 for (f=0; f<db->nfields; f++) { 1031 DataField df = db->field[f]; 1032 1033 sizeof_marker_contents += df->atomic_size; 1034 } 1035 1036 buffer = malloc(sizeof_marker_contents); 1037 memset(buffer,0,sizeof_marker_contents); 1038 1039 if (bytes) { *bytes = sizeof_marker_contents; } 1040 if (buf) { *buf = buffer; } 1041 } 1042 1043 void DataBucketDestroyPackedArray(DataBucket db,void **buf) 1044 { 1045 if (buf) { 1046 free(*buf); 1047 *buf = NULL; 1048 } 1049 } 1050 1051 void DataBucketFillPackedArray(DataBucket db,const PetscInt index,void *buf) 1052 { 1053 PetscInt f; 1054 void *data,*data_p; 1055 size_t asize,offset; 1056 1057 offset = 0; 1058 for (f=0; f<db->nfields; f++) { 1059 DataField df = db->field[f]; 1060 1061 asize = df->atomic_size; 1062 1063 data = (void*)( df->data ); 1064 data_p = (void*)( (char*)data + index*asize ); 1065 1066 memcpy( (void*)((char*)buf + offset), data_p, asize); 1067 offset = offset + asize; 1068 } 1069 } 1070 1071 void DataBucketInsertPackedArray(DataBucket db,const PetscInt idx,void *data) 1072 { 1073 PetscInt f; 1074 void *data_p; 1075 size_t offset; 1076 1077 offset = 0; 1078 for (f=0; f<db->nfields; f++) { 1079 DataField df = db->field[f]; 1080 1081 data_p = (void*)( (char*)data + offset ); 1082 1083 DataFieldInsertPoint(df, idx, (void*)data_p ); 1084 offset = offset + df->atomic_size; 1085 } 1086 } 1087