xref: /petsc/src/dm/impls/swarm/data_bucket.c (revision ea78f98c112368f404cd6d4fff6d4dfe73e5a1e7)
1279f676cSBarry Smith #include "../src/dm/impls/swarm/data_bucket.h"
252849c42SDave May 
352849c42SDave May /* string helpers */
477048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldStringInList(const char name[],const PetscInt N,const DMSwarmDataField gfield[],PetscBool *val)
552849c42SDave May {
65c18a9d6SDave May   PetscInt       i;
7521f74f9SMatthew G. Knepley   PetscErrorCode ierr;
852849c42SDave May 
9521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1052849c42SDave May   *val = PETSC_FALSE;
11521f74f9SMatthew G. Knepley   for (i = 0; i < N; ++i) {
12521f74f9SMatthew G. Knepley     PetscBool flg;
13521f74f9SMatthew G. Knepley     ierr = PetscStrcmp(name, gfield[i]->name, &flg);CHKERRQ(ierr);
14521f74f9SMatthew G. Knepley     if (flg) {
1552849c42SDave May       *val = PETSC_TRUE;
162eac95f8SDave May       PetscFunctionReturn(0);
1752849c42SDave May     }
1852849c42SDave May   }
192eac95f8SDave May   PetscFunctionReturn(0);
2052849c42SDave May }
2152849c42SDave May 
2277048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldStringFindInList(const char name[],const PetscInt N,const DMSwarmDataField gfield[],PetscInt *index)
2352849c42SDave May {
245c18a9d6SDave May   PetscInt       i;
25521f74f9SMatthew G. Knepley   PetscErrorCode ierr;
2652849c42SDave May 
27521f74f9SMatthew G. Knepley   PetscFunctionBegin;
2852849c42SDave May   *index = -1;
29521f74f9SMatthew G. Knepley   for (i = 0; i < N; ++i) {
30521f74f9SMatthew G. Knepley     PetscBool flg;
31521f74f9SMatthew G. Knepley     ierr = PetscStrcmp(name, gfield[i]->name, &flg);CHKERRQ(ierr);
32521f74f9SMatthew G. Knepley     if (flg) {
3352849c42SDave May       *index = i;
342eac95f8SDave May       PetscFunctionReturn(0);
3552849c42SDave May     }
3652849c42SDave May   }
372eac95f8SDave May   PetscFunctionReturn(0);
3852849c42SDave May }
3952849c42SDave May 
40ee71fbaeSPatrick Sanan PetscErrorCode DMSwarmDataFieldCreate(const char registration_function[],const char name[],const size_t size,const PetscInt L,DMSwarmDataField *DF)
4152849c42SDave May {
4277048351SPatrick Sanan   DMSwarmDataField df;
43521f74f9SMatthew G. Knepley   PetscErrorCode   ierr;
4452849c42SDave May 
45521f74f9SMatthew G. Knepley   PetscFunctionBegin;
4677048351SPatrick Sanan   ierr = PetscMalloc(sizeof(struct _p_DMSwarmDataField), &df);CHKERRQ(ierr);
4777048351SPatrick Sanan   ierr = PetscMemzero(df, sizeof(struct _p_DMSwarmDataField));CHKERRQ(ierr);
48ee71fbaeSPatrick Sanan   ierr = PetscStrallocpy(registration_function, &df->registration_function);CHKERRQ(ierr);
49521f74f9SMatthew G. Knepley   ierr = PetscStrallocpy(name, &df->name);CHKERRQ(ierr);
5052849c42SDave May   df->atomic_size = size;
5152849c42SDave May   df->L  = L;
5252c3ed93SDave May   df->bs = 1;
53521f74f9SMatthew G. Knepley   /* allocate something so we don't have to reallocate */
54521f74f9SMatthew G. Knepley   ierr = PetscMalloc(size * L, &df->data);CHKERRQ(ierr);
55521f74f9SMatthew G. Knepley   ierr = PetscMemzero(df->data, size * L);CHKERRQ(ierr);
5652849c42SDave May   *DF = df;
572eac95f8SDave May   PetscFunctionReturn(0);
5852849c42SDave May }
5952849c42SDave May 
6077048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldDestroy(DMSwarmDataField *DF)
6152849c42SDave May {
6277048351SPatrick Sanan   DMSwarmDataField df = *DF;
63521f74f9SMatthew G. Knepley   PetscErrorCode   ierr;
6452849c42SDave May 
65521f74f9SMatthew G. Knepley   PetscFunctionBegin;
66ee71fbaeSPatrick Sanan   ierr = PetscFree(df->registration_function);CHKERRQ(ierr);
67521f74f9SMatthew G. Knepley   ierr = PetscFree(df->name);CHKERRQ(ierr);
68521f74f9SMatthew G. Knepley   ierr = PetscFree(df->data);CHKERRQ(ierr);
69521f74f9SMatthew G. Knepley   ierr = PetscFree(df);CHKERRQ(ierr);
7052849c42SDave May   *DF  = NULL;
712eac95f8SDave May   PetscFunctionReturn(0);
7252849c42SDave May }
7352849c42SDave May 
7452849c42SDave May /* data bucket */
7577048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketCreate(DMSwarmDataBucket *DB)
7652849c42SDave May {
7777048351SPatrick Sanan   DMSwarmDataBucket db;
78521f74f9SMatthew G. Knepley   PetscErrorCode    ierr;
7952849c42SDave May 
80521f74f9SMatthew G. Knepley   PetscFunctionBegin;
8177048351SPatrick Sanan   ierr = PetscMalloc(sizeof(struct _p_DMSwarmDataBucket), &db);CHKERRQ(ierr);
8277048351SPatrick Sanan   ierr = PetscMemzero(db, sizeof(struct _p_DMSwarmDataBucket));CHKERRQ(ierr);
8352849c42SDave May 
8452849c42SDave May   db->finalised = PETSC_FALSE;
8552849c42SDave May   /* create empty spaces for fields */
863454631fSDave May   db->L         = -1;
8752849c42SDave May   db->buffer    = 1;
8852849c42SDave May   db->allocated = 1;
8952849c42SDave May   db->nfields   = 0;
90521f74f9SMatthew G. Knepley   ierr = PetscMalloc1(1, &db->field);CHKERRQ(ierr);
9152849c42SDave May   *DB  = db;
922eac95f8SDave May   PetscFunctionReturn(0);
9352849c42SDave May }
9452849c42SDave May 
9577048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketDestroy(DMSwarmDataBucket *DB)
9652849c42SDave May {
9777048351SPatrick Sanan   DMSwarmDataBucket db = *DB;
985c18a9d6SDave May   PetscInt          f;
99dbe06d34SDave May   PetscErrorCode    ierr;
10052849c42SDave May 
101521f74f9SMatthew G. Knepley   PetscFunctionBegin;
10252849c42SDave May   /* release fields */
103521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
10477048351SPatrick Sanan     ierr = DMSwarmDataFieldDestroy(&db->field[f]);CHKERRQ(ierr);
10552849c42SDave May   }
10652849c42SDave May   /* this will catch the initially allocated objects in the event that no fields are registered */
10752849c42SDave May   if (db->field != NULL) {
108521f74f9SMatthew G. Knepley     ierr = PetscFree(db->field);CHKERRQ(ierr);
10952849c42SDave May   }
110521f74f9SMatthew G. Knepley   ierr = PetscFree(db);CHKERRQ(ierr);
11152849c42SDave May   *DB = NULL;
1122eac95f8SDave May   PetscFunctionReturn(0);
11352849c42SDave May }
11452849c42SDave May 
11577048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketQueryForActiveFields(DMSwarmDataBucket db,PetscBool *any_active_fields)
1160cb17e37SDave May {
1170cb17e37SDave May   PetscInt f;
118521f74f9SMatthew G. Knepley 
119521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1200cb17e37SDave May   *any_active_fields = PETSC_FALSE;
121521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
1220cb17e37SDave May     if (db->field[f]->active) {
1230cb17e37SDave May       *any_active_fields = PETSC_TRUE;
1240cb17e37SDave May       PetscFunctionReturn(0);
1250cb17e37SDave May     }
1260cb17e37SDave May   }
1270cb17e37SDave May   PetscFunctionReturn(0);
1280cb17e37SDave May }
1290cb17e37SDave May 
13077048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketRegisterField(
13177048351SPatrick Sanan                               DMSwarmDataBucket db,
132ee71fbaeSPatrick Sanan                               const char registration_function[],
13352849c42SDave May                               const char field_name[],
13477048351SPatrick Sanan                               size_t atomic_size, DMSwarmDataField *_gfield)
13552849c42SDave May {
13652849c42SDave May   PetscBool val;
13777048351SPatrick Sanan   DMSwarmDataField fp;
138dbe06d34SDave May   PetscErrorCode ierr;
13952849c42SDave May 
140521f74f9SMatthew G. Knepley   PetscFunctionBegin;
14152849c42SDave May 	/* check we haven't finalised the registration of fields */
14252849c42SDave May 	/*
14352849c42SDave May    if(db->finalised==PETSC_TRUE) {
14477048351SPatrick Sanan    printf("ERROR: DMSwarmDataBucketFinalize() has been called. Cannot register more fields\n");
14552849c42SDave May    ERROR();
14652849c42SDave May    }
14752849c42SDave May    */
14852849c42SDave May   /* check for repeated name */
14977048351SPatrick Sanan   ierr = DMSwarmDataFieldStringInList(field_name, db->nfields, (const DMSwarmDataField*) db->field, &val);CHKERRQ(ierr);
1502eac95f8SDave May   if (val == PETSC_TRUE) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field %s already exists. Cannot add same field twice",field_name);
15152849c42SDave May   /* create new space for data */
15277048351SPatrick Sanan   ierr = PetscRealloc(sizeof(DMSwarmDataField)*(db->nfields+1), &db->field);CHKERRQ(ierr);
15352849c42SDave May   /* add field */
154ee71fbaeSPatrick Sanan   ierr = DMSwarmDataFieldCreate(registration_function, field_name, atomic_size, db->allocated, &fp);CHKERRQ(ierr);
15552849c42SDave May   db->field[db->nfields] = fp;
15652849c42SDave May   db->nfields++;
15752849c42SDave May   if (_gfield != NULL) {
15852849c42SDave May     *_gfield = fp;
15952849c42SDave May   }
1602eac95f8SDave May   PetscFunctionReturn(0);
16152849c42SDave May }
16252849c42SDave May 
16352849c42SDave May /*
16477048351SPatrick Sanan  #define DMSwarmDataBucketRegisterField(db,name,size,k) {\
16552849c42SDave May  char *location;\
16652849c42SDave May  asprintf(&location,"Registered by %s() at line %d within file %s", __FUNCTION__, __LINE__, __FILE__);\
16777048351SPatrick Sanan  _DMSwarmDataBucketRegisterField( (db), location, (name), (size), (k) );\
168521f74f9SMatthew G. Knepley  ierr = PetscFree(location);\
16952849c42SDave May  }
17052849c42SDave May  */
17152849c42SDave May 
17277048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketGetDMSwarmDataFieldByName(DMSwarmDataBucket db,const char name[],DMSwarmDataField *gfield)
17352849c42SDave May {
1745c18a9d6SDave May   PetscInt       idx;
17552849c42SDave May   PetscBool      found;
176521f74f9SMatthew G. Knepley   PetscErrorCode ierr;
17752849c42SDave May 
178521f74f9SMatthew G. Knepley   PetscFunctionBegin;
17977048351SPatrick Sanan   ierr = DMSwarmDataFieldStringInList(name,db->nfields,(const DMSwarmDataField*)db->field,&found);CHKERRQ(ierr);
18077048351SPatrick Sanan   if (!found) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot find DMSwarmDataField with name %s",name);
18177048351SPatrick Sanan   ierr = DMSwarmDataFieldStringFindInList(name,db->nfields,(const DMSwarmDataField*)db->field,&idx);CHKERRQ(ierr);
18252849c42SDave May   *gfield = db->field[idx];
1832eac95f8SDave May   PetscFunctionReturn(0);
18452849c42SDave May }
18552849c42SDave May 
18677048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketQueryDMSwarmDataFieldByName(DMSwarmDataBucket db,const char name[],PetscBool *found)
18752849c42SDave May {
1882635f519SDave May   PetscErrorCode ierr;
189521f74f9SMatthew G. Knepley 
190521f74f9SMatthew G. Knepley   PetscFunctionBegin;
19152849c42SDave May   *found = PETSC_FALSE;
19277048351SPatrick Sanan   ierr = DMSwarmDataFieldStringInList(name,db->nfields,(const DMSwarmDataField*)db->field,found);CHKERRQ(ierr);
1937fbf63aeSDave May   PetscFunctionReturn(0);
19452849c42SDave May }
19552849c42SDave May 
19677048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketFinalize(DMSwarmDataBucket db)
19752849c42SDave May {
198521f74f9SMatthew G. Knepley   PetscFunctionBegin;
19952849c42SDave May   db->finalised = PETSC_TRUE;
2007fbf63aeSDave May   PetscFunctionReturn(0);
20152849c42SDave May }
20252849c42SDave May 
20377048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldGetNumEntries(DMSwarmDataField df,PetscInt *sum)
20452849c42SDave May {
205521f74f9SMatthew G. Knepley   PetscFunctionBegin;
20652849c42SDave May   *sum = df->L;
2077fbf63aeSDave May   PetscFunctionReturn(0);
20852849c42SDave May }
20952849c42SDave May 
21077048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldSetBlockSize(DMSwarmDataField df,PetscInt blocksize)
21152c3ed93SDave May {
212521f74f9SMatthew G. Knepley   PetscFunctionBegin;
21352c3ed93SDave May   df->bs = blocksize;
21452c3ed93SDave May   PetscFunctionReturn(0);
21552c3ed93SDave May }
21652c3ed93SDave May 
21777048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldSetSize(DMSwarmDataField df,const PetscInt new_L)
21852849c42SDave May {
219521f74f9SMatthew G. Knepley   PetscErrorCode ierr;
22052849c42SDave May 
221521f74f9SMatthew G. Knepley   PetscFunctionBegin;
22277048351SPatrick Sanan   if (new_L < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot set size of DMSwarmDataField to be < 0");
2237fbf63aeSDave May   if (new_L == df->L) PetscFunctionReturn(0);
22452849c42SDave May   if (new_L > df->L) {
2254be7464cSMatthew G. Knepley     ierr = PetscRealloc(df->atomic_size * (new_L), &df->data);CHKERRQ(ierr);
22652849c42SDave May     /* init new contents */
227521f74f9SMatthew G. Knepley     ierr = PetscMemzero(( ((char*)df->data)+df->L*df->atomic_size), (new_L-df->L)*df->atomic_size);CHKERRQ(ierr);
2287fbf63aeSDave May   } else {
22952849c42SDave May     /* reallocate pointer list, add +1 in case new_L = 0 */
2304be7464cSMatthew G. Knepley     ierr = PetscRealloc(df->atomic_size * (new_L+1), &df->data);CHKERRQ(ierr);
23152849c42SDave May   }
23252849c42SDave May   df->L = new_L;
2337fbf63aeSDave May   PetscFunctionReturn(0);
23452849c42SDave May }
23552849c42SDave May 
23677048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldZeroBlock(DMSwarmDataField df,const PetscInt start,const PetscInt end)
23752849c42SDave May {
238521f74f9SMatthew G. Knepley   PetscErrorCode ierr;
239521f74f9SMatthew G. Knepley 
240521f74f9SMatthew G. Knepley   PetscFunctionBegin;
2417fbf63aeSDave May   if (start > end) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if start(%D) > end(%D)",start,end);
2427fbf63aeSDave May   if (start < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if start(%D) < 0",start);
243a233d522SDave May   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);
244521f74f9SMatthew G. Knepley   ierr = PetscMemzero((((char*)df->data)+start*df->atomic_size), (end-start)*df->atomic_size);CHKERRQ(ierr);
2457fbf63aeSDave May   PetscFunctionReturn(0);
24652849c42SDave May }
24752849c42SDave May 
24852849c42SDave May /*
24952849c42SDave May  A negative buffer value will simply be ignored and the old buffer value will be used.
25052849c42SDave May  */
25177048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketSetSizes(DMSwarmDataBucket db,const PetscInt L,const PetscInt buffer)
25252849c42SDave May {
2535c18a9d6SDave May   PetscInt       current_allocated,new_used,new_unused,new_buffer,new_allocated,f;
2540cb17e37SDave May   PetscBool      any_active_fields;
255dbe06d34SDave May   PetscErrorCode ierr;
25652849c42SDave May 
257521f74f9SMatthew G. Knepley   PetscFunctionBegin;
25877048351SPatrick Sanan   if (db->finalised == PETSC_FALSE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"You must call DMSwarmDataBucketFinalize() before DMSwarmDataBucketSetSizes()");
25977048351SPatrick Sanan   ierr = DMSwarmDataBucketQueryForActiveFields(db,&any_active_fields);CHKERRQ(ierr);
26077048351SPatrick Sanan   if (any_active_fields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot safely re-size as at least one DMSwarmDataField is currently being accessed");
2610cb17e37SDave May 
26252849c42SDave May   current_allocated = db->allocated;
26352849c42SDave May   new_used   = L;
26452849c42SDave May   new_unused = current_allocated - new_used;
26552849c42SDave May   new_buffer = db->buffer;
26652849c42SDave May   if (buffer >= 0) { /* update the buffer value */
26752849c42SDave May     new_buffer = buffer;
26852849c42SDave May   }
26952849c42SDave May   new_allocated = new_used + new_buffer;
27052849c42SDave May   /* action */
27152849c42SDave May   if (new_allocated > current_allocated) {
27252849c42SDave May     /* increase size to new_used + new_buffer */
27352849c42SDave May     for (f=0; f<db->nfields; f++) {
27477048351SPatrick Sanan       ierr = DMSwarmDataFieldSetSize(db->field[f], new_allocated);CHKERRQ(ierr);
27552849c42SDave May     }
27652849c42SDave May     db->L         = new_used;
27752849c42SDave May     db->buffer    = new_buffer;
27852849c42SDave May     db->allocated = new_used + new_buffer;
279521f74f9SMatthew G. Knepley   } else {
28052849c42SDave May     if (new_unused > 2 * new_buffer) {
28152849c42SDave May       /* shrink array to new_used + new_buffer */
282521f74f9SMatthew G. Knepley       for (f = 0; f < db->nfields; ++f) {
28377048351SPatrick Sanan         ierr = DMSwarmDataFieldSetSize(db->field[f], new_allocated);CHKERRQ(ierr);
28452849c42SDave May       }
28552849c42SDave May       db->L         = new_used;
28652849c42SDave May       db->buffer    = new_buffer;
28752849c42SDave May       db->allocated = new_used + new_buffer;
288521f74f9SMatthew G. Knepley     } else {
28952849c42SDave May       db->L      = new_used;
29052849c42SDave May       db->buffer = new_buffer;
29152849c42SDave May     }
29252849c42SDave May   }
29352849c42SDave May   /* zero all entries from db->L to db->allocated */
294521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
29577048351SPatrick Sanan     DMSwarmDataField field = db->field[f];
29677048351SPatrick Sanan     ierr = DMSwarmDataFieldZeroBlock(field, db->L,db->allocated);CHKERRQ(ierr);
29752849c42SDave May   }
2987fbf63aeSDave May   PetscFunctionReturn(0);
29952849c42SDave May }
30052849c42SDave May 
30177048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketSetInitialSizes(DMSwarmDataBucket db,const PetscInt L,const PetscInt buffer)
30252849c42SDave May {
3035c18a9d6SDave May   PetscInt       f;
304dbe06d34SDave May   PetscErrorCode ierr;
305dbe06d34SDave May 
306521f74f9SMatthew G. Knepley   PetscFunctionBegin;
30777048351SPatrick Sanan   ierr = DMSwarmDataBucketSetSizes(db,L,buffer);CHKERRQ(ierr);
308521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
30977048351SPatrick Sanan     DMSwarmDataField field = db->field[f];
31077048351SPatrick Sanan     ierr = DMSwarmDataFieldZeroBlock(field,0,db->allocated);CHKERRQ(ierr);
31152849c42SDave May   }
3127fbf63aeSDave May   PetscFunctionReturn(0);
31352849c42SDave May }
31452849c42SDave May 
31577048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketGetSizes(DMSwarmDataBucket db,PetscInt *L,PetscInt *buffer,PetscInt *allocated)
31652849c42SDave May {
317521f74f9SMatthew G. Knepley   PetscFunctionBegin;
31852849c42SDave May   if (L) {*L = db->L;}
31952849c42SDave May   if (buffer) {*buffer = db->buffer;}
32052849c42SDave May   if (allocated) {*allocated = db->allocated;}
3217fbf63aeSDave May   PetscFunctionReturn(0);
32252849c42SDave May }
32352849c42SDave May 
32477048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketGetGlobalSizes(MPI_Comm comm,DMSwarmDataBucket db,PetscInt *L,PetscInt *buffer,PetscInt *allocated)
32552849c42SDave May {
3265c18a9d6SDave May   PetscInt ierr;
32752849c42SDave May 
328521f74f9SMatthew G. Knepley   PetscFunctionBegin;
329df4bb746SBarry Smith   if (L) {         ierr = MPI_Allreduce(&db->L,L,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); }
330df4bb746SBarry Smith   if (buffer) {    ierr = MPI_Allreduce(&db->buffer,buffer,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); }
331df4bb746SBarry Smith   if (allocated) { ierr = MPI_Allreduce(&db->allocated,allocated,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); }
3327fbf63aeSDave May   PetscFunctionReturn(0);
33352849c42SDave May }
33452849c42SDave May 
33577048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketGetDMSwarmDataFields(DMSwarmDataBucket db,PetscInt *L,DMSwarmDataField *fields[])
33652849c42SDave May {
337521f74f9SMatthew G. Knepley   PetscFunctionBegin;
33852849c42SDave May   if (L)      {*L      = db->nfields;}
33952849c42SDave May   if (fields) {*fields = db->field;}
3407fbf63aeSDave May   PetscFunctionReturn(0);
34152849c42SDave May }
34252849c42SDave May 
34377048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldGetAccess(const DMSwarmDataField gfield)
34452849c42SDave May {
345521f74f9SMatthew G. Knepley   PetscFunctionBegin;
34677048351SPatrick Sanan   if (gfield->active) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is already active. You must call DMSwarmDataFieldRestoreAccess()",gfield->name);
34752849c42SDave May   gfield->active = PETSC_TRUE;
3487fbf63aeSDave May   PetscFunctionReturn(0);
34952849c42SDave May }
35052849c42SDave May 
35177048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldAccessPoint(const DMSwarmDataField gfield,const PetscInt pid,void **ctx_p)
35252849c42SDave May {
353521f74f9SMatthew G. Knepley   PetscFunctionBegin;
35472505a4dSJed Brown   *ctx_p = NULL;
355497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
35652849c42SDave May   /* debug mode */
35784bcda08SDave May   /* check point is valid */
3587fbf63aeSDave May   if (pid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
3597fbf63aeSDave May   if (pid >= gfield->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",gfield->L);
36077048351SPatrick Sanan   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);
36152849c42SDave May #endif
36277048351SPatrick Sanan   *ctx_p = DMSWARM_DATAFIELD_point_access(gfield->data,pid,gfield->atomic_size);
3637fbf63aeSDave May   PetscFunctionReturn(0);
36452849c42SDave May }
36552849c42SDave May 
36677048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldAccessPointOffset(const DMSwarmDataField gfield,const size_t offset,const PetscInt pid,void **ctx_p)
36752849c42SDave May {
368521f74f9SMatthew G. Knepley   PetscFunctionBegin;
369497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
37052849c42SDave May   /* debug mode */
37184bcda08SDave May   /* check point is valid */
372521f74f9SMatthew G. Knepley   /* if( offset < 0 ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be >= 0");*/
373521f74f9SMatthew G. Knepley   /* Note compiler realizes this can never happen with an unsigned PetscInt */
3747fbf63aeSDave May   if (offset >= gfield->atomic_size) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be < %zu",gfield->atomic_size);
37584bcda08SDave May   /* check point is valid */
3767fbf63aeSDave May   if (pid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
377a233d522SDave May   if (pid >= gfield->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",gfield->L);
37877048351SPatrick Sanan   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);
37952849c42SDave May #endif
38077048351SPatrick Sanan   *ctx_p = DMSWARM_DATAFIELD_point_access_offset(gfield->data,pid,gfield->atomic_size,offset);
3817fbf63aeSDave May   PetscFunctionReturn(0);
38252849c42SDave May }
38352849c42SDave May 
38477048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldRestoreAccess(DMSwarmDataField gfield)
38552849c42SDave May {
386521f74f9SMatthew G. Knepley   PetscFunctionBegin;
38777048351SPatrick Sanan   if (gfield->active == PETSC_FALSE) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is not active. You must call DMSwarmDataFieldGetAccess()", gfield->name);
38852849c42SDave May   gfield->active = PETSC_FALSE;
3897fbf63aeSDave May   PetscFunctionReturn(0);
39052849c42SDave May }
39152849c42SDave May 
39277048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldVerifyAccess(const DMSwarmDataField gfield,const size_t size)
39352849c42SDave May {
394521f74f9SMatthew G. Knepley   PetscFunctionBegin;
395497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
396521f74f9SMatthew G. Knepley   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 );
39752849c42SDave May #endif
3987fbf63aeSDave May   PetscFunctionReturn(0);
39952849c42SDave May }
40052849c42SDave May 
40177048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldGetAtomicSize(const DMSwarmDataField gfield,size_t *size)
40252849c42SDave May {
403521f74f9SMatthew G. Knepley   PetscFunctionBegin;
40452849c42SDave May   if (size) {*size = gfield->atomic_size;}
4057fbf63aeSDave May   PetscFunctionReturn(0);
40652849c42SDave May }
40752849c42SDave May 
40877048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldGetEntries(const DMSwarmDataField gfield,void **data)
40952849c42SDave May {
410521f74f9SMatthew G. Knepley   PetscFunctionBegin;
411521f74f9SMatthew G. Knepley   if (data) {*data = gfield->data;}
4127fbf63aeSDave May   PetscFunctionReturn(0);
41352849c42SDave May }
41452849c42SDave May 
41577048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldRestoreEntries(const DMSwarmDataField gfield,void **data)
41652849c42SDave May {
417521f74f9SMatthew G. Knepley   PetscFunctionBegin;
418521f74f9SMatthew G. Knepley   if (data) {*data = NULL;}
4197fbf63aeSDave May   PetscFunctionReturn(0);
42052849c42SDave May }
42152849c42SDave May 
42252849c42SDave May /* y = x */
42377048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketCopyPoint(const DMSwarmDataBucket xb,const PetscInt pid_x,
42477048351SPatrick Sanan                          const DMSwarmDataBucket yb,const PetscInt pid_y)
42552849c42SDave May {
4265c18a9d6SDave May   PetscInt f;
427dbe06d34SDave May   PetscErrorCode ierr;
428dbe06d34SDave May 
429521f74f9SMatthew G. Knepley   PetscFunctionBegin;
430521f74f9SMatthew G. Knepley   for (f = 0; f < xb->nfields; ++f) {
43152849c42SDave May     void *dest;
43252849c42SDave May     void *src;
43352849c42SDave May 
43477048351SPatrick Sanan     ierr = DMSwarmDataFieldGetAccess(xb->field[f]);CHKERRQ(ierr);
43577048351SPatrick Sanan     if (xb != yb) { ierr = DMSwarmDataFieldGetAccess( yb->field[f]);CHKERRQ(ierr); }
43677048351SPatrick Sanan     ierr = DMSwarmDataFieldAccessPoint(xb->field[f],pid_x, &src);CHKERRQ(ierr);
43777048351SPatrick Sanan     ierr = DMSwarmDataFieldAccessPoint(yb->field[f],pid_y, &dest);CHKERRQ(ierr);
438521f74f9SMatthew G. Knepley     ierr = PetscMemcpy(dest, src, xb->field[f]->atomic_size);CHKERRQ(ierr);
43977048351SPatrick Sanan     ierr = DMSwarmDataFieldRestoreAccess(xb->field[f]);CHKERRQ(ierr);
44077048351SPatrick Sanan     if (xb != yb) {ierr = DMSwarmDataFieldRestoreAccess(yb->field[f]);CHKERRQ(ierr);}
44152849c42SDave May   }
4427fbf63aeSDave May   PetscFunctionReturn(0);
44352849c42SDave May }
44452849c42SDave May 
44577048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketCreateFromSubset(DMSwarmDataBucket DBIn,const PetscInt N,const PetscInt list[],DMSwarmDataBucket *DB)
44652849c42SDave May {
4475c18a9d6SDave May   PetscInt nfields;
44877048351SPatrick Sanan   DMSwarmDataField *fields;
4495c18a9d6SDave May   PetscInt f,L,buffer,allocated,p;
450dbe06d34SDave May   PetscErrorCode ierr;
45152849c42SDave May 
452521f74f9SMatthew G. Knepley   PetscFunctionBegin;
45377048351SPatrick Sanan   ierr = DMSwarmDataBucketCreate(DB);CHKERRQ(ierr);
45452849c42SDave May   /* copy contents of DBIn */
45577048351SPatrick Sanan   ierr = DMSwarmDataBucketGetDMSwarmDataFields(DBIn,&nfields,&fields);CHKERRQ(ierr);
45677048351SPatrick Sanan   ierr = DMSwarmDataBucketGetSizes(DBIn,&L,&buffer,&allocated);CHKERRQ(ierr);
457521f74f9SMatthew G. Knepley   for (f = 0; f < nfields; ++f) {
45877048351SPatrick Sanan     ierr = DMSwarmDataBucketRegisterField(*DB,"DMSwarmDataBucketCreateFromSubset",fields[f]->name,fields[f]->atomic_size,NULL);CHKERRQ(ierr);
45952849c42SDave May   }
46077048351SPatrick Sanan   ierr = DMSwarmDataBucketFinalize(*DB);CHKERRQ(ierr);
46177048351SPatrick Sanan   ierr = DMSwarmDataBucketSetSizes(*DB,L,buffer);CHKERRQ(ierr);
46252849c42SDave May   /* now copy the desired guys from DBIn => DB */
463521f74f9SMatthew G. Knepley   for (p = 0; p < N; ++p) {
46477048351SPatrick Sanan     ierr = DMSwarmDataBucketCopyPoint(DBIn,list[p], *DB,p);CHKERRQ(ierr);
46552849c42SDave May   }
4667fbf63aeSDave May   PetscFunctionReturn(0);
46752849c42SDave May }
46852849c42SDave May 
469521f74f9SMatthew G. Knepley /* insert into an exisitng location */
47077048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldInsertPoint(const DMSwarmDataField field,const PetscInt index,const void *ctx)
47152849c42SDave May {
472521f74f9SMatthew G. Knepley   PetscErrorCode ierr;
47352849c42SDave May 
474521f74f9SMatthew G. Knepley   PetscFunctionBegin;
475497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
47684bcda08SDave May   /* check point is valid */
477a233d522SDave May   if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
478a233d522SDave May   if (index >= field->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",field->L);
47952849c42SDave May #endif
48077048351SPatrick Sanan   ierr = PetscMemcpy(DMSWARM_DATAFIELD_point_access(field->data,index,field->atomic_size), ctx, field->atomic_size);CHKERRQ(ierr);
4817fbf63aeSDave May   PetscFunctionReturn(0);
48252849c42SDave May }
48352849c42SDave May 
484521f74f9SMatthew G. Knepley /* remove data at index - replace with last point */
48577048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketRemovePointAtIndex(const DMSwarmDataBucket db,const PetscInt index)
48652849c42SDave May {
4875c18a9d6SDave May   PetscInt       f;
4880cb17e37SDave May   PetscBool      any_active_fields;
489dbe06d34SDave May   PetscErrorCode ierr;
49052849c42SDave May 
491521f74f9SMatthew G. Knepley   PetscFunctionBegin;
492497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
49384bcda08SDave May   /* check point is valid */
494a233d522SDave May   if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
495a233d522SDave May   if (index >= db->allocated) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",db->L+db->buffer);
49652849c42SDave May #endif
49777048351SPatrick Sanan   ierr = DMSwarmDataBucketQueryForActiveFields(db,&any_active_fields);CHKERRQ(ierr);
49877048351SPatrick Sanan   if (any_active_fields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot safely remove point as at least one DMSwarmDataField is currently being accessed");
499a233d522SDave May   if (index >= db->L) { /* this point is not in the list - no need to error, but I will anyway */
500a233d522SDave May     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);
50152849c42SDave May   }
502a233d522SDave May   if (index != db->L-1) { /* not last point in list */
503521f74f9SMatthew G. Knepley     for (f = 0; f < db->nfields; ++f) {
50477048351SPatrick Sanan       DMSwarmDataField field = db->field[f];
50552849c42SDave May 
50652849c42SDave May       /* copy then remove */
50777048351SPatrick Sanan       ierr = DMSwarmDataFieldCopyPoint(db->L-1, field, index, field);CHKERRQ(ierr);
50877048351SPatrick Sanan       /* DMSwarmDataFieldZeroPoint(field,index); */
50952849c42SDave May     }
51052849c42SDave May   }
51152849c42SDave May   /* decrement size */
51252849c42SDave May   /* this will zero out an crap at the end of the list */
51377048351SPatrick Sanan   ierr = DMSwarmDataBucketRemovePoint(db);CHKERRQ(ierr);
5147fbf63aeSDave May   PetscFunctionReturn(0);
51552849c42SDave May }
51652849c42SDave May 
51752849c42SDave May /* copy x into y */
51877048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldCopyPoint(const PetscInt pid_x,const DMSwarmDataField field_x,
51977048351SPatrick Sanan                         const PetscInt pid_y,const DMSwarmDataField field_y )
52052849c42SDave May {
521521f74f9SMatthew G. Knepley   PetscErrorCode ierr;
52252849c42SDave May 
523521f74f9SMatthew G. Knepley   PetscFunctionBegin;
524497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
52584bcda08SDave May   /* check point is valid */
526a233d522SDave May   if (pid_x < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"(IN) index must be >= 0");
527a233d522SDave May   if (pid_x >= field_x->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"(IN) index must be < %D",field_x->L);
528a233d522SDave May   if (pid_y < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"(OUT) index must be >= 0");
529a233d522SDave May   if (pid_y >= field_y->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"(OUT) index must be < %D",field_y->L);
530521f74f9SMatthew G. Knepley   if( field_y->atomic_size != field_x->atomic_size ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"atomic size must match");
53152849c42SDave May #endif
53277048351SPatrick Sanan   ierr = 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);CHKERRQ(ierr);
5337fbf63aeSDave May   PetscFunctionReturn(0);
53452849c42SDave May }
53552849c42SDave May 
53652849c42SDave May 
537521f74f9SMatthew G. Knepley /* zero only the datafield at this point */
53877048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldZeroPoint(const DMSwarmDataField field,const PetscInt index)
53952849c42SDave May {
540521f74f9SMatthew G. Knepley   PetscErrorCode ierr;
541521f74f9SMatthew G. Knepley 
542521f74f9SMatthew G. Knepley   PetscFunctionBegin;
543497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
54484bcda08SDave May   /* check point is valid */
545a233d522SDave May   if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
546a233d522SDave May   if (index >= field->L) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",field->L);
54752849c42SDave May #endif
54877048351SPatrick Sanan   ierr = PetscMemzero(DMSWARM_DATAFIELD_point_access(field->data,index,field->atomic_size), field->atomic_size);CHKERRQ(ierr);
5497fbf63aeSDave May   PetscFunctionReturn(0);
55052849c42SDave May }
55152849c42SDave May 
552521f74f9SMatthew G. Knepley /* zero ALL data for this point */
55377048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketZeroPoint(const DMSwarmDataBucket db,const PetscInt index)
55452849c42SDave May {
5555c18a9d6SDave May   PetscInt f;
556dbe06d34SDave May   PetscErrorCode ierr;
55752849c42SDave May 
558521f74f9SMatthew G. Knepley   PetscFunctionBegin;
55984bcda08SDave May   /* check point is valid */
560a233d522SDave May   if (index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
561a233d522SDave May   if (index >= db->allocated) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",db->allocated);
562521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
56377048351SPatrick Sanan     DMSwarmDataField field = db->field[f];
56477048351SPatrick Sanan     ierr = DMSwarmDataFieldZeroPoint(field,index);CHKERRQ(ierr);
56552849c42SDave May   }
5667fbf63aeSDave May   PetscFunctionReturn(0);
56752849c42SDave May }
56852849c42SDave May 
56952849c42SDave May /* increment */
57077048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketAddPoint(DMSwarmDataBucket db)
57152849c42SDave May {
572dbe06d34SDave May   PetscErrorCode ierr;
573dbe06d34SDave May 
574521f74f9SMatthew G. Knepley   PetscFunctionBegin;
57577048351SPatrick Sanan   ierr = DMSwarmDataBucketSetSizes(db,db->L+1,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr);
5767fbf63aeSDave May   PetscFunctionReturn(0);
57752849c42SDave May }
57852849c42SDave May 
5797fbf63aeSDave May /* decrement */
58077048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketRemovePoint(DMSwarmDataBucket db)
5817fbf63aeSDave May {
582dbe06d34SDave May   PetscErrorCode ierr;
583dbe06d34SDave May 
584521f74f9SMatthew G. Knepley   PetscFunctionBegin;
58577048351SPatrick Sanan   ierr = DMSwarmDataBucketSetSizes(db,db->L-1,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr);
5867fbf63aeSDave May   PetscFunctionReturn(0);
5877fbf63aeSDave May }
5887fbf63aeSDave May 
58977048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketView_stdout(MPI_Comm comm,DMSwarmDataBucket db)
590ab2be9a3SDave May {
591ab2be9a3SDave May   PetscInt       f;
592ab2be9a3SDave May   double         memory_usage_total,memory_usage_total_local = 0.0;
593ab2be9a3SDave May   PetscErrorCode ierr;
594ab2be9a3SDave May 
595ab2be9a3SDave May   PetscFunctionBegin;
59677048351SPatrick Sanan   ierr = PetscPrintf(comm,"DMSwarmDataBucketView: \n");CHKERRQ(ierr);
597ab2be9a3SDave May   ierr = PetscPrintf(comm,"  L                  = %D \n", db->L);CHKERRQ(ierr);
598ab2be9a3SDave May   ierr = PetscPrintf(comm,"  buffer             = %D \n", db->buffer);CHKERRQ(ierr);
599ab2be9a3SDave May   ierr = PetscPrintf(comm,"  allocated          = %D \n", db->allocated);CHKERRQ(ierr);
600ab2be9a3SDave May   ierr = PetscPrintf(comm,"  nfields registered = %D \n", db->nfields);CHKERRQ(ierr);
601ab2be9a3SDave May 
602ab2be9a3SDave May   for (f = 0; f < db->nfields; ++f) {
603ab2be9a3SDave May     double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
604ab2be9a3SDave May     memory_usage_total_local += memory_usage_f;
605ab2be9a3SDave May   }
606ab2be9a3SDave May   ierr = MPI_Allreduce(&memory_usage_total_local,&memory_usage_total,1,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
607ab2be9a3SDave May 
608ab2be9a3SDave May   for (f = 0; f < db->nfields; ++f) {
609ab2be9a3SDave May     double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
610ab2be9a3SDave May     ierr = PetscPrintf(comm,"    [%3D] %15s : Mem. usage       = %1.2e (MB) [rank0]\n", f, db->field[f]->name, memory_usage_f );CHKERRQ(ierr);
611ab2be9a3SDave May     ierr = PetscPrintf(comm,"                            blocksize        = %D \n", db->field[f]->bs);CHKERRQ(ierr);
612ab2be9a3SDave May     if (db->field[f]->bs != 1) {
613ab2be9a3SDave May       ierr = PetscPrintf(comm,"                            atomic size      = %zu [full block, bs=%D]\n", db->field[f]->atomic_size,db->field[f]->bs);CHKERRQ(ierr);
614ab2be9a3SDave May       ierr = PetscPrintf(comm,"                            atomic size/item = %zu \n", db->field[f]->atomic_size/db->field[f]->bs);CHKERRQ(ierr);
615ab2be9a3SDave May     } else {
616ab2be9a3SDave May       ierr = PetscPrintf(comm,"                            atomic size      = %zu \n", db->field[f]->atomic_size);CHKERRQ(ierr);
617ab2be9a3SDave May     }
618ab2be9a3SDave May   }
619ab2be9a3SDave May   ierr = PetscPrintf(comm,"  Total mem. usage                           = %1.2e (MB) (collective)\n", memory_usage_total);CHKERRQ(ierr);
620ab2be9a3SDave May   PetscFunctionReturn(0);
621ab2be9a3SDave May }
622ab2be9a3SDave May 
62377048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketView_SEQ(MPI_Comm comm,DMSwarmDataBucket db,const char filename[],DMSwarmDataBucketViewType type)
62452849c42SDave May {
625dbe06d34SDave May   PetscErrorCode ierr;
626dbe06d34SDave May 
627521f74f9SMatthew G. Knepley   PetscFunctionBegin;
62852849c42SDave May   switch (type) {
62952849c42SDave May   case DATABUCKET_VIEW_STDOUT:
63077048351SPatrick Sanan     ierr = DMSwarmDataBucketView_stdout(PETSC_COMM_SELF,db);CHKERRQ(ierr);
63152849c42SDave May     break;
63252849c42SDave May   case DATABUCKET_VIEW_ASCII:
633071ce369SDave May     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for ascii output");
63452849c42SDave May     break;
63552849c42SDave May   case DATABUCKET_VIEW_BINARY:
636071ce369SDave May     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for binary output");
63752849c42SDave May     break;
63852849c42SDave May   case DATABUCKET_VIEW_HDF5:
639071ce369SDave May     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for HDF5 output");
64052849c42SDave May     break;
641521f74f9SMatthew G. Knepley   default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer method requested");
64252849c42SDave May   }
6437fbf63aeSDave May   PetscFunctionReturn(0);
64452849c42SDave May }
64552849c42SDave May 
64677048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketView_MPI(MPI_Comm comm,DMSwarmDataBucket db,const char filename[],DMSwarmDataBucketViewType type)
64752849c42SDave May {
648997fa542SDave May   PetscErrorCode ierr;
649997fa542SDave May 
650521f74f9SMatthew G. Knepley   PetscFunctionBegin;
65152849c42SDave May   switch (type) {
65252849c42SDave May   case DATABUCKET_VIEW_STDOUT:
65377048351SPatrick Sanan     ierr = DMSwarmDataBucketView_stdout(comm,db);CHKERRQ(ierr);
65452849c42SDave May     break;
65552849c42SDave May   case DATABUCKET_VIEW_ASCII:
656071ce369SDave May     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for ascii output");
65752849c42SDave May     break;
65852849c42SDave May   case DATABUCKET_VIEW_BINARY:
659071ce369SDave May     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for binary output");
66052849c42SDave May     break;
66152849c42SDave May   case DATABUCKET_VIEW_HDF5:
662071ce369SDave May     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for HDF5 output");
66352849c42SDave May     break;
664521f74f9SMatthew G. Knepley   default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer method requested");
66552849c42SDave May   }
6667fbf63aeSDave May   PetscFunctionReturn(0);
66752849c42SDave May }
66852849c42SDave May 
66977048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketView(MPI_Comm comm,DMSwarmDataBucket db,const char filename[],DMSwarmDataBucketViewType type)
67052849c42SDave May {
671d7d19db6SBarry Smith   PetscMPIInt size;
672997fa542SDave May   PetscErrorCode ierr;
67352849c42SDave May 
674521f74f9SMatthew G. Knepley   PetscFunctionBegin;
675d7d19db6SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
676d7d19db6SBarry Smith   if (size == 1) {
67777048351SPatrick Sanan     ierr = DMSwarmDataBucketView_SEQ(comm,db,filename,type);CHKERRQ(ierr);
67852849c42SDave May   } else {
67977048351SPatrick Sanan     ierr = DMSwarmDataBucketView_MPI(comm,db,filename,type);CHKERRQ(ierr);
68052849c42SDave May   }
6817fbf63aeSDave May   PetscFunctionReturn(0);
68252849c42SDave May }
68352849c42SDave May 
68477048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketDuplicateFields(DMSwarmDataBucket dbA,DMSwarmDataBucket *dbB)
68552849c42SDave May {
68677048351SPatrick Sanan   DMSwarmDataBucket db2;
6875c18a9d6SDave May   PetscInt          f;
688dbe06d34SDave May   PetscErrorCode    ierr;
68952849c42SDave May 
690521f74f9SMatthew G. Knepley   PetscFunctionBegin;
69177048351SPatrick Sanan   ierr = DMSwarmDataBucketCreate(&db2);CHKERRQ(ierr);
69252849c42SDave May   /* copy contents from dbA into db2 */
693521f74f9SMatthew G. Knepley   for (f = 0; f < dbA->nfields; ++f) {
69477048351SPatrick Sanan     DMSwarmDataField field;
69552849c42SDave May     size_t    atomic_size;
69652849c42SDave May     char      *name;
69752849c42SDave May 
69852849c42SDave May     field = dbA->field[f];
69952849c42SDave May     atomic_size = field->atomic_size;
70052849c42SDave May     name        = field->name;
70177048351SPatrick Sanan     ierr = DMSwarmDataBucketRegisterField(db2,"DMSwarmDataBucketDuplicateFields",name,atomic_size,NULL);CHKERRQ(ierr);
70252849c42SDave May   }
70377048351SPatrick Sanan   ierr = DMSwarmDataBucketFinalize(db2);CHKERRQ(ierr);
70477048351SPatrick Sanan   ierr = DMSwarmDataBucketSetInitialSizes(db2,0,1000);CHKERRQ(ierr);
70552849c42SDave May   *dbB = db2;
7067fbf63aeSDave May   PetscFunctionReturn(0);
70752849c42SDave May }
70852849c42SDave May 
70952849c42SDave May /*
71052849c42SDave May  Insert points from db2 into db1
71152849c42SDave May  db1 <<== db2
71252849c42SDave May  */
71377048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketInsertValues(DMSwarmDataBucket db1,DMSwarmDataBucket db2)
71452849c42SDave May {
7155c18a9d6SDave May   PetscInt n_mp_points1,n_mp_points2;
7165c18a9d6SDave May   PetscInt n_mp_points1_new,p;
717dbe06d34SDave May   PetscErrorCode ierr;
71852849c42SDave May 
719521f74f9SMatthew G. Knepley   PetscFunctionBegin;
720*ea78f98cSLisandro Dalcin   ierr = DMSwarmDataBucketGetSizes(db1,&n_mp_points1,NULL,NULL);CHKERRQ(ierr);
721*ea78f98cSLisandro Dalcin   ierr = DMSwarmDataBucketGetSizes(db2,&n_mp_points2,NULL,NULL);CHKERRQ(ierr);
72252849c42SDave May   n_mp_points1_new = n_mp_points1 + n_mp_points2;
72377048351SPatrick Sanan   ierr = DMSwarmDataBucketSetSizes(db1,n_mp_points1_new,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr);
724521f74f9SMatthew G. Knepley   for (p = 0; p < n_mp_points2; ++p) {
725521f74f9SMatthew G. Knepley     /* db1 <<== db2 */
72677048351SPatrick Sanan     ierr = DMSwarmDataBucketCopyPoint(db2,p, db1,(n_mp_points1 + p));CHKERRQ(ierr);
72752849c42SDave May   }
7287fbf63aeSDave May   PetscFunctionReturn(0);
72952849c42SDave May }
73052849c42SDave May 
73152849c42SDave May /* helpers for parallel send/recv */
73277048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketCreatePackedArray(DMSwarmDataBucket db,size_t *bytes,void **buf)
73352849c42SDave May {
7345c18a9d6SDave May   PetscInt       f;
73552849c42SDave May   size_t         sizeof_marker_contents;
73652849c42SDave May   void          *buffer;
737521f74f9SMatthew G. Knepley   PetscErrorCode ierr;
73852849c42SDave May 
739521f74f9SMatthew G. Knepley   PetscFunctionBegin;
74052849c42SDave May   sizeof_marker_contents = 0;
741521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
74277048351SPatrick Sanan     DMSwarmDataField df = db->field[f];
74352849c42SDave May     sizeof_marker_contents += df->atomic_size;
74452849c42SDave May   }
745521f74f9SMatthew G. Knepley   ierr = PetscMalloc(sizeof_marker_contents, &buffer);CHKERRQ(ierr);
746521f74f9SMatthew G. Knepley   ierr = PetscMemzero(buffer, sizeof_marker_contents);CHKERRQ(ierr);
74752849c42SDave May   if (bytes) {*bytes = sizeof_marker_contents;}
74852849c42SDave May   if (buf)   {*buf   = buffer;}
7497fbf63aeSDave May   PetscFunctionReturn(0);
75052849c42SDave May }
75152849c42SDave May 
75277048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketDestroyPackedArray(DMSwarmDataBucket db,void **buf)
75352849c42SDave May {
754521f74f9SMatthew G. Knepley   PetscErrorCode ierr;
755521f74f9SMatthew G. Knepley 
756521f74f9SMatthew G. Knepley   PetscFunctionBegin;
75752849c42SDave May   if (buf) {
758521f74f9SMatthew G. Knepley     ierr = PetscFree(*buf);CHKERRQ(ierr);
75952849c42SDave May     *buf = NULL;
76052849c42SDave May   }
7617fbf63aeSDave May   PetscFunctionReturn(0);
76252849c42SDave May }
76352849c42SDave May 
76477048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketFillPackedArray(DMSwarmDataBucket db,const PetscInt index,void *buf)
76552849c42SDave May {
7665c18a9d6SDave May   PetscInt       f;
76752849c42SDave May   void          *data, *data_p;
76852849c42SDave May   size_t         asize, offset;
769521f74f9SMatthew G. Knepley   PetscErrorCode ierr;
77052849c42SDave May 
771521f74f9SMatthew G. Knepley   PetscFunctionBegin;
77252849c42SDave May   offset = 0;
773521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
77477048351SPatrick Sanan     DMSwarmDataField df = db->field[f];
77552849c42SDave May 
77652849c42SDave May     asize = df->atomic_size;
77752849c42SDave May     data = (void*)( df->data );
77852849c42SDave May     data_p = (void*)( (char*)data + index*asize );
779521f74f9SMatthew G. Knepley     ierr = PetscMemcpy((void*)((char*)buf + offset), data_p, asize);CHKERRQ(ierr);
78052849c42SDave May     offset = offset + asize;
78152849c42SDave May   }
7827fbf63aeSDave May   PetscFunctionReturn(0);
78352849c42SDave May }
78452849c42SDave May 
78577048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketInsertPackedArray(DMSwarmDataBucket db,const PetscInt idx,void *data)
78652849c42SDave May {
7875c18a9d6SDave May   PetscInt f;
78852849c42SDave May   void *data_p;
78952849c42SDave May   size_t offset;
790dbe06d34SDave May   PetscErrorCode ierr;
79152849c42SDave May 
792521f74f9SMatthew G. Knepley   PetscFunctionBegin;
79352849c42SDave May   offset = 0;
794521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
79577048351SPatrick Sanan     DMSwarmDataField df = db->field[f];
79652849c42SDave May 
79752849c42SDave May     data_p = (void*)( (char*)data + offset );
79877048351SPatrick Sanan     ierr = DMSwarmDataFieldInsertPoint(df, idx, (void*)data_p);CHKERRQ(ierr);
79952849c42SDave May     offset = offset + df->atomic_size;
80052849c42SDave May   }
8017fbf63aeSDave May   PetscFunctionReturn(0);
80252849c42SDave May }
803