Actual source code: data_bucket.c

petsc-3.11.4 2019-09-28
Report Typos and Errors
  1: #include "../src/dm/impls/swarm/data_bucket.h"

  3: /* string helpers */
  4: PetscErrorCode DMSwarmDataFieldStringInList(const char name[],const PetscInt N,const DMSwarmDataField gfield[],PetscBool *val)
  5: {
  6:   PetscInt       i;

 10:   *val = PETSC_FALSE;
 11:   for (i = 0; i < N; ++i) {
 12:     PetscBool flg;
 13:     PetscStrcmp(name, gfield[i]->name, &flg);
 14:     if (flg) {
 15:       *val = PETSC_TRUE;
 16:       return(0);
 17:     }
 18:   }
 19:   return(0);
 20: }

 22: PetscErrorCode DMSwarmDataFieldStringFindInList(const char name[],const PetscInt N,const DMSwarmDataField gfield[],PetscInt *index)
 23: {
 24:   PetscInt       i;

 28:   *index = -1;
 29:   for (i = 0; i < N; ++i) {
 30:     PetscBool flg;
 31:     PetscStrcmp(name, gfield[i]->name, &flg);
 32:     if (flg) {
 33:       *index = i;
 34:       return(0);
 35:     }
 36:   }
 37:   return(0);
 38: }

 40: PetscErrorCode DMSwarmDataFieldCreate(const char registration_function[],const char name[],const size_t size,const PetscInt L,DMSwarmDataField *DF)
 41: {
 42:   DMSwarmDataField df;
 43:   PetscErrorCode   ierr;

 46:   PetscMalloc(sizeof(struct _p_DMSwarmDataField), &df);
 47:   PetscMemzero(df, sizeof(struct _p_DMSwarmDataField));
 48:   PetscStrallocpy(registration_function, &df->registration_function);
 49:   PetscStrallocpy(name, &df->name);
 50:   df->atomic_size = size;
 51:   df->L  = L;
 52:   df->bs = 1;
 53:   /* allocate something so we don't have to reallocate */
 54:   PetscMalloc(size * L, &df->data);
 55:   PetscMemzero(df->data, size * L);
 56:   *DF = df;
 57:   return(0);
 58: }

 60: PetscErrorCode DMSwarmDataFieldDestroy(DMSwarmDataField *DF)
 61: {
 62:   DMSwarmDataField df = *DF;
 63:   PetscErrorCode   ierr;

 66:   PetscFree(df->registration_function);
 67:   PetscFree(df->name);
 68:   PetscFree(df->data);
 69:   PetscFree(df);
 70:   *DF  = NULL;
 71:   return(0);
 72: }

 74: /* data bucket */
 75: PetscErrorCode DMSwarmDataBucketCreate(DMSwarmDataBucket *DB)
 76: {
 77:   DMSwarmDataBucket db;
 78:   PetscErrorCode    ierr;

 81:   PetscMalloc(sizeof(struct _p_DMSwarmDataBucket), &db);
 82:   PetscMemzero(db, sizeof(struct _p_DMSwarmDataBucket));

 84:   db->finalised = PETSC_FALSE;
 85:   /* create empty spaces for fields */
 86:   db->L         = -1;
 87:   db->buffer    = 1;
 88:   db->allocated = 1;
 89:   db->nfields   = 0;
 90:   PetscMalloc1(1, &db->field);
 91:   *DB  = db;
 92:   return(0);
 93: }

 95: PetscErrorCode DMSwarmDataBucketDestroy(DMSwarmDataBucket *DB)
 96: {
 97:   DMSwarmDataBucket db = *DB;
 98:   PetscInt          f;
 99:   PetscErrorCode    ierr;

102:   /* release fields */
103:   for (f = 0; f < db->nfields; ++f) {
104:     DMSwarmDataFieldDestroy(&db->field[f]);
105:   }
106:   /* this will catch the initially allocated objects in the event that no fields are registered */
107:   if (db->field != NULL) {
108:     PetscFree(db->field);
109:   }
110:   PetscFree(db);
111:   *DB = NULL;
112:   return(0);
113: }

115: PetscErrorCode DMSwarmDataBucketQueryForActiveFields(DMSwarmDataBucket db,PetscBool *any_active_fields)
116: {
117:   PetscInt f;

120:   *any_active_fields = PETSC_FALSE;
121:   for (f = 0; f < db->nfields; ++f) {
122:     if (db->field[f]->active) {
123:       *any_active_fields = PETSC_TRUE;
124:       return(0);
125:     }
126:   }
127:   return(0);
128: }

130: PetscErrorCode DMSwarmDataBucketRegisterField(
131:                               DMSwarmDataBucket db,
132:                               const char registration_function[],
133:                               const char field_name[],
134:                               size_t atomic_size, DMSwarmDataField *_gfield)
135: {
136:   PetscBool val;
137:   DMSwarmDataField fp;

141:         /* check we haven't finalised the registration of fields */
142:         /*
143:    if(db->finalised==PETSC_TRUE) {
144:    printf("ERROR: DMSwarmDataBucketFinalize() has been called. Cannot register more fields\n");
145:    ERROR();
146:    }
147:    */
148:   /* check for repeated name */
149:   DMSwarmDataFieldStringInList(field_name, db->nfields, (const DMSwarmDataField*) db->field, &val);
150:   if (val == PETSC_TRUE) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field %s already exists. Cannot add same field twice",field_name);
151:   /* create new space for data */
152:   PetscRealloc(sizeof(DMSwarmDataField)*(db->nfields+1), &db->field);
153:   /* add field */
154:   DMSwarmDataFieldCreate(registration_function, field_name, atomic_size, db->allocated, &fp);
155:   db->field[db->nfields] = fp;
156:   db->nfields++;
157:   if (_gfield != NULL) {
158:     *_gfield = fp;
159:   }
160:   return(0);
161: }

163: /*
164:  #define DMSwarmDataBucketRegisterField(db,name,size,k) {\
165:  char *location;\
166:  asprintf(&location,"Registered by %s() at line %d within file %s", __FUNCTION__, __LINE__, __FILE__);\
167:  _DMSwarmDataBucketRegisterField( (db), location, (name), (size), (k) );\
168:  PetscFree(location);\
169:  }
170:  */

172: PetscErrorCode DMSwarmDataBucketGetDMSwarmDataFieldByName(DMSwarmDataBucket db,const char name[],DMSwarmDataField *gfield)
173: {
174:   PetscInt       idx;
175:   PetscBool      found;

179:   DMSwarmDataFieldStringInList(name,db->nfields,(const DMSwarmDataField*)db->field,&found);
180:   if (!found) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot find DMSwarmDataField with name %s",name);
181:   DMSwarmDataFieldStringFindInList(name,db->nfields,(const DMSwarmDataField*)db->field,&idx);
182:   *gfield = db->field[idx];
183:   return(0);
184: }

186: PetscErrorCode DMSwarmDataBucketQueryDMSwarmDataFieldByName(DMSwarmDataBucket db,const char name[],PetscBool *found)
187: {

191:   *found = PETSC_FALSE;
192:   DMSwarmDataFieldStringInList(name,db->nfields,(const DMSwarmDataField*)db->field,found);
193:   return(0);
194: }

196: PetscErrorCode DMSwarmDataBucketFinalize(DMSwarmDataBucket db)
197: {
199:   db->finalised = PETSC_TRUE;
200:   return(0);
201: }

203: PetscErrorCode DMSwarmDataFieldGetNumEntries(DMSwarmDataField df,PetscInt *sum)
204: {
206:   *sum = df->L;
207:   return(0);
208: }

210: PetscErrorCode DMSwarmDataFieldSetBlockSize(DMSwarmDataField df,PetscInt blocksize)
211: {
213:   df->bs = blocksize;
214:   return(0);
215: }

217: PetscErrorCode DMSwarmDataFieldSetSize(DMSwarmDataField df,const PetscInt new_L)
218: {

222:   if (new_L < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot set size of DMSwarmDataField to be < 0");
223:   if (new_L == df->L) return(0);
224:   if (new_L > df->L) {
225:     PetscRealloc(df->atomic_size * (new_L), &df->data);
226:     /* init new contents */
227:     PetscMemzero(( ((char*)df->data)+df->L*df->atomic_size), (new_L-df->L)*df->atomic_size);
228:   } else {
229:     /* reallocate pointer list, add +1 in case new_L = 0 */
230:     PetscRealloc(df->atomic_size * (new_L+1), &df->data);
231:   }
232:   df->L = new_L;
233:   return(0);
234: }

236: PetscErrorCode DMSwarmDataFieldZeroBlock(DMSwarmDataField df,const PetscInt start,const PetscInt end)
237: {

241:   if (start > end) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if start(%D) > end(%D)",start,end);
242:   if (start < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if start(%D) < 0",start);
243:   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);
244:   PetscMemzero((((char*)df->data)+start*df->atomic_size), (end-start)*df->atomic_size);
245:   return(0);
246: }

248: /*
249:  A negative buffer value will simply be ignored and the old buffer value will be used.
250:  */
251: PetscErrorCode DMSwarmDataBucketSetSizes(DMSwarmDataBucket db,const PetscInt L,const PetscInt buffer)
252: {
253:   PetscInt       current_allocated,new_used,new_unused,new_buffer,new_allocated,f;
254:   PetscBool      any_active_fields;

258:   if (db->finalised == PETSC_FALSE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"You must call DMSwarmDataBucketFinalize() before DMSwarmDataBucketSetSizes()");
259:   DMSwarmDataBucketQueryForActiveFields(db,&any_active_fields);
260:   if (any_active_fields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot safely re-size as at least one DMSwarmDataField is currently being accessed");

262:   current_allocated = db->allocated;
263:   new_used   = L;
264:   new_unused = current_allocated - new_used;
265:   new_buffer = db->buffer;
266:   if (buffer >= 0) { /* update the buffer value */
267:     new_buffer = buffer;
268:   }
269:   new_allocated = new_used + new_buffer;
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:       DMSwarmDataFieldSetSize(db->field[f], new_allocated);
275:     }
276:     db->L         = new_used;
277:     db->buffer    = new_buffer;
278:     db->allocated = new_used + new_buffer;
279:   } else {
280:     if (new_unused > 2 * new_buffer) {
281:       /* shrink array to new_used + new_buffer */
282:       for (f = 0; f < db->nfields; ++f) {
283:         DMSwarmDataFieldSetSize(db->field[f], new_allocated);
284:       }
285:       db->L         = new_used;
286:       db->buffer    = new_buffer;
287:       db->allocated = new_used + new_buffer;
288:     } else {
289:       db->L      = new_used;
290:       db->buffer = new_buffer;
291:     }
292:   }
293:   /* zero all entries from db->L to db->allocated */
294:   for (f = 0; f < db->nfields; ++f) {
295:     DMSwarmDataField field = db->field[f];
296:     DMSwarmDataFieldZeroBlock(field, db->L,db->allocated);
297:   }
298:   return(0);
299: }

301: PetscErrorCode DMSwarmDataBucketSetInitialSizes(DMSwarmDataBucket db,const PetscInt L,const PetscInt buffer)
302: {
303:   PetscInt       f;

307:   DMSwarmDataBucketSetSizes(db,L,buffer);
308:   for (f = 0; f < db->nfields; ++f) {
309:     DMSwarmDataField field = db->field[f];
310:     DMSwarmDataFieldZeroBlock(field,0,db->allocated);
311:   }
312:   return(0);
313: }

315: PetscErrorCode DMSwarmDataBucketGetSizes(DMSwarmDataBucket db,PetscInt *L,PetscInt *buffer,PetscInt *allocated)
316: {
318:   if (L) {*L = db->L;}
319:   if (buffer) {*buffer = db->buffer;}
320:   if (allocated) {*allocated = db->allocated;}
321:   return(0);
322: }

324: PetscErrorCode DMSwarmDataBucketGetGlobalSizes(MPI_Comm comm,DMSwarmDataBucket db,PetscInt *L,PetscInt *buffer,PetscInt *allocated)
325: {
326:   PetscInt ierr;

329:   if (L) {         MPI_Allreduce(&db->L,L,1,MPIU_INT,MPI_SUM,comm); }
330:   if (buffer) {    MPI_Allreduce(&db->buffer,buffer,1,MPIU_INT,MPI_SUM,comm); }
331:   if (allocated) { MPI_Allreduce(&db->allocated,allocated,1,MPIU_INT,MPI_SUM,comm); }
332:   return(0);
333: }

335: PetscErrorCode DMSwarmDataBucketGetDMSwarmDataFields(DMSwarmDataBucket db,PetscInt *L,DMSwarmDataField *fields[])
336: {
338:   if (L)      {*L      = db->nfields;}
339:   if (fields) {*fields = db->field;}
340:   return(0);
341: }

343: PetscErrorCode DMSwarmDataFieldGetAccess(const DMSwarmDataField gfield)
344: {
346:   if (gfield->active) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is already active. You must call DMSwarmDataFieldRestoreAccess()",gfield->name);
347:   gfield->active = PETSC_TRUE;
348:   return(0);
349: }

351: PetscErrorCode DMSwarmDataFieldAccessPoint(const DMSwarmDataField gfield,const PetscInt pid,void **ctx_p)
352: {
354:   *ctx_p = NULL;
355: #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
356:   /* debug mode */
357:   /* check point is valid */
358:   if (pid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
359:   if (pid >= gfield->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",gfield->L);
360:   if (gfield->active == PETSC_FALSE) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is not active. You must call DMSwarmDataFieldGetAccess() before point data can be retrivied",gfield->name);
361: #endif
362:   *ctx_p = DMSWARM_DATAFIELD_point_access(gfield->data,pid,gfield->atomic_size);
363:   return(0);
364: }

366: PetscErrorCode DMSwarmDataFieldAccessPointOffset(const DMSwarmDataField gfield,const size_t offset,const PetscInt pid,void **ctx_p)
367: {
369: #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
370:   /* debug mode */
371:   /* check point is valid */
372:   /* if( offset < 0 ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be >= 0");*/
373:   /* Note compiler realizes this can never happen with an unsigned PetscInt */
374:   if (offset >= gfield->atomic_size) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be < %zu",gfield->atomic_size);
375:   /* check point is valid */
376:   if (pid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
377:   if (pid >= gfield->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",gfield->L);
378:   if (gfield->active == PETSC_FALSE) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is not active. You must call DMSwarmDataFieldGetAccess() before point data can be retrivied",gfield->name);
379: #endif
380:   *ctx_p = DMSWARM_DATAFIELD_point_access_offset(gfield->data,pid,gfield->atomic_size,offset);
381:   return(0);
382: }

384: PetscErrorCode DMSwarmDataFieldRestoreAccess(DMSwarmDataField gfield)
385: {
387:   if (gfield->active == PETSC_FALSE) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is not active. You must call DMSwarmDataFieldGetAccess()", gfield->name);
388:   gfield->active = PETSC_FALSE;
389:   return(0);
390: }

392: PetscErrorCode DMSwarmDataFieldVerifyAccess(const DMSwarmDataField gfield,const size_t size)
393: {
395: #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
396:   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 );
397: #endif
398:   return(0);
399: }

401: PetscErrorCode DMSwarmDataFieldGetAtomicSize(const DMSwarmDataField gfield,size_t *size)
402: {
404:   if (size) {*size = gfield->atomic_size;}
405:   return(0);
406: }

408: PetscErrorCode DMSwarmDataFieldGetEntries(const DMSwarmDataField gfield,void **data)
409: {
411:   if (data) {*data = gfield->data;}
412:   return(0);
413: }

415: PetscErrorCode DMSwarmDataFieldRestoreEntries(const DMSwarmDataField gfield,void **data)
416: {
418:   if (data) {*data = NULL;}
419:   return(0);
420: }

422: /* y = x */
423: PetscErrorCode DMSwarmDataBucketCopyPoint(const DMSwarmDataBucket xb,const PetscInt pid_x,
424:                          const DMSwarmDataBucket yb,const PetscInt pid_y)
425: {
426:   PetscInt f;

430:   for (f = 0; f < xb->nfields; ++f) {
431:     void *dest;
432:     void *src;

434:     DMSwarmDataFieldGetAccess(xb->field[f]);
435:     if (xb != yb) { DMSwarmDataFieldGetAccess( yb->field[f]); }
436:     DMSwarmDataFieldAccessPoint(xb->field[f],pid_x, &src);
437:     DMSwarmDataFieldAccessPoint(yb->field[f],pid_y, &dest);
438:     PetscMemcpy(dest, src, xb->field[f]->atomic_size);
439:     DMSwarmDataFieldRestoreAccess(xb->field[f]);
440:     if (xb != yb) {DMSwarmDataFieldRestoreAccess(yb->field[f]);}
441:   }
442:   return(0);
443: }

445: PetscErrorCode DMSwarmDataBucketCreateFromSubset(DMSwarmDataBucket DBIn,const PetscInt N,const PetscInt list[],DMSwarmDataBucket *DB)
446: {
447:   PetscInt nfields;
448:   DMSwarmDataField *fields;
449:   PetscInt f,L,buffer,allocated,p;

453:   DMSwarmDataBucketCreate(DB);
454:   /* copy contents of DBIn */
455:   DMSwarmDataBucketGetDMSwarmDataFields(DBIn,&nfields,&fields);
456:   DMSwarmDataBucketGetSizes(DBIn,&L,&buffer,&allocated);
457:   for (f = 0; f < nfields; ++f) {
458:     DMSwarmDataBucketRegisterField(*DB,"DMSwarmDataBucketCreateFromSubset",fields[f]->name,fields[f]->atomic_size,NULL);
459:   }
460:   DMSwarmDataBucketFinalize(*DB);
461:   DMSwarmDataBucketSetSizes(*DB,L,buffer);
462:   /* now copy the desired guys from DBIn => DB */
463:   for (p = 0; p < N; ++p) {
464:     DMSwarmDataBucketCopyPoint(DBIn,list[p], *DB,p);
465:   }
466:   return(0);
467: }

469: /* insert into an exisitng location */
470: PetscErrorCode DMSwarmDataFieldInsertPoint(const DMSwarmDataField field,const PetscInt index,const void *ctx)
471: {

475: #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
476:   /* check point is valid */
477:   if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
478:   if (index >= field->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",field->L);
479: #endif
480:   PetscMemcpy(DMSWARM_DATAFIELD_point_access(field->data,index,field->atomic_size), ctx, field->atomic_size);
481:   return(0);
482: }

484: /* remove data at index - replace with last point */
485: PetscErrorCode DMSwarmDataBucketRemovePointAtIndex(const DMSwarmDataBucket db,const PetscInt index)
486: {
487:   PetscInt       f;
488:   PetscBool      any_active_fields;

492: #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
493:   /* check point is valid */
494:   if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
495:   if (index >= db->allocated) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",db->L+db->buffer);
496: #endif
497:   DMSwarmDataBucketQueryForActiveFields(db,&any_active_fields);
498:   if (any_active_fields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot safely remove point as at least one DMSwarmDataField is currently being accessed");
499:   if (index >= db->L) { /* this point is not in the list - no need to error, but I will anyway */
500:     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);
501:   }
502:   if (index != db->L-1) { /* not last point in list */
503:     for (f = 0; f < db->nfields; ++f) {
504:       DMSwarmDataField field = db->field[f];

506:       /* copy then remove */
507:       DMSwarmDataFieldCopyPoint(db->L-1, field, index, field);
508:       /* DMSwarmDataFieldZeroPoint(field,index); */
509:     }
510:   }
511:   /* decrement size */
512:   /* this will zero out an crap at the end of the list */
513:   DMSwarmDataBucketRemovePoint(db);
514:   return(0);
515: }

517: /* copy x into y */
518: PetscErrorCode DMSwarmDataFieldCopyPoint(const PetscInt pid_x,const DMSwarmDataField field_x,
519:                         const PetscInt pid_y,const DMSwarmDataField field_y )
520: {

524: #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
525:   /* check point is valid */
526:   if (pid_x < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"(IN) index must be >= 0");
527:   if (pid_x >= field_x->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"(IN) index must be < %D",field_x->L);
528:   if (pid_y < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"(OUT) index must be >= 0");
529:   if (pid_y >= field_y->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"(OUT) index must be < %D",field_y->L);
530:   if( field_y->atomic_size != field_x->atomic_size ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"atomic size must match");
531: #endif
532:   PetscMemcpy(DMSWARM_DATAFIELD_point_access(field_y->data,pid_y,field_y->atomic_size),DMSWARM_DATAFIELD_point_access(field_x->data,pid_x,field_x->atomic_size),field_y->atomic_size);
533:   return(0);
534: }


537: /* zero only the datafield at this point */
538: PetscErrorCode DMSwarmDataFieldZeroPoint(const DMSwarmDataField field,const PetscInt index)
539: {

543: #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
544:   /* check point is valid */
545:   if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
546:   if (index >= field->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",field->L);
547: #endif
548:   PetscMemzero(DMSWARM_DATAFIELD_point_access(field->data,index,field->atomic_size), field->atomic_size);
549:   return(0);
550: }

552: /* zero ALL data for this point */
553: PetscErrorCode DMSwarmDataBucketZeroPoint(const DMSwarmDataBucket db,const PetscInt index)
554: {
555:   PetscInt f;

559:   /* check point is valid */
560:   if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
561:   if (index >= db->allocated) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",db->allocated);
562:   for (f = 0; f < db->nfields; ++f) {
563:     DMSwarmDataField field = db->field[f];
564:     DMSwarmDataFieldZeroPoint(field,index);
565:   }
566:   return(0);
567: }

569: /* increment */
570: PetscErrorCode DMSwarmDataBucketAddPoint(DMSwarmDataBucket db)
571: {

575:   DMSwarmDataBucketSetSizes(db,db->L+1,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT);
576:   return(0);
577: }

579: /* decrement */
580: PetscErrorCode DMSwarmDataBucketRemovePoint(DMSwarmDataBucket db)
581: {

585:   DMSwarmDataBucketSetSizes(db,db->L-1,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT);
586:   return(0);
587: }

589: PetscErrorCode DMSwarmDataBucketView_stdout(MPI_Comm comm,DMSwarmDataBucket db)
590: {
591:   PetscInt       f;
592:   double         memory_usage_total,memory_usage_total_local = 0.0;
594: 
596:   PetscPrintf(comm,"DMSwarmDataBucketView: \n");
597:   PetscPrintf(comm,"  L                  = %D \n", db->L);
598:   PetscPrintf(comm,"  buffer             = %D \n", db->buffer);
599:   PetscPrintf(comm,"  allocated          = %D \n", db->allocated);
600:   PetscPrintf(comm,"  nfields registered = %D \n", db->nfields);
601: 
602:   for (f = 0; f < db->nfields; ++f) {
603:     double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
604:     memory_usage_total_local += memory_usage_f;
605:   }
606:   MPI_Allreduce(&memory_usage_total_local,&memory_usage_total,1,MPI_DOUBLE,MPI_SUM,comm);
607: 
608:   for (f = 0; f < db->nfields; ++f) {
609:     double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
610:     PetscPrintf(comm,"    [%3D] %15s : Mem. usage       = %1.2e (MB) [rank0]\n", f, db->field[f]->name, memory_usage_f );
611:     PetscPrintf(comm,"                            blocksize        = %D \n", db->field[f]->bs);
612:     if (db->field[f]->bs != 1) {
613:       PetscPrintf(comm,"                            atomic size      = %zu [full block, bs=%D]\n", db->field[f]->atomic_size,db->field[f]->bs);
614:       PetscPrintf(comm,"                            atomic size/item = %zu \n", db->field[f]->atomic_size/db->field[f]->bs);
615:     } else {
616:       PetscPrintf(comm,"                            atomic size      = %zu \n", db->field[f]->atomic_size);
617:     }
618:   }
619:   PetscPrintf(comm,"  Total mem. usage                           = %1.2e (MB) (collective)\n", memory_usage_total);
620:   return(0);
621: }

623: PetscErrorCode DMSwarmDataBucketView_SEQ(MPI_Comm comm,DMSwarmDataBucket db,const char filename[],DMSwarmDataBucketViewType type)
624: {

628:   switch (type) {
629:   case DATABUCKET_VIEW_STDOUT:
630:     DMSwarmDataBucketView_stdout(PETSC_COMM_SELF,db);
631:     break;
632:   case DATABUCKET_VIEW_ASCII:
633:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for ascii output");
634:     break;
635:   case DATABUCKET_VIEW_BINARY:
636:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for binary output");
637:     break;
638:   case DATABUCKET_VIEW_HDF5:
639:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for HDF5 output");
640:     break;
641:   default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer method requested");
642:   }
643:   return(0);
644: }

646: PetscErrorCode DMSwarmDataBucketView_MPI(MPI_Comm comm,DMSwarmDataBucket db,const char filename[],DMSwarmDataBucketViewType type)
647: {

651:   switch (type) {
652:   case DATABUCKET_VIEW_STDOUT:
653:     DMSwarmDataBucketView_stdout(comm,db);
654:     break;
655:   case DATABUCKET_VIEW_ASCII:
656:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for ascii output");
657:     break;
658:   case DATABUCKET_VIEW_BINARY:
659:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for binary output");
660:     break;
661:   case DATABUCKET_VIEW_HDF5:
662:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for HDF5 output");
663:     break;
664:   default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer method requested");
665:   }
666:   return(0);
667: }

669: PetscErrorCode DMSwarmDataBucketView(MPI_Comm comm,DMSwarmDataBucket db,const char filename[],DMSwarmDataBucketViewType type)
670: {
671:   PetscMPIInt size;

675:   MPI_Comm_size(comm,&size);
676:   if (size == 1) {
677:     DMSwarmDataBucketView_SEQ(comm,db,filename,type);
678:   } else {
679:     DMSwarmDataBucketView_MPI(comm,db,filename,type);
680:   }
681:   return(0);
682: }

684: PetscErrorCode DMSwarmDataBucketDuplicateFields(DMSwarmDataBucket dbA,DMSwarmDataBucket *dbB)
685: {
686:   DMSwarmDataBucket db2;
687:   PetscInt          f;
688:   PetscErrorCode    ierr;

691:   DMSwarmDataBucketCreate(&db2);
692:   /* copy contents from dbA into db2 */
693:   for (f = 0; f < dbA->nfields; ++f) {
694:     DMSwarmDataField field;
695:     size_t    atomic_size;
696:     char      *name;

698:     field = dbA->field[f];
699:     atomic_size = field->atomic_size;
700:     name        = field->name;
701:     DMSwarmDataBucketRegisterField(db2,"DMSwarmDataBucketDuplicateFields",name,atomic_size,NULL);
702:   }
703:   DMSwarmDataBucketFinalize(db2);
704:   DMSwarmDataBucketSetInitialSizes(db2,0,1000);
705:   *dbB = db2;
706:   return(0);
707: }

709: /*
710:  Insert points from db2 into db1
711:  db1 <<== db2
712:  */
713: PetscErrorCode DMSwarmDataBucketInsertValues(DMSwarmDataBucket db1,DMSwarmDataBucket db2)
714: {
715:   PetscInt n_mp_points1,n_mp_points2;
716:   PetscInt n_mp_points1_new,p;

720:   DMSwarmDataBucketGetSizes(db1,&n_mp_points1,0,0);
721:   DMSwarmDataBucketGetSizes(db2,&n_mp_points2,0,0);
722:   n_mp_points1_new = n_mp_points1 + n_mp_points2;
723:   DMSwarmDataBucketSetSizes(db1,n_mp_points1_new,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT);
724:   for (p = 0; p < n_mp_points2; ++p) {
725:     /* db1 <<== db2 */
726:     DMSwarmDataBucketCopyPoint(db2,p, db1,(n_mp_points1 + p));
727:   }
728:   return(0);
729: }

731: /* helpers for parallel send/recv */
732: PetscErrorCode DMSwarmDataBucketCreatePackedArray(DMSwarmDataBucket db,size_t *bytes,void **buf)
733: {
734:   PetscInt       f;
735:   size_t         sizeof_marker_contents;
736:   void          *buffer;

740:   sizeof_marker_contents = 0;
741:   for (f = 0; f < db->nfields; ++f) {
742:     DMSwarmDataField df = db->field[f];
743:     sizeof_marker_contents += df->atomic_size;
744:   }
745:   PetscMalloc(sizeof_marker_contents, &buffer);
746:   PetscMemzero(buffer, sizeof_marker_contents);
747:   if (bytes) {*bytes = sizeof_marker_contents;}
748:   if (buf)   {*buf   = buffer;}
749:   return(0);
750: }

752: PetscErrorCode DMSwarmDataBucketDestroyPackedArray(DMSwarmDataBucket db,void **buf)
753: {

757:   if (buf) {
758:     PetscFree(*buf);
759:     *buf = NULL;
760:   }
761:   return(0);
762: }

764: PetscErrorCode DMSwarmDataBucketFillPackedArray(DMSwarmDataBucket db,const PetscInt index,void *buf)
765: {
766:   PetscInt       f;
767:   void          *data, *data_p;
768:   size_t         asize, offset;

772:   offset = 0;
773:   for (f = 0; f < db->nfields; ++f) {
774:     DMSwarmDataField df = db->field[f];

776:     asize = df->atomic_size;
777:     data = (void*)( df->data );
778:     data_p = (void*)( (char*)data + index*asize );
779:     PetscMemcpy((void*)((char*)buf + offset), data_p, asize);
780:     offset = offset + asize;
781:   }
782:   return(0);
783: }

785: PetscErrorCode DMSwarmDataBucketInsertPackedArray(DMSwarmDataBucket db,const PetscInt idx,void *data)
786: {
787:   PetscInt f;
788:   void *data_p;
789:   size_t offset;

793:   offset = 0;
794:   for (f = 0; f < db->nfields; ++f) {
795:     DMSwarmDataField df = db->field[f];

797:     data_p = (void*)( (char*)data + offset );
798:     DMSwarmDataFieldInsertPoint(df, idx, (void*)data_p);
799:     offset = offset + df->atomic_size;
800:   }
801:   return(0);
802: }