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