xref: /petsc/src/dm/impls/swarm/data_bucket.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
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;
752849c42SDave May 
8521f74f9SMatthew G. Knepley   PetscFunctionBegin;
952849c42SDave May   *val = PETSC_FALSE;
10521f74f9SMatthew G. Knepley   for (i = 0; i < N; ++i) {
11521f74f9SMatthew G. Knepley     PetscBool flg;
125f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscStrcmp(name, gfield[i]->name, &flg));
13521f74f9SMatthew G. Knepley     if (flg) {
1452849c42SDave May       *val = PETSC_TRUE;
152eac95f8SDave May       PetscFunctionReturn(0);
1652849c42SDave May     }
1752849c42SDave May   }
182eac95f8SDave May   PetscFunctionReturn(0);
1952849c42SDave May }
2052849c42SDave May 
2177048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldStringFindInList(const char name[],const PetscInt N,const DMSwarmDataField gfield[],PetscInt *index)
2252849c42SDave May {
235c18a9d6SDave May   PetscInt       i;
2452849c42SDave May 
25521f74f9SMatthew G. Knepley   PetscFunctionBegin;
2652849c42SDave May   *index = -1;
27521f74f9SMatthew G. Knepley   for (i = 0; i < N; ++i) {
28521f74f9SMatthew G. Knepley     PetscBool flg;
295f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscStrcmp(name, gfield[i]->name, &flg));
30521f74f9SMatthew G. Knepley     if (flg) {
3152849c42SDave May       *index = i;
322eac95f8SDave May       PetscFunctionReturn(0);
3352849c42SDave May     }
3452849c42SDave May   }
352eac95f8SDave May   PetscFunctionReturn(0);
3652849c42SDave May }
3752849c42SDave May 
38ee71fbaeSPatrick Sanan PetscErrorCode DMSwarmDataFieldCreate(const char registration_function[],const char name[],const size_t size,const PetscInt L,DMSwarmDataField *DF)
3952849c42SDave May {
4077048351SPatrick Sanan   DMSwarmDataField df;
4152849c42SDave May 
42521f74f9SMatthew G. Knepley   PetscFunctionBegin;
435f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNew(&df));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrallocpy(registration_function, &df->registration_function));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrallocpy(name, &df->name));
4652849c42SDave May   df->atomic_size = size;
4752849c42SDave May   df->L  = L;
4852c3ed93SDave May   df->bs = 1;
49521f74f9SMatthew G. Knepley   /* allocate something so we don't have to reallocate */
505f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc(size * L, &df->data));
515f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMemzero(df->data, size * L));
5252849c42SDave May   *DF = df;
532eac95f8SDave May   PetscFunctionReturn(0);
5452849c42SDave May }
5552849c42SDave May 
5677048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldDestroy(DMSwarmDataField *DF)
5752849c42SDave May {
5877048351SPatrick Sanan   DMSwarmDataField df = *DF;
5952849c42SDave May 
60521f74f9SMatthew G. Knepley   PetscFunctionBegin;
615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(df->registration_function));
625f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(df->name));
635f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(df->data));
645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(df));
6552849c42SDave May   *DF  = NULL;
662eac95f8SDave May   PetscFunctionReturn(0);
6752849c42SDave May }
6852849c42SDave May 
6952849c42SDave May /* data bucket */
7077048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketCreate(DMSwarmDataBucket *DB)
7152849c42SDave May {
7277048351SPatrick Sanan   DMSwarmDataBucket db;
7352849c42SDave May 
74521f74f9SMatthew G. Knepley   PetscFunctionBegin;
755f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNew(&db));
7652849c42SDave May 
7752849c42SDave May   db->finalised = PETSC_FALSE;
7852849c42SDave May   /* create empty spaces for fields */
793454631fSDave May   db->L         = -1;
8052849c42SDave May   db->buffer    = 1;
8152849c42SDave May   db->allocated = 1;
8252849c42SDave May   db->nfields   = 0;
835f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(1, &db->field));
8452849c42SDave May   *DB  = db;
852eac95f8SDave May   PetscFunctionReturn(0);
8652849c42SDave May }
8752849c42SDave May 
8877048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketDestroy(DMSwarmDataBucket *DB)
8952849c42SDave May {
9077048351SPatrick Sanan   DMSwarmDataBucket db = *DB;
915c18a9d6SDave May   PetscInt          f;
9252849c42SDave May 
93521f74f9SMatthew G. Knepley   PetscFunctionBegin;
9452849c42SDave May   /* release fields */
95521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
965f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataFieldDestroy(&db->field[f]));
9752849c42SDave May   }
9852849c42SDave May   /* this will catch the initially allocated objects in the event that no fields are registered */
9952849c42SDave May   if (db->field != NULL) {
1005f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(db->field));
10152849c42SDave May   }
1025f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(db));
10352849c42SDave May   *DB = NULL;
1042eac95f8SDave May   PetscFunctionReturn(0);
10552849c42SDave May }
10652849c42SDave May 
10777048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketQueryForActiveFields(DMSwarmDataBucket db,PetscBool *any_active_fields)
1080cb17e37SDave May {
1090cb17e37SDave May   PetscInt f;
110521f74f9SMatthew G. Knepley 
111521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1120cb17e37SDave May   *any_active_fields = PETSC_FALSE;
113521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
1140cb17e37SDave May     if (db->field[f]->active) {
1150cb17e37SDave May       *any_active_fields = PETSC_TRUE;
1160cb17e37SDave May       PetscFunctionReturn(0);
1170cb17e37SDave May     }
1180cb17e37SDave May   }
1190cb17e37SDave May   PetscFunctionReturn(0);
1200cb17e37SDave May }
1210cb17e37SDave May 
1225627991aSBarry Smith PetscErrorCode DMSwarmDataBucketRegisterField(DMSwarmDataBucket db,const char registration_function[],const char field_name[],size_t atomic_size, DMSwarmDataField *_gfield)
12352849c42SDave May {
12452849c42SDave May   PetscBool        val;
12577048351SPatrick Sanan   DMSwarmDataField fp;
12652849c42SDave May 
127521f74f9SMatthew G. Knepley   PetscFunctionBegin;
12852849c42SDave May   /* check we haven't finalised the registration of fields */
12952849c42SDave May         /*
13052849c42SDave May    if (db->finalised==PETSC_TRUE) {
13177048351SPatrick Sanan    printf("ERROR: DMSwarmDataBucketFinalize() has been called. Cannot register more fields\n");
13252849c42SDave May    ERROR();
13352849c42SDave May    }
13452849c42SDave May   */
13552849c42SDave May   /* check for repeated name */
1365f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataFieldStringInList(field_name, db->nfields, (const DMSwarmDataField*) db->field, &val));
1372c71b3e2SJacob Faibussowitsch   PetscCheckFalse(val == PETSC_TRUE,PETSC_COMM_SELF,PETSC_ERR_USER,"Field %s already exists. Cannot add same field twice",field_name);
13852849c42SDave May   /* create new space for data */
1395f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRealloc(sizeof(DMSwarmDataField)*(db->nfields+1), &db->field));
14052849c42SDave May   /* add field */
1415f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataFieldCreate(registration_function, field_name, atomic_size, db->allocated, &fp));
14252849c42SDave May   db->field[db->nfields] = fp;
14352849c42SDave May   db->nfields++;
14452849c42SDave May   if (_gfield != NULL) {
14552849c42SDave May     *_gfield = fp;
14652849c42SDave May   }
1472eac95f8SDave May   PetscFunctionReturn(0);
14852849c42SDave May }
14952849c42SDave May 
15052849c42SDave May /*
15177048351SPatrick Sanan  #define DMSwarmDataBucketRegisterField(db,name,size,k) {\
15252849c42SDave May  char *location;\
15352849c42SDave May  asprintf(&location,"Registered by %s() at line %d within file %s", __FUNCTION__, __LINE__, __FILE__);\
15477048351SPatrick Sanan  _DMSwarmDataBucketRegisterField( (db), location, (name), (size), (k));\
155521f74f9SMatthew G. Knepley  ierr = PetscFree(location);\
15652849c42SDave May  }
15752849c42SDave May  */
15852849c42SDave May 
15977048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketGetDMSwarmDataFieldByName(DMSwarmDataBucket db,const char name[],DMSwarmDataField *gfield)
16052849c42SDave May {
1615c18a9d6SDave May   PetscInt       idx;
16252849c42SDave May   PetscBool      found;
16352849c42SDave May 
164521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1655f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataFieldStringInList(name,db->nfields,(const DMSwarmDataField*)db->field,&found));
166*28b400f6SJacob Faibussowitsch   PetscCheck(found,PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot find DMSwarmDataField with name %s",name);
1675f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataFieldStringFindInList(name,db->nfields,(const DMSwarmDataField*)db->field,&idx));
16852849c42SDave May   *gfield = db->field[idx];
1692eac95f8SDave May   PetscFunctionReturn(0);
17052849c42SDave May }
17152849c42SDave May 
17277048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketQueryDMSwarmDataFieldByName(DMSwarmDataBucket db,const char name[],PetscBool *found)
17352849c42SDave May {
174521f74f9SMatthew G. Knepley   PetscFunctionBegin;
17552849c42SDave May   *found = PETSC_FALSE;
1765f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataFieldStringInList(name,db->nfields,(const DMSwarmDataField*)db->field,found));
1777fbf63aeSDave May   PetscFunctionReturn(0);
17852849c42SDave May }
17952849c42SDave May 
18077048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketFinalize(DMSwarmDataBucket db)
18152849c42SDave May {
182521f74f9SMatthew G. Knepley   PetscFunctionBegin;
18352849c42SDave May   db->finalised = PETSC_TRUE;
1847fbf63aeSDave May   PetscFunctionReturn(0);
18552849c42SDave May }
18652849c42SDave May 
18777048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldGetNumEntries(DMSwarmDataField df,PetscInt *sum)
18852849c42SDave May {
189521f74f9SMatthew G. Knepley   PetscFunctionBegin;
19052849c42SDave May   *sum = df->L;
1917fbf63aeSDave May   PetscFunctionReturn(0);
19252849c42SDave May }
19352849c42SDave May 
19477048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldSetBlockSize(DMSwarmDataField df,PetscInt blocksize)
19552c3ed93SDave May {
196521f74f9SMatthew G. Knepley   PetscFunctionBegin;
19752c3ed93SDave May   df->bs = blocksize;
19852c3ed93SDave May   PetscFunctionReturn(0);
19952c3ed93SDave May }
20052c3ed93SDave May 
20177048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldSetSize(DMSwarmDataField df,const PetscInt new_L)
20252849c42SDave May {
203521f74f9SMatthew G. Knepley   PetscFunctionBegin;
2042c71b3e2SJacob Faibussowitsch   PetscCheckFalse(new_L < 0,PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot set size of DMSwarmDataField to be < 0");
2057fbf63aeSDave May   if (new_L == df->L) PetscFunctionReturn(0);
20652849c42SDave May   if (new_L > df->L) {
2075f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRealloc(df->atomic_size * (new_L), &df->data));
20852849c42SDave May     /* init new contents */
2095f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMemzero(( ((char*)df->data)+df->L*df->atomic_size), (new_L-df->L)*df->atomic_size));
2107fbf63aeSDave May   } else {
21152849c42SDave May     /* reallocate pointer list, add +1 in case new_L = 0 */
2125f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRealloc(df->atomic_size * (new_L+1), &df->data));
21352849c42SDave May   }
21452849c42SDave May   df->L = new_L;
2157fbf63aeSDave May   PetscFunctionReturn(0);
21652849c42SDave May }
21752849c42SDave May 
21877048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldZeroBlock(DMSwarmDataField df,const PetscInt start,const PetscInt end)
21952849c42SDave May {
220521f74f9SMatthew G. Knepley   PetscFunctionBegin;
2212c71b3e2SJacob Faibussowitsch   PetscCheckFalse(start > end,PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if start(%D) > end(%D)",start,end);
2222c71b3e2SJacob Faibussowitsch   PetscCheckFalse(start < 0,PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if start(%D) < 0",start);
2232c71b3e2SJacob Faibussowitsch   PetscCheckFalse(end > df->L,PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot zero a block of entries if end(%D) >= array size(%D)",end,df->L);
2245f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMemzero((((char*)df->data)+start*df->atomic_size), (end-start)*df->atomic_size));
2257fbf63aeSDave May   PetscFunctionReturn(0);
22652849c42SDave May }
22752849c42SDave May 
22852849c42SDave May /*
22952849c42SDave May  A negative buffer value will simply be ignored and the old buffer value will be used.
23052849c42SDave May  */
23177048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketSetSizes(DMSwarmDataBucket db,const PetscInt L,const PetscInt buffer)
23252849c42SDave May {
2335c18a9d6SDave May   PetscInt       current_allocated,new_used,new_unused,new_buffer,new_allocated,f;
2340cb17e37SDave May   PetscBool      any_active_fields;
23552849c42SDave May 
236521f74f9SMatthew G. Knepley   PetscFunctionBegin;
2372c71b3e2SJacob Faibussowitsch   PetscCheckFalse(db->finalised == PETSC_FALSE,PETSC_COMM_SELF,PETSC_ERR_USER,"You must call DMSwarmDataBucketFinalize() before DMSwarmDataBucketSetSizes()");
2385f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketQueryForActiveFields(db,&any_active_fields));
239*28b400f6SJacob Faibussowitsch   PetscCheck(!any_active_fields,PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot safely re-size as at least one DMSwarmDataField is currently being accessed");
2400cb17e37SDave May 
24152849c42SDave May   current_allocated = db->allocated;
24252849c42SDave May   new_used   = L;
24352849c42SDave May   new_unused = current_allocated - new_used;
24452849c42SDave May   new_buffer = db->buffer;
24552849c42SDave May   if (buffer >= 0) { /* update the buffer value */
24652849c42SDave May     new_buffer = buffer;
24752849c42SDave May   }
24852849c42SDave May   new_allocated = new_used + new_buffer;
24952849c42SDave May   /* action */
25052849c42SDave May   if (new_allocated > current_allocated) {
25152849c42SDave May     /* increase size to new_used + new_buffer */
25252849c42SDave May     for (f=0; f<db->nfields; f++) {
2535f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSwarmDataFieldSetSize(db->field[f], new_allocated));
25452849c42SDave May     }
25552849c42SDave May     db->L         = new_used;
25652849c42SDave May     db->buffer    = new_buffer;
25752849c42SDave May     db->allocated = new_used + new_buffer;
258521f74f9SMatthew G. Knepley   } else {
25952849c42SDave May     if (new_unused > 2 * new_buffer) {
26052849c42SDave May       /* shrink array to new_used + new_buffer */
261521f74f9SMatthew G. Knepley       for (f = 0; f < db->nfields; ++f) {
2625f80ce2aSJacob Faibussowitsch         CHKERRQ(DMSwarmDataFieldSetSize(db->field[f], new_allocated));
26352849c42SDave May       }
26452849c42SDave May       db->L         = new_used;
26552849c42SDave May       db->buffer    = new_buffer;
26652849c42SDave May       db->allocated = new_used + new_buffer;
267521f74f9SMatthew G. Knepley     } else {
26852849c42SDave May       db->L      = new_used;
26952849c42SDave May       db->buffer = new_buffer;
27052849c42SDave May     }
27152849c42SDave May   }
27252849c42SDave May   /* zero all entries from db->L to db->allocated */
273521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
27477048351SPatrick Sanan     DMSwarmDataField field = db->field[f];
2755f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataFieldZeroBlock(field, db->L,db->allocated));
27652849c42SDave May   }
2777fbf63aeSDave May   PetscFunctionReturn(0);
27852849c42SDave May }
27952849c42SDave May 
28077048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketSetInitialSizes(DMSwarmDataBucket db,const PetscInt L,const PetscInt buffer)
28152849c42SDave May {
2825c18a9d6SDave May   PetscInt       f;
283dbe06d34SDave May 
284521f74f9SMatthew G. Knepley   PetscFunctionBegin;
2855f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketSetSizes(db,L,buffer));
286521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
28777048351SPatrick Sanan     DMSwarmDataField field = db->field[f];
2885f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataFieldZeroBlock(field,0,db->allocated));
28952849c42SDave May   }
2907fbf63aeSDave May   PetscFunctionReturn(0);
29152849c42SDave May }
29252849c42SDave May 
29377048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketGetSizes(DMSwarmDataBucket db,PetscInt *L,PetscInt *buffer,PetscInt *allocated)
29452849c42SDave May {
295521f74f9SMatthew G. Knepley   PetscFunctionBegin;
29652849c42SDave May   if (L) {*L = db->L;}
29752849c42SDave May   if (buffer) {*buffer = db->buffer;}
29852849c42SDave May   if (allocated) {*allocated = db->allocated;}
2997fbf63aeSDave May   PetscFunctionReturn(0);
30052849c42SDave May }
30152849c42SDave May 
30277048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketGetGlobalSizes(MPI_Comm comm,DMSwarmDataBucket db,PetscInt *L,PetscInt *buffer,PetscInt *allocated)
30352849c42SDave May {
304521f74f9SMatthew G. Knepley   PetscFunctionBegin;
3055f80ce2aSJacob Faibussowitsch   if (L) CHKERRMPI(MPI_Allreduce(&db->L,L,1,MPIU_INT,MPI_SUM,comm));
3065f80ce2aSJacob Faibussowitsch   if (buffer) CHKERRMPI(MPI_Allreduce(&db->buffer,buffer,1,MPIU_INT,MPI_SUM,comm));
3075f80ce2aSJacob Faibussowitsch   if (allocated) CHKERRMPI(MPI_Allreduce(&db->allocated,allocated,1,MPIU_INT,MPI_SUM,comm));
3087fbf63aeSDave May   PetscFunctionReturn(0);
30952849c42SDave May }
31052849c42SDave May 
31177048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketGetDMSwarmDataFields(DMSwarmDataBucket db,PetscInt *L,DMSwarmDataField *fields[])
31252849c42SDave May {
313521f74f9SMatthew G. Knepley   PetscFunctionBegin;
31452849c42SDave May   if (L)      {*L      = db->nfields;}
31552849c42SDave May   if (fields) {*fields = db->field;}
3167fbf63aeSDave May   PetscFunctionReturn(0);
31752849c42SDave May }
31852849c42SDave May 
31977048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldGetAccess(const DMSwarmDataField gfield)
32052849c42SDave May {
321521f74f9SMatthew G. Knepley   PetscFunctionBegin;
322*28b400f6SJacob Faibussowitsch   PetscCheck(!gfield->active,PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is already active. You must call DMSwarmDataFieldRestoreAccess()",gfield->name);
32352849c42SDave May   gfield->active = PETSC_TRUE;
3247fbf63aeSDave May   PetscFunctionReturn(0);
32552849c42SDave May }
32652849c42SDave May 
32777048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldAccessPoint(const DMSwarmDataField gfield,const PetscInt pid,void **ctx_p)
32852849c42SDave May {
329521f74f9SMatthew G. Knepley   PetscFunctionBegin;
33072505a4dSJed Brown   *ctx_p = NULL;
331497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
33252849c42SDave May   /* debug mode */
33384bcda08SDave May   /* check point is valid */
3342c71b3e2SJacob Faibussowitsch   PetscCheckFalse(pid < 0,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
3352c71b3e2SJacob Faibussowitsch   PetscCheckFalse(pid >= gfield->L,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",gfield->L);
3362c71b3e2SJacob Faibussowitsch   PetscCheckFalse(gfield->active == PETSC_FALSE,PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is not active. You must call DMSwarmDataFieldGetAccess() before point data can be retrivied",gfield->name);
33752849c42SDave May #endif
33877048351SPatrick Sanan   *ctx_p = DMSWARM_DATAFIELD_point_access(gfield->data,pid,gfield->atomic_size);
3397fbf63aeSDave May   PetscFunctionReturn(0);
34052849c42SDave May }
34152849c42SDave May 
34277048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldAccessPointOffset(const DMSwarmDataField gfield,const size_t offset,const PetscInt pid,void **ctx_p)
34352849c42SDave May {
344521f74f9SMatthew G. Knepley   PetscFunctionBegin;
345497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
34652849c42SDave May   /* debug mode */
34784bcda08SDave May   /* check point is valid */
3482c71b3e2SJacob Faibussowitsch   /* PetscCheckFalse(offset < 0,PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be >= 0");*/
349521f74f9SMatthew G. Knepley   /* Note compiler realizes this can never happen with an unsigned PetscInt */
3502c71b3e2SJacob Faibussowitsch   PetscCheckFalse(offset >= gfield->atomic_size,PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be < %zu",gfield->atomic_size);
35184bcda08SDave May   /* check point is valid */
3522c71b3e2SJacob Faibussowitsch   PetscCheckFalse(pid < 0,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
3532c71b3e2SJacob Faibussowitsch   PetscCheckFalse(pid >= gfield->L,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",gfield->L);
3542c71b3e2SJacob Faibussowitsch   PetscCheckFalse(gfield->active == PETSC_FALSE,PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is not active. You must call DMSwarmDataFieldGetAccess() before point data can be retrivied",gfield->name);
35552849c42SDave May #endif
35677048351SPatrick Sanan   *ctx_p = DMSWARM_DATAFIELD_point_access_offset(gfield->data,pid,gfield->atomic_size,offset);
3577fbf63aeSDave May   PetscFunctionReturn(0);
35852849c42SDave May }
35952849c42SDave May 
36077048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldRestoreAccess(DMSwarmDataField gfield)
36152849c42SDave May {
362521f74f9SMatthew G. Knepley   PetscFunctionBegin;
3632c71b3e2SJacob Faibussowitsch   PetscCheckFalse(gfield->active == PETSC_FALSE,PETSC_COMM_SELF,PETSC_ERR_USER,"Field \"%s\" is not active. You must call DMSwarmDataFieldGetAccess()", gfield->name);
36452849c42SDave May   gfield->active = PETSC_FALSE;
3657fbf63aeSDave May   PetscFunctionReturn(0);
36652849c42SDave May }
36752849c42SDave May 
36877048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldVerifyAccess(const DMSwarmDataField gfield,const size_t size)
36952849c42SDave May {
370521f74f9SMatthew G. Knepley   PetscFunctionBegin;
371497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
3722c71b3e2SJacob Faibussowitsch   PetscCheckFalse(gfield->atomic_size != size,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);
37352849c42SDave May #endif
3747fbf63aeSDave May   PetscFunctionReturn(0);
37552849c42SDave May }
37652849c42SDave May 
37777048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldGetAtomicSize(const DMSwarmDataField gfield,size_t *size)
37852849c42SDave May {
379521f74f9SMatthew G. Knepley   PetscFunctionBegin;
38052849c42SDave May   if (size) {*size = gfield->atomic_size;}
3817fbf63aeSDave May   PetscFunctionReturn(0);
38252849c42SDave May }
38352849c42SDave May 
38477048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldGetEntries(const DMSwarmDataField gfield,void **data)
38552849c42SDave May {
386521f74f9SMatthew G. Knepley   PetscFunctionBegin;
387521f74f9SMatthew G. Knepley   if (data) {*data = gfield->data;}
3887fbf63aeSDave May   PetscFunctionReturn(0);
38952849c42SDave May }
39052849c42SDave May 
39177048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldRestoreEntries(const DMSwarmDataField gfield,void **data)
39252849c42SDave May {
393521f74f9SMatthew G. Knepley   PetscFunctionBegin;
394521f74f9SMatthew G. Knepley   if (data) {*data = NULL;}
3957fbf63aeSDave May   PetscFunctionReturn(0);
39652849c42SDave May }
39752849c42SDave May 
39852849c42SDave May /* y = x */
3995627991aSBarry Smith PetscErrorCode DMSwarmDataBucketCopyPoint(const DMSwarmDataBucket xb,const PetscInt pid_x,const DMSwarmDataBucket yb,const PetscInt pid_y)
40052849c42SDave May {
4015c18a9d6SDave May   PetscInt       f;
402dbe06d34SDave May 
403521f74f9SMatthew G. Knepley   PetscFunctionBegin;
404521f74f9SMatthew G. Knepley   for (f = 0; f < xb->nfields; ++f) {
40552849c42SDave May     void *dest;
40652849c42SDave May     void *src;
40752849c42SDave May 
4085f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataFieldGetAccess(xb->field[f]));
4095f80ce2aSJacob Faibussowitsch     if (xb != yb) CHKERRQ(DMSwarmDataFieldGetAccess( yb->field[f]));
4105f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataFieldAccessPoint(xb->field[f],pid_x, &src));
4115f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataFieldAccessPoint(yb->field[f],pid_y, &dest));
4125f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMemcpy(dest, src, xb->field[f]->atomic_size));
4135f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataFieldRestoreAccess(xb->field[f]));
4145f80ce2aSJacob Faibussowitsch     if (xb != yb) CHKERRQ(DMSwarmDataFieldRestoreAccess(yb->field[f]));
41552849c42SDave May   }
4167fbf63aeSDave May   PetscFunctionReturn(0);
41752849c42SDave May }
41852849c42SDave May 
41977048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketCreateFromSubset(DMSwarmDataBucket DBIn,const PetscInt N,const PetscInt list[],DMSwarmDataBucket *DB)
42052849c42SDave May {
4215c18a9d6SDave May   PetscInt         nfields;
42277048351SPatrick Sanan   DMSwarmDataField *fields;
4235c18a9d6SDave May   PetscInt         f,L,buffer,allocated,p;
42452849c42SDave May 
425521f74f9SMatthew G. Knepley   PetscFunctionBegin;
4265f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketCreate(DB));
42752849c42SDave May   /* copy contents of DBIn */
4285f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketGetDMSwarmDataFields(DBIn,&nfields,&fields));
4295f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketGetSizes(DBIn,&L,&buffer,&allocated));
430521f74f9SMatthew G. Knepley   for (f = 0; f < nfields; ++f) {
4315f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataBucketRegisterField(*DB,"DMSwarmDataBucketCreateFromSubset",fields[f]->name,fields[f]->atomic_size,NULL));
43252849c42SDave May   }
4335f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketFinalize(*DB));
4345f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketSetSizes(*DB,L,buffer));
43552849c42SDave May   /* now copy the desired guys from DBIn => DB */
436521f74f9SMatthew G. Knepley   for (p = 0; p < N; ++p) {
4375f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataBucketCopyPoint(DBIn,list[p], *DB,list[p]));
43852849c42SDave May   }
4397fbf63aeSDave May   PetscFunctionReturn(0);
44052849c42SDave May }
44152849c42SDave May 
442521f74f9SMatthew G. Knepley /* insert into an exisitng location */
44377048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldInsertPoint(const DMSwarmDataField field,const PetscInt index,const void *ctx)
44452849c42SDave May {
445521f74f9SMatthew G. Knepley   PetscFunctionBegin;
446497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
44784bcda08SDave May   /* check point is valid */
4482c71b3e2SJacob Faibussowitsch   PetscCheckFalse(index < 0,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
4492c71b3e2SJacob Faibussowitsch   PetscCheckFalse(index >= field->L,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",field->L);
45052849c42SDave May #endif
4515f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMemcpy(DMSWARM_DATAFIELD_point_access(field->data,index,field->atomic_size), ctx, field->atomic_size));
4527fbf63aeSDave May   PetscFunctionReturn(0);
45352849c42SDave May }
45452849c42SDave May 
455521f74f9SMatthew G. Knepley /* remove data at index - replace with last point */
45677048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketRemovePointAtIndex(const DMSwarmDataBucket db,const PetscInt index)
45752849c42SDave May {
4585c18a9d6SDave May   PetscInt       f;
4590cb17e37SDave May   PetscBool      any_active_fields;
46052849c42SDave May 
461521f74f9SMatthew G. Knepley   PetscFunctionBegin;
462497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
46384bcda08SDave May   /* check point is valid */
4642c71b3e2SJacob Faibussowitsch   PetscCheckFalse(index < 0,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
4652c71b3e2SJacob Faibussowitsch   PetscCheckFalse(index >= db->allocated,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",db->L+db->buffer);
46652849c42SDave May #endif
4675f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketQueryForActiveFields(db,&any_active_fields));
468*28b400f6SJacob Faibussowitsch   PetscCheck(!any_active_fields,PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot safely remove point as at least one DMSwarmDataField is currently being accessed");
469a233d522SDave May   if (index >= db->L) { /* this point is not in the list - no need to error, but I will anyway */
47098921bdaSJacob Faibussowitsch     SETERRQ(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);
47152849c42SDave May   }
472a233d522SDave May   if (index != db->L-1) { /* not last point in list */
473521f74f9SMatthew G. Knepley     for (f = 0; f < db->nfields; ++f) {
47477048351SPatrick Sanan       DMSwarmDataField field = db->field[f];
47552849c42SDave May 
47652849c42SDave May       /* copy then remove */
4775f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSwarmDataFieldCopyPoint(db->L-1, field, index, field));
47877048351SPatrick Sanan       /* DMSwarmDataFieldZeroPoint(field,index); */
47952849c42SDave May     }
48052849c42SDave May   }
48152849c42SDave May   /* decrement size */
48252849c42SDave May   /* this will zero out an crap at the end of the list */
4835f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketRemovePoint(db));
4847fbf63aeSDave May   PetscFunctionReturn(0);
48552849c42SDave May }
48652849c42SDave May 
48752849c42SDave May /* copy x into y */
4885627991aSBarry Smith PetscErrorCode DMSwarmDataFieldCopyPoint(const PetscInt pid_x,const DMSwarmDataField field_x,const PetscInt pid_y,const DMSwarmDataField field_y)
48952849c42SDave May {
490521f74f9SMatthew G. Knepley   PetscFunctionBegin;
491497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
49284bcda08SDave May   /* check point is valid */
4932c71b3e2SJacob Faibussowitsch   PetscCheckFalse(pid_x < 0,PETSC_COMM_SELF,PETSC_ERR_USER,"(IN) index must be >= 0");
4942c71b3e2SJacob Faibussowitsch   PetscCheckFalse(pid_x >= field_x->L,PETSC_COMM_SELF,PETSC_ERR_USER,"(IN) index must be < %D",field_x->L);
4952c71b3e2SJacob Faibussowitsch   PetscCheckFalse(pid_y < 0,PETSC_COMM_SELF,PETSC_ERR_USER,"(OUT) index must be >= 0");
4962c71b3e2SJacob Faibussowitsch   PetscCheckFalse(pid_y >= field_y->L,PETSC_COMM_SELF,PETSC_ERR_USER,"(OUT) index must be < %D",field_y->L);
4972c71b3e2SJacob Faibussowitsch   PetscCheckFalse(field_y->atomic_size != field_x->atomic_size,PETSC_COMM_SELF,PETSC_ERR_USER,"atomic size must match");
49852849c42SDave May #endif
4995f80ce2aSJacob Faibussowitsch   CHKERRQ(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));
5007fbf63aeSDave May   PetscFunctionReturn(0);
50152849c42SDave May }
50252849c42SDave May 
503521f74f9SMatthew G. Knepley /* zero only the datafield at this point */
50477048351SPatrick Sanan PetscErrorCode DMSwarmDataFieldZeroPoint(const DMSwarmDataField field,const PetscInt index)
50552849c42SDave May {
506521f74f9SMatthew G. Knepley   PetscFunctionBegin;
507497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
50884bcda08SDave May   /* check point is valid */
5092c71b3e2SJacob Faibussowitsch   PetscCheckFalse(index < 0,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
5102c71b3e2SJacob Faibussowitsch   PetscCheckFalse(index >= field->L,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",field->L);
51152849c42SDave May #endif
5125f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMemzero(DMSWARM_DATAFIELD_point_access(field->data,index,field->atomic_size), field->atomic_size));
5137fbf63aeSDave May   PetscFunctionReturn(0);
51452849c42SDave May }
51552849c42SDave May 
516521f74f9SMatthew G. Knepley /* zero ALL data for this point */
51777048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketZeroPoint(const DMSwarmDataBucket db,const PetscInt index)
51852849c42SDave May {
5195c18a9d6SDave May   PetscInt       f;
52052849c42SDave May 
521521f74f9SMatthew G. Knepley   PetscFunctionBegin;
52284bcda08SDave May   /* check point is valid */
5232c71b3e2SJacob Faibussowitsch   PetscCheckFalse(index < 0,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be >= 0");
5242c71b3e2SJacob Faibussowitsch   PetscCheckFalse(index >= db->allocated,PETSC_COMM_SELF,PETSC_ERR_USER,"index must be < %D",db->allocated);
525521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
52677048351SPatrick Sanan     DMSwarmDataField field = db->field[f];
5275f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataFieldZeroPoint(field,index));
52852849c42SDave May   }
5297fbf63aeSDave May   PetscFunctionReturn(0);
53052849c42SDave May }
53152849c42SDave May 
53252849c42SDave May /* increment */
53377048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketAddPoint(DMSwarmDataBucket db)
53452849c42SDave May {
535521f74f9SMatthew G. Knepley   PetscFunctionBegin;
5365f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketSetSizes(db,db->L+1,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
5377fbf63aeSDave May   PetscFunctionReturn(0);
53852849c42SDave May }
53952849c42SDave May 
5407fbf63aeSDave May /* decrement */
54177048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketRemovePoint(DMSwarmDataBucket db)
5427fbf63aeSDave May {
543521f74f9SMatthew G. Knepley   PetscFunctionBegin;
5445f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketSetSizes(db,db->L-1,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
5457fbf63aeSDave May   PetscFunctionReturn(0);
5467fbf63aeSDave May }
5477fbf63aeSDave May 
5485627991aSBarry Smith /*  Should be redone to user PetscViewer */
54977048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketView_stdout(MPI_Comm comm,DMSwarmDataBucket db)
550ab2be9a3SDave May {
551ab2be9a3SDave May   PetscInt       f;
552ab2be9a3SDave May   double         memory_usage_total,memory_usage_total_local = 0.0;
553ab2be9a3SDave May 
554ab2be9a3SDave May   PetscFunctionBegin;
5555f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(comm,"DMSwarmDataBucketView: \n"));
5565f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(comm,"  L                  = %D \n", db->L));
5575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(comm,"  buffer             = %D \n", db->buffer));
5585f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(comm,"  allocated          = %D \n", db->allocated));
5595f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(comm,"  nfields registered = %D \n", db->nfields));
560ab2be9a3SDave May 
561ab2be9a3SDave May   for (f = 0; f < db->nfields; ++f) {
562ab2be9a3SDave May     double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
563ab2be9a3SDave May     memory_usage_total_local += memory_usage_f;
564ab2be9a3SDave May   }
5655f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Allreduce(&memory_usage_total_local,&memory_usage_total,1,MPI_DOUBLE,MPI_SUM,comm));
566ab2be9a3SDave May 
567ab2be9a3SDave May   for (f = 0; f < db->nfields; ++f) {
568ab2be9a3SDave May     double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
5695f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(comm,"    [%3D] %15s : Mem. usage       = %1.2e (MB) [rank0]\n", f, db->field[f]->name, memory_usage_f));
5705f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(comm,"                            blocksize        = %D \n", db->field[f]->bs));
571ab2be9a3SDave May     if (db->field[f]->bs != 1) {
5725f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(comm,"                            atomic size      = %zu [full block, bs=%D]\n", db->field[f]->atomic_size,db->field[f]->bs));
5735f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(comm,"                            atomic size/item = %zu \n", db->field[f]->atomic_size/db->field[f]->bs));
574ab2be9a3SDave May     } else {
5755f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(comm,"                            atomic size      = %zu \n", db->field[f]->atomic_size));
576ab2be9a3SDave May     }
577ab2be9a3SDave May   }
5785f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(comm,"  Total mem. usage                           = %1.2e (MB) (collective)\n", memory_usage_total));
579ab2be9a3SDave May   PetscFunctionReturn(0);
580ab2be9a3SDave May }
581ab2be9a3SDave May 
5825627991aSBarry Smith PetscErrorCode DMSwarmDataBucketView_Seq(MPI_Comm comm,DMSwarmDataBucket db,const char filename[],DMSwarmDataBucketViewType type)
58352849c42SDave May {
584521f74f9SMatthew G. Knepley   PetscFunctionBegin;
58552849c42SDave May   switch (type) {
58652849c42SDave May   case DATABUCKET_VIEW_STDOUT:
5875f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataBucketView_stdout(PETSC_COMM_SELF,db));
58852849c42SDave May     break;
58952849c42SDave May   case DATABUCKET_VIEW_ASCII:
590071ce369SDave May     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for ascii output");
59152849c42SDave May   case DATABUCKET_VIEW_BINARY:
592071ce369SDave May     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for binary output");
59352849c42SDave May   case DATABUCKET_VIEW_HDF5:
594071ce369SDave May     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for HDF5 output");
595521f74f9SMatthew G. Knepley   default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer method requested");
59652849c42SDave May   }
5977fbf63aeSDave May   PetscFunctionReturn(0);
59852849c42SDave May }
59952849c42SDave May 
60077048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketView_MPI(MPI_Comm comm,DMSwarmDataBucket db,const char filename[],DMSwarmDataBucketViewType type)
60152849c42SDave May {
602521f74f9SMatthew G. Knepley   PetscFunctionBegin;
60352849c42SDave May   switch (type) {
60452849c42SDave May   case DATABUCKET_VIEW_STDOUT:
6055f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataBucketView_stdout(comm,db));
60652849c42SDave May     break;
60752849c42SDave May   case DATABUCKET_VIEW_ASCII:
608071ce369SDave May     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for ascii output");
60952849c42SDave May   case DATABUCKET_VIEW_BINARY:
610071ce369SDave May     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for binary output");
61152849c42SDave May   case DATABUCKET_VIEW_HDF5:
612071ce369SDave May     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for HDF5 output");
613521f74f9SMatthew G. Knepley   default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer method requested");
61452849c42SDave May   }
6157fbf63aeSDave May   PetscFunctionReturn(0);
61652849c42SDave May }
61752849c42SDave May 
61877048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketView(MPI_Comm comm,DMSwarmDataBucket db,const char filename[],DMSwarmDataBucketViewType type)
61952849c42SDave May {
620d7d19db6SBarry Smith   PetscMPIInt    size;
62152849c42SDave May 
622521f74f9SMatthew G. Knepley   PetscFunctionBegin;
6235f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(comm,&size));
624d7d19db6SBarry Smith   if (size == 1) {
6255f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataBucketView_Seq(comm,db,filename,type));
62652849c42SDave May   } else {
6275f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataBucketView_MPI(comm,db,filename,type));
62852849c42SDave May   }
6297fbf63aeSDave May   PetscFunctionReturn(0);
63052849c42SDave May }
63152849c42SDave May 
63277048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketDuplicateFields(DMSwarmDataBucket dbA,DMSwarmDataBucket *dbB)
63352849c42SDave May {
63477048351SPatrick Sanan   DMSwarmDataBucket db2;
6355c18a9d6SDave May   PetscInt          f;
63652849c42SDave May 
637521f74f9SMatthew G. Knepley   PetscFunctionBegin;
6385f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketCreate(&db2));
63952849c42SDave May   /* copy contents from dbA into db2 */
640521f74f9SMatthew G. Knepley   for (f = 0; f < dbA->nfields; ++f) {
64177048351SPatrick Sanan     DMSwarmDataField field;
64252849c42SDave May     size_t           atomic_size;
64352849c42SDave May     char             *name;
64452849c42SDave May 
64552849c42SDave May     field = dbA->field[f];
64652849c42SDave May     atomic_size = field->atomic_size;
64752849c42SDave May     name        = field->name;
6485f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataBucketRegisterField(db2,"DMSwarmDataBucketDuplicateFields",name,atomic_size,NULL));
64952849c42SDave May   }
6505f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketFinalize(db2));
6515f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketSetInitialSizes(db2,0,1000));
65252849c42SDave May   *dbB = db2;
6537fbf63aeSDave May   PetscFunctionReturn(0);
65452849c42SDave May }
65552849c42SDave May 
65652849c42SDave May /*
65752849c42SDave May  Insert points from db2 into db1
65852849c42SDave May  db1 <<== db2
65952849c42SDave May  */
66077048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketInsertValues(DMSwarmDataBucket db1,DMSwarmDataBucket db2)
66152849c42SDave May {
6625c18a9d6SDave May   PetscInt       n_mp_points1,n_mp_points2;
6635c18a9d6SDave May   PetscInt       n_mp_points1_new,p;
66452849c42SDave May 
665521f74f9SMatthew G. Knepley   PetscFunctionBegin;
6665f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketGetSizes(db1,&n_mp_points1,NULL,NULL));
6675f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketGetSizes(db2,&n_mp_points2,NULL,NULL));
66852849c42SDave May   n_mp_points1_new = n_mp_points1 + n_mp_points2;
6695f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDataBucketSetSizes(db1,n_mp_points1_new,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
670521f74f9SMatthew G. Knepley   for (p = 0; p < n_mp_points2; ++p) {
671521f74f9SMatthew G. Knepley     /* db1 <<== db2 */
6725f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataBucketCopyPoint(db2,p, db1,(n_mp_points1 + p)));
67352849c42SDave May   }
6747fbf63aeSDave May   PetscFunctionReturn(0);
67552849c42SDave May }
67652849c42SDave May 
67752849c42SDave May /* helpers for parallel send/recv */
67877048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketCreatePackedArray(DMSwarmDataBucket db,size_t *bytes,void **buf)
67952849c42SDave May {
6805c18a9d6SDave May   PetscInt       f;
68152849c42SDave May   size_t         sizeof_marker_contents;
68252849c42SDave May   void          *buffer;
68352849c42SDave May 
684521f74f9SMatthew G. Knepley   PetscFunctionBegin;
68552849c42SDave May   sizeof_marker_contents = 0;
686521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
68777048351SPatrick Sanan     DMSwarmDataField df = db->field[f];
68852849c42SDave May     sizeof_marker_contents += df->atomic_size;
68952849c42SDave May   }
6905f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc(sizeof_marker_contents, &buffer));
6915f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMemzero(buffer, sizeof_marker_contents));
69252849c42SDave May   if (bytes) {*bytes = sizeof_marker_contents;}
69352849c42SDave May   if (buf)   {*buf   = buffer;}
6947fbf63aeSDave May   PetscFunctionReturn(0);
69552849c42SDave May }
69652849c42SDave May 
69777048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketDestroyPackedArray(DMSwarmDataBucket db,void **buf)
69852849c42SDave May {
699521f74f9SMatthew G. Knepley   PetscFunctionBegin;
70052849c42SDave May   if (buf) {
7015f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(*buf));
70252849c42SDave May     *buf = NULL;
70352849c42SDave May   }
7047fbf63aeSDave May   PetscFunctionReturn(0);
70552849c42SDave May }
70652849c42SDave May 
70777048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketFillPackedArray(DMSwarmDataBucket db,const PetscInt index,void *buf)
70852849c42SDave May {
7095c18a9d6SDave May   PetscInt       f;
71052849c42SDave May   void          *data, *data_p;
71152849c42SDave May   size_t         asize, offset;
71252849c42SDave May 
713521f74f9SMatthew G. Knepley   PetscFunctionBegin;
71452849c42SDave May   offset = 0;
715521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
71677048351SPatrick Sanan     DMSwarmDataField df = db->field[f];
71752849c42SDave May 
71852849c42SDave May     asize = df->atomic_size;
71952849c42SDave May     data = (void*)( df->data);
72052849c42SDave May     data_p = (void*)( (char*)data + index*asize);
7215f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMemcpy((void*)((char*)buf + offset), data_p, asize));
72252849c42SDave May     offset = offset + asize;
72352849c42SDave May   }
7247fbf63aeSDave May   PetscFunctionReturn(0);
72552849c42SDave May }
72652849c42SDave May 
72777048351SPatrick Sanan PetscErrorCode DMSwarmDataBucketInsertPackedArray(DMSwarmDataBucket db,const PetscInt idx,void *data)
72852849c42SDave May {
7295c18a9d6SDave May   PetscInt       f;
73052849c42SDave May   void           *data_p;
73152849c42SDave May   size_t         offset;
73252849c42SDave May 
733521f74f9SMatthew G. Knepley   PetscFunctionBegin;
73452849c42SDave May   offset = 0;
735521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
73677048351SPatrick Sanan     DMSwarmDataField df = db->field[f];
73752849c42SDave May 
73852849c42SDave May     data_p = (void*)( (char*)data + offset);
7395f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmDataFieldInsertPoint(df, idx, (void*)data_p));
74052849c42SDave May     offset = offset + df->atomic_size;
74152849c42SDave May   }
7427fbf63aeSDave May   PetscFunctionReturn(0);
74352849c42SDave May }
744