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