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