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