Actual source code: data_bucket.c

petsc-3.9.4 2018-09-11
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 _L,_buffer,_allocated;
327:   PetscInt ierr;

330:   _L = db->L;
331:   _buffer = db->buffer;
332:   _allocated = db->allocated;

334:   if (L) {         MPI_Allreduce(&_L,L,1,MPIU_INT,MPI_SUM,comm); }
335:   if (buffer) {    MPI_Allreduce(&_buffer,buffer,1,MPIU_INT,MPI_SUM,comm); }
336:   if (allocated) { MPI_Allreduce(&_allocated,allocated,1,MPIU_INT,MPI_SUM,comm); }
337:   return(0);
338: }

340: PetscErrorCode DMSwarmDataBucketGetDMSwarmDataFields(DMSwarmDataBucket db,PetscInt *L,DMSwarmDataField *fields[])
341: {
343:   if (L)      {*L      = db->nfields;}
344:   if (fields) {*fields = db->field;}
345:   return(0);
346: }

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

356: PetscErrorCode DMSwarmDataFieldAccessPoint(const DMSwarmDataField gfield,const PetscInt pid,void **ctx_p)
357: {
359:   *ctx_p = NULL;
360: #ifdef DMSWARM_DATAFIELD_POINT_ACCESS_GUARD
361:   /* debug mode */
362:   /* check point is valid */
363:   if (pid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
364:   if (pid >= gfield->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",gfield->L);
365:   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);
366: #endif
367:   *ctx_p = DMSWARM_DATAFIELD_point_access(gfield->data,pid,gfield->atomic_size);
368:   return(0);
369: }

371: PetscErrorCode DMSwarmDataFieldAccessPointOffset(const DMSwarmDataField gfield,const size_t offset,const PetscInt pid,void **ctx_p)
372: {
374: #ifdef DMSWARM_DATAFIELD_POINT_ACCESS_GUARD
375:   /* debug mode */
376:   /* check point is valid */
377:   /* if( offset < 0 ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be >= 0");*/
378:   /* Note compiler realizes this can never happen with an unsigned PetscInt */
379:   if (offset >= gfield->atomic_size) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be < %zu",gfield->atomic_size);
380:   /* check point is valid */
381:   if (pid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
382:   if (pid >= gfield->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",gfield->L);
383:   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);
384: #endif
385:   *ctx_p = DMSWARM_DATAFIELD_point_access_offset(gfield->data,pid,gfield->atomic_size,offset);
386:   return(0);
387: }

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

397: PetscErrorCode DMSwarmDataFieldVerifyAccess(const DMSwarmDataField gfield,const size_t size)
398: {
400: #ifdef DMSWARM_DATAFIELD_POINT_ACCESS_GUARD
401:   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 );
402: #endif
403:   return(0);
404: }

406: PetscErrorCode DMSwarmDataFieldGetAtomicSize(const DMSwarmDataField gfield,size_t *size)
407: {
409:   if (size) {*size = gfield->atomic_size;}
410:   return(0);
411: }

413: PetscErrorCode DMSwarmDataFieldGetEntries(const DMSwarmDataField gfield,void **data)
414: {
416:   if (data) {*data = gfield->data;}
417:   return(0);
418: }

420: PetscErrorCode DMSwarmDataFieldRestoreEntries(const DMSwarmDataField gfield,void **data)
421: {
423:   if (data) {*data = NULL;}
424:   return(0);
425: }

427: /* y = x */
428: PetscErrorCode DMSwarmDataBucketCopyPoint(const DMSwarmDataBucket xb,const PetscInt pid_x,
429:                          const DMSwarmDataBucket yb,const PetscInt pid_y)
430: {
431:   PetscInt f;

435:   for (f = 0; f < xb->nfields; ++f) {
436:     void *dest;
437:     void *src;

439:     DMSwarmDataFieldGetAccess(xb->field[f]);
440:     if (xb != yb) { DMSwarmDataFieldGetAccess( yb->field[f]); }
441:     DMSwarmDataFieldAccessPoint(xb->field[f],pid_x, &src);
442:     DMSwarmDataFieldAccessPoint(yb->field[f],pid_y, &dest);
443:     PetscMemcpy(dest, src, xb->field[f]->atomic_size);
444:     DMSwarmDataFieldRestoreAccess(xb->field[f]);
445:     if (xb != yb) {DMSwarmDataFieldRestoreAccess(yb->field[f]);}
446:   }
447:   return(0);
448: }

450: PetscErrorCode DMSwarmDataBucketCreateFromSubset(DMSwarmDataBucket DBIn,const PetscInt N,const PetscInt list[],DMSwarmDataBucket *DB)
451: {
452:   PetscInt nfields;
453:   DMSwarmDataField *fields;
454:   PetscInt f,L,buffer,allocated,p;

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

474: /* insert into an exisitng location */
475: PetscErrorCode DMSwarmDataFieldInsertPoint(const DMSwarmDataField field,const PetscInt index,const void *ctx)
476: {

480: #ifdef DMSWARM_DATAFIELD_POINT_ACCESS_GUARD
481:   /* check point is valid */
482:   if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
483:   if (index >= field->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",field->L);
484: #endif
485:   PetscMemcpy(DMSWARM_DATAFIELD_point_access(field->data,index,field->atomic_size), ctx, field->atomic_size);
486:   return(0);
487: }

489: /* remove data at index - replace with last point */
490: PetscErrorCode DMSwarmDataBucketRemovePointAtIndex(const DMSwarmDataBucket db,const PetscInt index)
491: {
492:   PetscInt       f;
493:   PetscBool      any_active_fields;

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

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

522: /* copy x into y */
523: PetscErrorCode DMSwarmDataFieldCopyPoint(const PetscInt pid_x,const DMSwarmDataField field_x,
524:                         const PetscInt pid_y,const DMSwarmDataField field_y )
525: {

529: #ifdef DMSWARM_DATAFIELD_POINT_ACCESS_GUARD
530:   /* check point is valid */
531:   if (pid_x < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"(IN) index must be >= 0");
532:   if (pid_x >= field_x->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"(IN) index must be < %D",field_x->L);
533:   if (pid_y < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"(OUT) index must be >= 0");
534:   if (pid_y >= field_y->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"(OUT) index must be < %D",field_y->L);
535:   if( field_y->atomic_size != field_x->atomic_size ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"atomic size must match");
536: #endif
537:   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);
538:   return(0);
539: }


542: /* zero only the datafield at this point */
543: PetscErrorCode DMSwarmDataFieldZeroPoint(const DMSwarmDataField field,const PetscInt index)
544: {

548: #ifdef DMSWARM_DATAFIELD_POINT_ACCESS_GUARD
549:   /* check point is valid */
550:   if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
551:   if (index >= field->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",field->L);
552: #endif
553:   PetscMemzero(DMSWARM_DATAFIELD_point_access(field->data,index,field->atomic_size), field->atomic_size);
554:   return(0);
555: }

557: /* zero ALL data for this point */
558: PetscErrorCode DMSwarmDataBucketZeroPoint(const DMSwarmDataBucket db,const PetscInt index)
559: {
560:   PetscInt f;

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

574: /* increment */
575: PetscErrorCode DMSwarmDataBucketAddPoint(DMSwarmDataBucket db)
576: {

580:   DMSwarmDataBucketSetSizes(db,db->L+1,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT);
581:   return(0);
582: }

584: /* decrement */
585: PetscErrorCode DMSwarmDataBucketRemovePoint(DMSwarmDataBucket db)
586: {

590:   DMSwarmDataBucketSetSizes(db,db->L-1,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT);
591:   return(0);
592: }

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

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

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

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

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

674: PetscErrorCode DMSwarmDataBucketView(MPI_Comm comm,DMSwarmDataBucket db,const char filename[],DMSwarmDataBucketViewType type)
675: {
676:   PetscMPIInt size;

680:   MPI_Comm_size(comm,&size);
681:   if (size == 1) {
682:     DMSwarmDataBucketView_SEQ(comm,db,filename,type);
683:   } else {
684:     DMSwarmDataBucketView_MPI(comm,db,filename,type);
685:   }
686:   return(0);
687: }

689: PetscErrorCode DMSwarmDataBucketDuplicateFields(DMSwarmDataBucket dbA,DMSwarmDataBucket *dbB)
690: {
691:   DMSwarmDataBucket db2;
692:   PetscInt          f;
693:   PetscErrorCode    ierr;

696:   DMSwarmDataBucketCreate(&db2);
697:   /* copy contents from dbA into db2 */
698:   for (f = 0; f < dbA->nfields; ++f) {
699:     DMSwarmDataField field;
700:     size_t    atomic_size;
701:     char      *name;

703:     field = dbA->field[f];
704:     atomic_size = field->atomic_size;
705:     name        = field->name;
706:     DMSwarmDataBucketRegisterField(db2,"DMSwarmDataBucketDuplicateFields",name,atomic_size,NULL);
707:   }
708:   DMSwarmDataBucketFinalize(db2);
709:   DMSwarmDataBucketSetInitialSizes(db2,0,1000);
710:   *dbB = db2;
711:   return(0);
712: }

714: /*
715:  Insert points from db2 into db1
716:  db1 <<== db2
717:  */
718: PetscErrorCode DMSwarmDataBucketInsertValues(DMSwarmDataBucket db1,DMSwarmDataBucket db2)
719: {
720:   PetscInt n_mp_points1,n_mp_points2;
721:   PetscInt n_mp_points1_new,p;

725:   DMSwarmDataBucketGetSizes(db1,&n_mp_points1,0,0);
726:   DMSwarmDataBucketGetSizes(db2,&n_mp_points2,0,0);
727:   n_mp_points1_new = n_mp_points1 + n_mp_points2;
728:   DMSwarmDataBucketSetSizes(db1,n_mp_points1_new,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT);
729:   for (p = 0; p < n_mp_points2; ++p) {
730:     /* db1 <<== db2 */
731:     DMSwarmDataBucketCopyPoint(db2,p, db1,(n_mp_points1 + p));
732:   }
733:   return(0);
734: }

736: /* helpers for parallel send/recv */
737: PetscErrorCode DMSwarmDataBucketCreatePackedArray(DMSwarmDataBucket db,size_t *bytes,void **buf)
738: {
739:   PetscInt       f;
740:   size_t         sizeof_marker_contents;
741:   void          *buffer;

745:   sizeof_marker_contents = 0;
746:   for (f = 0; f < db->nfields; ++f) {
747:     DMSwarmDataField df = db->field[f];
748:     sizeof_marker_contents += df->atomic_size;
749:   }
750:   PetscMalloc(sizeof_marker_contents, &buffer);
751:   PetscMemzero(buffer, sizeof_marker_contents);
752:   if (bytes) {*bytes = sizeof_marker_contents;}
753:   if (buf)   {*buf   = buffer;}
754:   return(0);
755: }

757: PetscErrorCode DMSwarmDataBucketDestroyPackedArray(DMSwarmDataBucket db,void **buf)
758: {

762:   if (buf) {
763:     PetscFree(*buf);
764:     *buf = NULL;
765:   }
766:   return(0);
767: }

769: PetscErrorCode DMSwarmDataBucketFillPackedArray(DMSwarmDataBucket db,const PetscInt index,void *buf)
770: {
771:   PetscInt       f;
772:   void          *data, *data_p;
773:   size_t         asize, offset;

777:   offset = 0;
778:   for (f = 0; f < db->nfields; ++f) {
779:     DMSwarmDataField df = db->field[f];

781:     asize = df->atomic_size;
782:     data = (void*)( df->data );
783:     data_p = (void*)( (char*)data + index*asize );
784:     PetscMemcpy((void*)((char*)buf + offset), data_p, asize);
785:     offset = offset + asize;
786:   }
787:   return(0);
788: }

790: PetscErrorCode DMSwarmDataBucketInsertPackedArray(DMSwarmDataBucket db,const PetscInt idx,void *data)
791: {
792:   PetscInt f;
793:   void *data_p;
794:   size_t offset;

798:   offset = 0;
799:   for (f = 0; f < db->nfields; ++f) {
800:     DMSwarmDataField df = db->field[f];

802:     data_p = (void*)( (char*)data + offset );
803:     DMSwarmDataFieldInsertPoint(df, idx, (void*)data_p);
804:     offset = offset + df->atomic_size;
805:   }
806:   return(0);
807: }