xref: /petsc/src/dm/impls/swarm/data_bucket.c (revision 670292f542dd698a111da7ff6b3eb1295e84523d)
1279f676cSBarry Smith #include "../src/dm/impls/swarm/data_bucket.h"
252849c42SDave May 
352849c42SDave May /* string helpers */
4d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldStringInList(const char name[], const PetscInt N, const DMSwarmDataField gfield[], PetscBool *val)
5d71ae5a4SJacob Faibussowitsch {
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;
129566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(name, gfield[i]->name, &flg));
13521f74f9SMatthew G. Knepley     if (flg) {
1452849c42SDave May       *val = PETSC_TRUE;
153ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
1652849c42SDave May     }
1752849c42SDave May   }
183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1952849c42SDave May }
2052849c42SDave May 
21d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldStringFindInList(const char name[], const PetscInt N, const DMSwarmDataField gfield[], PetscInt *index)
22d71ae5a4SJacob Faibussowitsch {
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;
299566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(name, gfield[i]->name, &flg));
30521f74f9SMatthew G. Knepley     if (flg) {
3152849c42SDave May       *index = i;
323ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
3352849c42SDave May     }
3452849c42SDave May   }
353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3652849c42SDave May }
3752849c42SDave May 
38d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldCreate(const char registration_function[], const char name[], const size_t size, const PetscInt L, DMSwarmDataField *DF)
39d71ae5a4SJacob Faibussowitsch {
4077048351SPatrick Sanan   DMSwarmDataField df;
4152849c42SDave May 
42521f74f9SMatthew G. Knepley   PetscFunctionBegin;
439566063dSJacob Faibussowitsch   PetscCall(PetscNew(&df));
449566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(registration_function, &df->registration_function));
459566063dSJacob Faibussowitsch   PetscCall(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 */
509566063dSJacob Faibussowitsch   PetscCall(PetscMalloc(size * L, &df->data));
519566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(df->data, size * L));
5252849c42SDave May   *DF = df;
533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5452849c42SDave May }
5552849c42SDave May 
56d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldDestroy(DMSwarmDataField *DF)
57d71ae5a4SJacob Faibussowitsch {
5877048351SPatrick Sanan   DMSwarmDataField df = *DF;
5952849c42SDave May 
60521f74f9SMatthew G. Knepley   PetscFunctionBegin;
619566063dSJacob Faibussowitsch   PetscCall(PetscFree(df->registration_function));
629566063dSJacob Faibussowitsch   PetscCall(PetscFree(df->name));
639566063dSJacob Faibussowitsch   PetscCall(PetscFree(df->data));
649566063dSJacob Faibussowitsch   PetscCall(PetscFree(df));
6552849c42SDave May   *DF = NULL;
663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6752849c42SDave May }
6852849c42SDave May 
6952849c42SDave May /* data bucket */
70d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketCreate(DMSwarmDataBucket *DB)
71d71ae5a4SJacob Faibussowitsch {
7277048351SPatrick Sanan   DMSwarmDataBucket db;
7352849c42SDave May 
74521f74f9SMatthew G. Knepley   PetscFunctionBegin;
759566063dSJacob Faibussowitsch   PetscCall(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;
839566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(1, &db->field));
8452849c42SDave May   *DB = db;
853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8652849c42SDave May }
8752849c42SDave May 
88d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketDestroy(DMSwarmDataBucket *DB)
89d71ae5a4SJacob Faibussowitsch {
9077048351SPatrick Sanan   DMSwarmDataBucket db = *DB;
915c18a9d6SDave May   PetscInt          f;
9252849c42SDave May 
93521f74f9SMatthew G. Knepley   PetscFunctionBegin;
9452849c42SDave May   /* release fields */
9548a46eb9SPierre Jolivet   for (f = 0; f < db->nfields; ++f) PetscCall(DMSwarmDataFieldDestroy(&db->field[f]));
9652849c42SDave May   /* this will catch the initially allocated objects in the event that no fields are registered */
9748a46eb9SPierre Jolivet   if (db->field != NULL) PetscCall(PetscFree(db->field));
989566063dSJacob Faibussowitsch   PetscCall(PetscFree(db));
9952849c42SDave May   *DB = NULL;
1003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10152849c42SDave May }
10252849c42SDave May 
103d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketQueryForActiveFields(DMSwarmDataBucket db, PetscBool *any_active_fields)
104d71ae5a4SJacob Faibussowitsch {
1050cb17e37SDave May   PetscInt f;
106521f74f9SMatthew G. Knepley 
107521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1080cb17e37SDave May   *any_active_fields = PETSC_FALSE;
109521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
1100cb17e37SDave May     if (db->field[f]->active) {
1110cb17e37SDave May       *any_active_fields = PETSC_TRUE;
1123ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
1130cb17e37SDave May     }
1140cb17e37SDave May   }
1153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1160cb17e37SDave May }
1170cb17e37SDave May 
118d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketRegisterField(DMSwarmDataBucket db, const char registration_function[], const char field_name[], size_t atomic_size, DMSwarmDataField *_gfield)
119d71ae5a4SJacob Faibussowitsch {
12052849c42SDave May   PetscBool        val;
12177048351SPatrick Sanan   DMSwarmDataField fp;
12252849c42SDave May 
123521f74f9SMatthew G. Knepley   PetscFunctionBegin;
12452849c42SDave May   /* check we haven't finalised the registration of fields */
12552849c42SDave May   /*
12652849c42SDave May    if (db->finalised==PETSC_TRUE) {
12777048351SPatrick Sanan    printf("ERROR: DMSwarmDataBucketFinalize() has been called. Cannot register more fields\n");
12852849c42SDave May    ERROR();
12952849c42SDave May    }
13052849c42SDave May   */
13152849c42SDave May   /* check for repeated name */
1329566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataFieldStringInList(field_name, db->nfields, (const DMSwarmDataField *)db->field, &val));
13308401ef6SPierre Jolivet   PetscCheck(val != PETSC_TRUE, PETSC_COMM_SELF, PETSC_ERR_USER, "Field %s already exists. Cannot add same field twice", field_name);
13452849c42SDave May   /* create new space for data */
1359566063dSJacob Faibussowitsch   PetscCall(PetscRealloc(sizeof(DMSwarmDataField) * (db->nfields + 1), &db->field));
13652849c42SDave May   /* add field */
1379566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataFieldCreate(registration_function, field_name, atomic_size, db->allocated, &fp));
13852849c42SDave May   db->field[db->nfields] = fp;
13952849c42SDave May   db->nfields++;
140ad540459SPierre Jolivet   if (_gfield != NULL) *_gfield = fp;
1413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14252849c42SDave May }
14352849c42SDave May 
14452849c42SDave May /*
14577048351SPatrick Sanan  #define DMSwarmDataBucketRegisterField(db,name,size,k) {\
14652849c42SDave May  char *location;\
14752849c42SDave May  asprintf(&location,"Registered by %s() at line %d within file %s", __FUNCTION__, __LINE__, __FILE__);\
14877048351SPatrick Sanan  _DMSwarmDataBucketRegisterField( (db), location, (name), (size), (k));\
149521f74f9SMatthew G. Knepley  ierr = PetscFree(location);\
15052849c42SDave May  }
15152849c42SDave May  */
15252849c42SDave May 
153d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketGetDMSwarmDataFieldIdByName(DMSwarmDataBucket db, const char name[], PetscInt *idx)
154d71ae5a4SJacob Faibussowitsch {
1552e956fe4SStefano Zampini   PetscBool found;
1562e956fe4SStefano Zampini 
1572e956fe4SStefano Zampini   PetscFunctionBegin;
1582e956fe4SStefano Zampini   *idx = -1;
1592e956fe4SStefano Zampini   PetscCall(DMSwarmDataFieldStringInList(name, db->nfields, (const DMSwarmDataField *)db->field, &found));
1602e956fe4SStefano Zampini   PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot find DMSwarmDataField with name %s", name);
1612e956fe4SStefano Zampini   PetscCall(DMSwarmDataFieldStringFindInList(name, db->nfields, (const DMSwarmDataField *)db->field, idx));
1623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1632e956fe4SStefano Zampini }
1642e956fe4SStefano Zampini 
165d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketGetDMSwarmDataFieldByName(DMSwarmDataBucket db, const char name[], DMSwarmDataField *gfield)
166d71ae5a4SJacob Faibussowitsch {
1675c18a9d6SDave May   PetscInt  idx;
16852849c42SDave May   PetscBool found;
16952849c42SDave May 
170521f74f9SMatthew G. Knepley   PetscFunctionBegin;
1719566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataFieldStringInList(name, db->nfields, (const DMSwarmDataField *)db->field, &found));
17228b400f6SJacob Faibussowitsch   PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot find DMSwarmDataField with name %s", name);
1739566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataFieldStringFindInList(name, db->nfields, (const DMSwarmDataField *)db->field, &idx));
17452849c42SDave May   *gfield = db->field[idx];
1753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17652849c42SDave May }
17752849c42SDave May 
178d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketQueryDMSwarmDataFieldByName(DMSwarmDataBucket db, const char name[], PetscBool *found)
179d71ae5a4SJacob Faibussowitsch {
180521f74f9SMatthew G. Knepley   PetscFunctionBegin;
18152849c42SDave May   *found = PETSC_FALSE;
1829566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataFieldStringInList(name, db->nfields, (const DMSwarmDataField *)db->field, found));
1833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18452849c42SDave May }
18552849c42SDave May 
186d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketFinalize(DMSwarmDataBucket db)
187d71ae5a4SJacob Faibussowitsch {
188521f74f9SMatthew G. Knepley   PetscFunctionBegin;
18952849c42SDave May   db->finalised = PETSC_TRUE;
1903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19152849c42SDave May }
19252849c42SDave May 
193d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldGetNumEntries(DMSwarmDataField df, PetscInt *sum)
194d71ae5a4SJacob Faibussowitsch {
195521f74f9SMatthew G. Knepley   PetscFunctionBegin;
19652849c42SDave May   *sum = df->L;
1973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19852849c42SDave May }
19952849c42SDave May 
200d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldSetBlockSize(DMSwarmDataField df, PetscInt blocksize)
201d71ae5a4SJacob Faibussowitsch {
202521f74f9SMatthew G. Knepley   PetscFunctionBegin;
20352c3ed93SDave May   df->bs = blocksize;
2043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
20552c3ed93SDave May }
20652c3ed93SDave May 
207d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldSetSize(DMSwarmDataField df, const PetscInt new_L)
208d71ae5a4SJacob Faibussowitsch {
209521f74f9SMatthew G. Knepley   PetscFunctionBegin;
21008401ef6SPierre Jolivet   PetscCheck(new_L >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot set size of DMSwarmDataField to be < 0");
2113ba16761SJacob Faibussowitsch   if (new_L == df->L) PetscFunctionReturn(PETSC_SUCCESS);
21252849c42SDave May   if (new_L > df->L) {
2139566063dSJacob Faibussowitsch     PetscCall(PetscRealloc(df->atomic_size * (new_L), &df->data));
21452849c42SDave May     /* init new contents */
21557508eceSPierre Jolivet     PetscCall(PetscMemzero(((char *)df->data) + df->L * df->atomic_size, (new_L - df->L) * df->atomic_size));
2167fbf63aeSDave May   } else {
21752849c42SDave May     /* reallocate pointer list, add +1 in case new_L = 0 */
2189566063dSJacob Faibussowitsch     PetscCall(PetscRealloc(df->atomic_size * (new_L + 1), &df->data));
21952849c42SDave May   }
22052849c42SDave May   df->L = new_L;
2213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
22252849c42SDave May }
22352849c42SDave May 
224d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldZeroBlock(DMSwarmDataField df, const PetscInt start, const PetscInt end)
225d71ae5a4SJacob Faibussowitsch {
226521f74f9SMatthew G. Knepley   PetscFunctionBegin;
22763a3b9bcSJacob Faibussowitsch   PetscCheck(start <= end, PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot zero a block of entries if start(%" PetscInt_FMT ") > end(%" PetscInt_FMT ")", start, end);
22863a3b9bcSJacob Faibussowitsch   PetscCheck(start >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot zero a block of entries if start(%" PetscInt_FMT ") < 0", start);
22963a3b9bcSJacob Faibussowitsch   PetscCheck(end <= df->L, PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot zero a block of entries if end(%" PetscInt_FMT ") >= array size(%" PetscInt_FMT ")", end, df->L);
23057508eceSPierre Jolivet   PetscCall(PetscMemzero(((char *)df->data) + start * df->atomic_size, (end - start) * df->atomic_size));
2313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23252849c42SDave May }
23352849c42SDave May 
23452849c42SDave May /*
23552849c42SDave May  A negative buffer value will simply be ignored and the old buffer value will be used.
23652849c42SDave May  */
237d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketSetSizes(DMSwarmDataBucket db, const PetscInt L, const PetscInt buffer)
238d71ae5a4SJacob Faibussowitsch {
2395c18a9d6SDave May   PetscInt  current_allocated, new_used, new_unused, new_buffer, new_allocated, f;
2400cb17e37SDave May   PetscBool any_active_fields;
24152849c42SDave May 
242521f74f9SMatthew G. Knepley   PetscFunctionBegin;
24308401ef6SPierre Jolivet   PetscCheck(db->finalised != PETSC_FALSE, PETSC_COMM_SELF, PETSC_ERR_USER, "You must call DMSwarmDataBucketFinalize() before DMSwarmDataBucketSetSizes()");
2449566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketQueryForActiveFields(db, &any_active_fields));
24528b400f6SJacob Faibussowitsch   PetscCheck(!any_active_fields, PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot safely re-size as at least one DMSwarmDataField is currently being accessed");
2460cb17e37SDave May 
24752849c42SDave May   current_allocated = db->allocated;
24852849c42SDave May   new_used          = L;
24952849c42SDave May   new_unused        = current_allocated - new_used;
25052849c42SDave May   new_buffer        = db->buffer;
25152849c42SDave May   if (buffer >= 0) { /* update the buffer value */
25252849c42SDave May     new_buffer = buffer;
25352849c42SDave May   }
25452849c42SDave May   new_allocated = new_used + new_buffer;
25552849c42SDave May   /* action */
25652849c42SDave May   if (new_allocated > current_allocated) {
25752849c42SDave May     /* increase size to new_used + new_buffer */
25848a46eb9SPierre Jolivet     for (f = 0; f < db->nfields; f++) PetscCall(DMSwarmDataFieldSetSize(db->field[f], new_allocated));
25952849c42SDave May     db->L         = new_used;
26052849c42SDave May     db->buffer    = new_buffer;
26152849c42SDave May     db->allocated = new_used + new_buffer;
262521f74f9SMatthew G. Knepley   } else {
26352849c42SDave May     if (new_unused > 2 * new_buffer) {
26452849c42SDave May       /* shrink array to new_used + new_buffer */
26548a46eb9SPierre Jolivet       for (f = 0; f < db->nfields; ++f) PetscCall(DMSwarmDataFieldSetSize(db->field[f], new_allocated));
26652849c42SDave May       db->L         = new_used;
26752849c42SDave May       db->buffer    = new_buffer;
26852849c42SDave May       db->allocated = new_used + new_buffer;
269521f74f9SMatthew G. Knepley     } else {
27052849c42SDave May       db->L      = new_used;
27152849c42SDave May       db->buffer = new_buffer;
27252849c42SDave May     }
27352849c42SDave May   }
27452849c42SDave May   /* zero all entries from db->L to db->allocated */
275521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
27677048351SPatrick Sanan     DMSwarmDataField field = db->field[f];
2779566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataFieldZeroBlock(field, db->L, db->allocated));
27852849c42SDave May   }
2793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
28052849c42SDave May }
28152849c42SDave May 
282d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketSetInitialSizes(DMSwarmDataBucket db, const PetscInt L, const PetscInt buffer)
283d71ae5a4SJacob Faibussowitsch {
2845c18a9d6SDave May   PetscInt f;
285dbe06d34SDave May 
286521f74f9SMatthew G. Knepley   PetscFunctionBegin;
2879566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketSetSizes(db, L, buffer));
288521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
28977048351SPatrick Sanan     DMSwarmDataField field = db->field[f];
2909566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataFieldZeroBlock(field, 0, db->allocated));
29152849c42SDave May   }
2923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
29352849c42SDave May }
29452849c42SDave May 
295d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketGetSizes(DMSwarmDataBucket db, PetscInt *L, PetscInt *buffer, PetscInt *allocated)
296d71ae5a4SJacob Faibussowitsch {
297521f74f9SMatthew G. Knepley   PetscFunctionBegin;
298ad540459SPierre Jolivet   if (L) *L = db->L;
299ad540459SPierre Jolivet   if (buffer) *buffer = db->buffer;
300ad540459SPierre Jolivet   if (allocated) *allocated = db->allocated;
3013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
30252849c42SDave May }
30352849c42SDave May 
304d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketGetGlobalSizes(MPI_Comm comm, DMSwarmDataBucket db, PetscInt *L, PetscInt *buffer, PetscInt *allocated)
305d71ae5a4SJacob Faibussowitsch {
306521f74f9SMatthew G. Knepley   PetscFunctionBegin;
307462c564dSBarry Smith   if (L) PetscCallMPI(MPIU_Allreduce(&db->L, L, 1, MPIU_INT, MPI_SUM, comm));
308462c564dSBarry Smith   if (buffer) PetscCallMPI(MPIU_Allreduce(&db->buffer, buffer, 1, MPIU_INT, MPI_SUM, comm));
309462c564dSBarry Smith   if (allocated) PetscCallMPI(MPIU_Allreduce(&db->allocated, allocated, 1, MPIU_INT, MPI_SUM, comm));
3103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
31152849c42SDave May }
31252849c42SDave May 
313d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketGetDMSwarmDataFields(DMSwarmDataBucket db, PetscInt *L, DMSwarmDataField *fields[])
314d71ae5a4SJacob Faibussowitsch {
315521f74f9SMatthew G. Knepley   PetscFunctionBegin;
316ad540459SPierre Jolivet   if (L) *L = db->nfields;
317ad540459SPierre Jolivet   if (fields) *fields = db->field;
3183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
31952849c42SDave May }
32052849c42SDave May 
321d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldGetAccess(const DMSwarmDataField gfield)
322d71ae5a4SJacob Faibussowitsch {
323521f74f9SMatthew G. Knepley   PetscFunctionBegin;
32428b400f6SJacob Faibussowitsch   PetscCheck(!gfield->active, PETSC_COMM_SELF, PETSC_ERR_USER, "Field \"%s\" is already active. You must call DMSwarmDataFieldRestoreAccess()", gfield->name);
32552849c42SDave May   gfield->active = PETSC_TRUE;
3263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
32752849c42SDave May }
32852849c42SDave May 
329d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldAccessPoint(const DMSwarmDataField gfield, const PetscInt pid, void **ctx_p)
330d71ae5a4SJacob Faibussowitsch {
331521f74f9SMatthew G. Knepley   PetscFunctionBegin;
33272505a4dSJed Brown   *ctx_p = NULL;
333497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
33452849c42SDave May   /* debug mode */
33584bcda08SDave May   /* check point is valid */
33608401ef6SPierre Jolivet   PetscCheck(pid >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be >= 0");
33763a3b9bcSJacob Faibussowitsch   PetscCheck(pid < gfield->L, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be < %" PetscInt_FMT, gfield->L);
33808401ef6SPierre Jolivet   PetscCheck(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);
33952849c42SDave May #endif
34077048351SPatrick Sanan   *ctx_p = DMSWARM_DATAFIELD_point_access(gfield->data, pid, gfield->atomic_size);
3413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
34252849c42SDave May }
34352849c42SDave May 
344d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldAccessPointOffset(const DMSwarmDataField gfield, const size_t offset, const PetscInt pid, void **ctx_p)
345d71ae5a4SJacob Faibussowitsch {
346521f74f9SMatthew G. Knepley   PetscFunctionBegin;
347497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
34852849c42SDave May   /* debug mode */
34984bcda08SDave May   /* check point is valid */
35008401ef6SPierre Jolivet   /* PetscCheck(offset >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"offset must be >= 0");*/
351521f74f9SMatthew G. Knepley   /* Note compiler realizes this can never happen with an unsigned PetscInt */
35208401ef6SPierre Jolivet   PetscCheck(offset < gfield->atomic_size, PETSC_COMM_SELF, PETSC_ERR_USER, "offset must be < %zu", gfield->atomic_size);
35384bcda08SDave May   /* check point is valid */
35408401ef6SPierre Jolivet   PetscCheck(pid >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be >= 0");
35563a3b9bcSJacob Faibussowitsch   PetscCheck(pid < gfield->L, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be < %" PetscInt_FMT, gfield->L);
35608401ef6SPierre Jolivet   PetscCheck(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);
35752849c42SDave May #endif
35877048351SPatrick Sanan   *ctx_p = DMSWARM_DATAFIELD_point_access_offset(gfield->data, pid, gfield->atomic_size, offset);
3593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
36052849c42SDave May }
36152849c42SDave May 
362d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldRestoreAccess(DMSwarmDataField gfield)
363d71ae5a4SJacob Faibussowitsch {
364521f74f9SMatthew G. Knepley   PetscFunctionBegin;
36508401ef6SPierre Jolivet   PetscCheck(gfield->active != PETSC_FALSE, PETSC_COMM_SELF, PETSC_ERR_USER, "Field \"%s\" is not active. You must call DMSwarmDataFieldGetAccess()", gfield->name);
36652849c42SDave May   gfield->active = PETSC_FALSE;
3673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
36852849c42SDave May }
36952849c42SDave May 
370d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldVerifyAccess(const DMSwarmDataField gfield, const size_t size)
371d71ae5a4SJacob Faibussowitsch {
372521f74f9SMatthew G. Knepley   PetscFunctionBegin;
373497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
37408401ef6SPierre Jolivet   PetscCheck(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);
37552849c42SDave May #endif
3763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
37752849c42SDave May }
37852849c42SDave May 
379d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldGetAtomicSize(const DMSwarmDataField gfield, size_t *size)
380d71ae5a4SJacob Faibussowitsch {
381521f74f9SMatthew G. Knepley   PetscFunctionBegin;
382ad540459SPierre Jolivet   if (size) *size = gfield->atomic_size;
3833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
38452849c42SDave May }
38552849c42SDave May 
386d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldGetEntries(const DMSwarmDataField gfield, void **data)
387d71ae5a4SJacob Faibussowitsch {
388521f74f9SMatthew G. Knepley   PetscFunctionBegin;
389ad540459SPierre Jolivet   if (data) *data = gfield->data;
3903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
39152849c42SDave May }
39252849c42SDave May 
393d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldRestoreEntries(const DMSwarmDataField gfield, void **data)
394d71ae5a4SJacob Faibussowitsch {
395521f74f9SMatthew G. Knepley   PetscFunctionBegin;
396ad540459SPierre Jolivet   if (data) *data = NULL;
3973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
39852849c42SDave May }
39952849c42SDave May 
40052849c42SDave May /* y = x */
401d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketCopyPoint(const DMSwarmDataBucket xb, const PetscInt pid_x, const DMSwarmDataBucket yb, const PetscInt pid_y)
402d71ae5a4SJacob Faibussowitsch {
4035c18a9d6SDave May   PetscInt f;
404dbe06d34SDave May 
405521f74f9SMatthew G. Knepley   PetscFunctionBegin;
406521f74f9SMatthew G. Knepley   for (f = 0; f < xb->nfields; ++f) {
40752849c42SDave May     void *dest;
40852849c42SDave May     void *src;
40952849c42SDave May 
4109566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataFieldGetAccess(xb->field[f]));
4119566063dSJacob Faibussowitsch     if (xb != yb) PetscCall(DMSwarmDataFieldGetAccess(yb->field[f]));
4129566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataFieldAccessPoint(xb->field[f], pid_x, &src));
4139566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataFieldAccessPoint(yb->field[f], pid_y, &dest));
4149566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(dest, src, xb->field[f]->atomic_size));
4159566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataFieldRestoreAccess(xb->field[f]));
4169566063dSJacob Faibussowitsch     if (xb != yb) PetscCall(DMSwarmDataFieldRestoreAccess(yb->field[f]));
41752849c42SDave May   }
4183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
41952849c42SDave May }
42052849c42SDave May 
421d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketCreateFromSubset(DMSwarmDataBucket DBIn, const PetscInt N, const PetscInt list[], DMSwarmDataBucket *DB)
422d71ae5a4SJacob Faibussowitsch {
4235c18a9d6SDave May   PetscInt          nfields;
42477048351SPatrick Sanan   DMSwarmDataField *fields;
4255c18a9d6SDave May   PetscInt          f, L, buffer, allocated, p;
42652849c42SDave May 
427521f74f9SMatthew G. Knepley   PetscFunctionBegin;
4289566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCreate(DB));
42952849c42SDave May   /* copy contents of DBIn */
4309566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetDMSwarmDataFields(DBIn, &nfields, &fields));
4319566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(DBIn, &L, &buffer, &allocated));
43248a46eb9SPierre Jolivet   for (f = 0; f < nfields; ++f) PetscCall(DMSwarmDataBucketRegisterField(*DB, "DMSwarmDataBucketCreateFromSubset", fields[f]->name, fields[f]->atomic_size, NULL));
4339566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketFinalize(*DB));
4349566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketSetSizes(*DB, L, buffer));
43519307e5cSMatthew G. Knepley   for (f = 0; f < nfields; ++f) {
43619307e5cSMatthew G. Knepley     DMSwarmDataField gfield;
43719307e5cSMatthew G. Knepley 
43819307e5cSMatthew G. Knepley     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(*DB, fields[f]->name, &gfield));
43919307e5cSMatthew G. Knepley     PetscCall(DMSwarmDataFieldSetBlockSize(gfield, fields[f]->bs));
44019307e5cSMatthew G. Knepley     gfield->petsc_type = fields[f]->petsc_type;
44119307e5cSMatthew G. Knepley   }
44252849c42SDave May   /* now copy the desired guys from DBIn => DB */
44348a46eb9SPierre Jolivet   for (p = 0; p < N; ++p) PetscCall(DMSwarmDataBucketCopyPoint(DBIn, list[p], *DB, list[p]));
4443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
44552849c42SDave May }
44652849c42SDave May 
447d5b43468SJose E. Roman /* insert into an existing location */
448d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldInsertPoint(const DMSwarmDataField field, const PetscInt index, const void *ctx)
449d71ae5a4SJacob Faibussowitsch {
450521f74f9SMatthew G. Knepley   PetscFunctionBegin;
451497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
45284bcda08SDave May   /* check point is valid */
45308401ef6SPierre Jolivet   PetscCheck(index >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be >= 0");
45463a3b9bcSJacob Faibussowitsch   PetscCheck(index < field->L, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be < %" PetscInt_FMT, field->L);
45552849c42SDave May #endif
4569566063dSJacob Faibussowitsch   PetscCall(PetscMemcpy(DMSWARM_DATAFIELD_point_access(field->data, index, field->atomic_size), ctx, field->atomic_size));
4573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
45852849c42SDave May }
45952849c42SDave May 
460521f74f9SMatthew G. Knepley /* remove data at index - replace with last point */
461d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketRemovePointAtIndex(const DMSwarmDataBucket db, const PetscInt index)
462d71ae5a4SJacob Faibussowitsch {
4635c18a9d6SDave May   PetscInt  f;
4640cb17e37SDave May   PetscBool any_active_fields;
46552849c42SDave May 
466521f74f9SMatthew G. Knepley   PetscFunctionBegin;
467497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
46884bcda08SDave May   /* check point is valid */
46908401ef6SPierre Jolivet   PetscCheck(index >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be >= 0");
47063a3b9bcSJacob Faibussowitsch   PetscCheck(index < db->allocated, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be < %" PetscInt_FMT, db->L + db->buffer);
47152849c42SDave May #endif
4729566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketQueryForActiveFields(db, &any_active_fields));
47328b400f6SJacob Faibussowitsch   PetscCheck(!any_active_fields, PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot safely remove point as at least one DMSwarmDataField is currently being accessed");
474a233d522SDave May   if (index >= db->L) { /* this point is not in the list - no need to error, but I will anyway */
47563a3b9bcSJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "You should not be trying to remove point at index=%" PetscInt_FMT " since it's < db->L = %" PetscInt_FMT, index, db->L);
47652849c42SDave May   }
477a233d522SDave May   if (index != db->L - 1) { /* not last point in list */
478521f74f9SMatthew G. Knepley     for (f = 0; f < db->nfields; ++f) {
47977048351SPatrick Sanan       DMSwarmDataField field = db->field[f];
48052849c42SDave May 
48152849c42SDave May       /* copy then remove */
4829566063dSJacob Faibussowitsch       PetscCall(DMSwarmDataFieldCopyPoint(db->L - 1, field, index, field));
48377048351SPatrick Sanan       /* DMSwarmDataFieldZeroPoint(field,index); */
48452849c42SDave May     }
48552849c42SDave May   }
48652849c42SDave May   /* decrement size */
48752849c42SDave May   /* this will zero out an crap at the end of the list */
4889566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketRemovePoint(db));
4893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
49052849c42SDave May }
49152849c42SDave May 
49252849c42SDave May /* copy x into y */
493d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldCopyPoint(const PetscInt pid_x, const DMSwarmDataField field_x, const PetscInt pid_y, const DMSwarmDataField field_y)
494d71ae5a4SJacob Faibussowitsch {
495521f74f9SMatthew G. Knepley   PetscFunctionBegin;
496497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
49784bcda08SDave May   /* check point is valid */
49808401ef6SPierre Jolivet   PetscCheck(pid_x >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "(IN) index must be >= 0");
49963a3b9bcSJacob Faibussowitsch   PetscCheck(pid_x < field_x->L, PETSC_COMM_SELF, PETSC_ERR_USER, "(IN) index must be < %" PetscInt_FMT, field_x->L);
50008401ef6SPierre Jolivet   PetscCheck(pid_y >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "(OUT) index must be >= 0");
50163a3b9bcSJacob Faibussowitsch   PetscCheck(pid_y < field_y->L, PETSC_COMM_SELF, PETSC_ERR_USER, "(OUT) index must be < %" PetscInt_FMT, field_y->L);
50208401ef6SPierre Jolivet   PetscCheck(field_y->atomic_size == field_x->atomic_size, PETSC_COMM_SELF, PETSC_ERR_USER, "atomic size must match");
50352849c42SDave May #endif
5049566063dSJacob Faibussowitsch   PetscCall(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));
5053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
50652849c42SDave May }
50752849c42SDave May 
508521f74f9SMatthew G. Knepley /* zero only the datafield at this point */
509d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataFieldZeroPoint(const DMSwarmDataField field, const PetscInt index)
510d71ae5a4SJacob Faibussowitsch {
511521f74f9SMatthew G. Knepley   PetscFunctionBegin;
512497880caSRichard Tran Mills #if defined(DMSWARM_DATAFIELD_POINT_ACCESS_GUARD)
51384bcda08SDave May   /* check point is valid */
51408401ef6SPierre Jolivet   PetscCheck(index >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be >= 0");
51563a3b9bcSJacob Faibussowitsch   PetscCheck(index < field->L, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be < %" PetscInt_FMT, field->L);
51652849c42SDave May #endif
5179566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(DMSWARM_DATAFIELD_point_access(field->data, index, field->atomic_size), field->atomic_size));
5183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
51952849c42SDave May }
52052849c42SDave May 
521521f74f9SMatthew G. Knepley /* zero ALL data for this point */
522d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketZeroPoint(const DMSwarmDataBucket db, const PetscInt index)
523d71ae5a4SJacob Faibussowitsch {
5245c18a9d6SDave May   PetscInt f;
52552849c42SDave May 
526521f74f9SMatthew G. Knepley   PetscFunctionBegin;
52784bcda08SDave May   /* check point is valid */
52808401ef6SPierre Jolivet   PetscCheck(index >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be >= 0");
52963a3b9bcSJacob Faibussowitsch   PetscCheck(index < db->allocated, PETSC_COMM_SELF, PETSC_ERR_USER, "index must be < %" PetscInt_FMT, db->allocated);
530521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
53177048351SPatrick Sanan     DMSwarmDataField field = db->field[f];
5329566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataFieldZeroPoint(field, index));
53352849c42SDave May   }
5343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
53552849c42SDave May }
53652849c42SDave May 
53752849c42SDave May /* increment */
538d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketAddPoint(DMSwarmDataBucket db)
539d71ae5a4SJacob Faibussowitsch {
540521f74f9SMatthew G. Knepley   PetscFunctionBegin;
541*670292f5SMatthew Knepley   PetscCall(DMSwarmDataBucketSetSizes(db, PetscMax(db->L, 0) + 1, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
5423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
54352849c42SDave May }
54452849c42SDave May 
5457fbf63aeSDave May /* decrement */
546d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketRemovePoint(DMSwarmDataBucket db)
547d71ae5a4SJacob Faibussowitsch {
548521f74f9SMatthew G. Knepley   PetscFunctionBegin;
549*670292f5SMatthew Knepley   PetscCheck(db->L > 0, PetscObjectComm((PetscObject)db), PETSC_ERR_ARG_WRONG, "Swarm has no points to be removed");
5509566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketSetSizes(db, db->L - 1, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
5513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5527fbf63aeSDave May }
5537fbf63aeSDave May 
5545627991aSBarry Smith /*  Should be redone to user PetscViewer */
55566976f2fSJacob Faibussowitsch static PetscErrorCode DMSwarmDataBucketView_stdout(MPI_Comm comm, DMSwarmDataBucket db)
556d71ae5a4SJacob Faibussowitsch {
557ab2be9a3SDave May   PetscInt f;
558ab2be9a3SDave May   double   memory_usage_total, memory_usage_total_local = 0.0;
559ab2be9a3SDave May 
560ab2be9a3SDave May   PetscFunctionBegin;
5619566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "DMSwarmDataBucketView: \n"));
56263a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "  L                  = %" PetscInt_FMT " \n", db->L));
56363a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "  buffer             = %" PetscInt_FMT " \n", db->buffer));
56463a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "  allocated          = %" PetscInt_FMT " \n", db->allocated));
56563a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "  nfields registered = %" PetscInt_FMT " \n", db->nfields));
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;
569ab2be9a3SDave May     memory_usage_total_local += memory_usage_f;
570ab2be9a3SDave May   }
5716a210b70SBarry Smith   PetscCallMPI(MPIU_Allreduce(&memory_usage_total_local, &memory_usage_total, 1, MPI_DOUBLE, MPI_SUM, comm));
572ab2be9a3SDave May 
573ab2be9a3SDave May   for (f = 0; f < db->nfields; ++f) {
574ab2be9a3SDave May     double memory_usage_f = (double)(db->field[f]->atomic_size * db->allocated) * 1.0e-6;
57563a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "    [%3" PetscInt_FMT "] %15s : Mem. usage       = %1.2e (MB) [rank0]\n", f, db->field[f]->name, memory_usage_f));
57663a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "                            blocksize        = %" PetscInt_FMT " \n", db->field[f]->bs));
577ab2be9a3SDave May     if (db->field[f]->bs != 1) {
57863a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "                            atomic size      = %zu [full block, bs=%" PetscInt_FMT "]\n", db->field[f]->atomic_size, db->field[f]->bs));
57963a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "                            atomic size/item = %zu \n", (size_t)(db->field[f]->atomic_size / db->field[f]->bs)));
580ab2be9a3SDave May     } else {
5819566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "                            atomic size      = %zu \n", db->field[f]->atomic_size));
582ab2be9a3SDave May     }
583ab2be9a3SDave May   }
5849566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "  Total mem. usage                           = %1.2e (MB) (collective)\n", memory_usage_total));
5853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
586ab2be9a3SDave May }
587ab2be9a3SDave May 
58866976f2fSJacob Faibussowitsch static PetscErrorCode DMSwarmDataBucketView_Seq(MPI_Comm comm, DMSwarmDataBucket db, const char filename[], DMSwarmDataBucketViewType type)
589d71ae5a4SJacob Faibussowitsch {
590521f74f9SMatthew G. Knepley   PetscFunctionBegin;
59152849c42SDave May   switch (type) {
592d71ae5a4SJacob Faibussowitsch   case DATABUCKET_VIEW_STDOUT:
593d71ae5a4SJacob Faibussowitsch     PetscCall(DMSwarmDataBucketView_stdout(PETSC_COMM_SELF, db));
594d71ae5a4SJacob Faibussowitsch     break;
595d71ae5a4SJacob Faibussowitsch   case DATABUCKET_VIEW_ASCII:
596d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for ascii output");
597d71ae5a4SJacob Faibussowitsch   case DATABUCKET_VIEW_BINARY:
598d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for binary output");
599d71ae5a4SJacob Faibussowitsch   case DATABUCKET_VIEW_HDF5:
600d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for HDF5 output");
601d71ae5a4SJacob Faibussowitsch   default:
602d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown viewer method requested");
60352849c42SDave May   }
6043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
60552849c42SDave May }
60652849c42SDave May 
60766976f2fSJacob Faibussowitsch static PetscErrorCode DMSwarmDataBucketView_MPI(MPI_Comm comm, DMSwarmDataBucket db, const char filename[], DMSwarmDataBucketViewType type)
608d71ae5a4SJacob Faibussowitsch {
609521f74f9SMatthew G. Knepley   PetscFunctionBegin;
61052849c42SDave May   switch (type) {
611d71ae5a4SJacob Faibussowitsch   case DATABUCKET_VIEW_STDOUT:
612d71ae5a4SJacob Faibussowitsch     PetscCall(DMSwarmDataBucketView_stdout(comm, db));
613d71ae5a4SJacob Faibussowitsch     break;
614d71ae5a4SJacob Faibussowitsch   case DATABUCKET_VIEW_ASCII:
615d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for ascii output");
616d71ae5a4SJacob Faibussowitsch   case DATABUCKET_VIEW_BINARY:
617d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for binary output");
618d71ae5a4SJacob Faibussowitsch   case DATABUCKET_VIEW_HDF5:
619d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for HDF5 output");
620d71ae5a4SJacob Faibussowitsch   default:
621d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown viewer method requested");
62252849c42SDave May   }
6233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
62452849c42SDave May }
62552849c42SDave May 
626d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketView(MPI_Comm comm, DMSwarmDataBucket db, const char filename[], DMSwarmDataBucketViewType type)
627d71ae5a4SJacob Faibussowitsch {
628d7d19db6SBarry Smith   PetscMPIInt size;
62952849c42SDave May 
630521f74f9SMatthew G. Knepley   PetscFunctionBegin;
6319566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
632d7d19db6SBarry Smith   if (size == 1) {
6339566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketView_Seq(comm, db, filename, type));
63452849c42SDave May   } else {
6359566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketView_MPI(comm, db, filename, type));
63652849c42SDave May   }
6373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
63852849c42SDave May }
63952849c42SDave May 
640d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketDuplicateFields(DMSwarmDataBucket dbA, DMSwarmDataBucket *dbB)
641d71ae5a4SJacob Faibussowitsch {
64277048351SPatrick Sanan   DMSwarmDataBucket db2;
6435c18a9d6SDave May   PetscInt          f;
64452849c42SDave May 
645521f74f9SMatthew G. Knepley   PetscFunctionBegin;
6469566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketCreate(&db2));
64752849c42SDave May   /* copy contents from dbA into db2 */
648521f74f9SMatthew G. Knepley   for (f = 0; f < dbA->nfields; ++f) {
64977048351SPatrick Sanan     DMSwarmDataField field;
65052849c42SDave May     size_t           atomic_size;
65152849c42SDave May     char            *name;
65252849c42SDave May 
65352849c42SDave May     field       = dbA->field[f];
65452849c42SDave May     atomic_size = field->atomic_size;
65552849c42SDave May     name        = field->name;
6569566063dSJacob Faibussowitsch     PetscCall(DMSwarmDataBucketRegisterField(db2, "DMSwarmDataBucketDuplicateFields", name, atomic_size, NULL));
65752849c42SDave May   }
6589566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketFinalize(db2));
6599566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketSetInitialSizes(db2, 0, 1000));
66052849c42SDave May   *dbB = db2;
6613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
66252849c42SDave May }
66352849c42SDave May 
66452849c42SDave May /*
66552849c42SDave May  Insert points from db2 into db1
66652849c42SDave May  db1 <<== db2
66752849c42SDave May  */
668d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketInsertValues(DMSwarmDataBucket db1, DMSwarmDataBucket db2)
669d71ae5a4SJacob Faibussowitsch {
6705c18a9d6SDave May   PetscInt n_mp_points1, n_mp_points2;
6715c18a9d6SDave May   PetscInt n_mp_points1_new, p;
67252849c42SDave May 
673521f74f9SMatthew G. Knepley   PetscFunctionBegin;
6749566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(db1, &n_mp_points1, NULL, NULL));
6759566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketGetSizes(db2, &n_mp_points2, NULL, NULL));
67652849c42SDave May   n_mp_points1_new = n_mp_points1 + n_mp_points2;
6779566063dSJacob Faibussowitsch   PetscCall(DMSwarmDataBucketSetSizes(db1, n_mp_points1_new, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
678521f74f9SMatthew G. Knepley   for (p = 0; p < n_mp_points2; ++p) {
679521f74f9SMatthew G. Knepley     /* db1 <<== db2 */
68057508eceSPierre Jolivet     PetscCall(DMSwarmDataBucketCopyPoint(db2, p, db1, n_mp_points1 + p));
68152849c42SDave May   }
6823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
68352849c42SDave May }
68452849c42SDave May 
68552849c42SDave May /* helpers for parallel send/recv */
686d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketCreatePackedArray(DMSwarmDataBucket db, size_t *bytes, void **buf)
687d71ae5a4SJacob Faibussowitsch {
6885c18a9d6SDave May   PetscInt f;
68952849c42SDave May   size_t   sizeof_marker_contents;
69052849c42SDave May   void    *buffer;
69152849c42SDave May 
692521f74f9SMatthew G. Knepley   PetscFunctionBegin;
69352849c42SDave May   sizeof_marker_contents = 0;
694521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
69577048351SPatrick Sanan     DMSwarmDataField df = db->field[f];
69652849c42SDave May     sizeof_marker_contents += df->atomic_size;
69752849c42SDave May   }
6989566063dSJacob Faibussowitsch   PetscCall(PetscMalloc(sizeof_marker_contents, &buffer));
6999566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(buffer, sizeof_marker_contents));
700ad540459SPierre Jolivet   if (bytes) *bytes = sizeof_marker_contents;
701ad540459SPierre Jolivet   if (buf) *buf = buffer;
7023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
70352849c42SDave May }
70452849c42SDave May 
705d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketDestroyPackedArray(DMSwarmDataBucket db, void **buf)
706d71ae5a4SJacob Faibussowitsch {
707521f74f9SMatthew G. Knepley   PetscFunctionBegin;
70852849c42SDave May   if (buf) {
7099566063dSJacob Faibussowitsch     PetscCall(PetscFree(*buf));
71052849c42SDave May     *buf = NULL;
71152849c42SDave May   }
7123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
71352849c42SDave May }
71452849c42SDave May 
715d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketFillPackedArray(DMSwarmDataBucket db, const PetscInt index, void *buf)
716d71ae5a4SJacob Faibussowitsch {
7175c18a9d6SDave May   PetscInt f;
71852849c42SDave May   void    *data, *data_p;
71952849c42SDave May   size_t   asize, offset;
72052849c42SDave May 
721521f74f9SMatthew G. Knepley   PetscFunctionBegin;
72252849c42SDave May   offset = 0;
723521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
72477048351SPatrick Sanan     DMSwarmDataField df = db->field[f];
72552849c42SDave May 
72652849c42SDave May     asize  = df->atomic_size;
727835f2295SStefano Zampini     data   = df->data;
72852849c42SDave May     data_p = (void *)((char *)data + index * asize);
7299566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy((void *)((char *)buf + offset), data_p, asize));
73052849c42SDave May     offset = offset + asize;
73152849c42SDave May   }
7323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
73352849c42SDave May }
73452849c42SDave May 
735d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmDataBucketInsertPackedArray(DMSwarmDataBucket db, const PetscInt idx, void *data)
736d71ae5a4SJacob Faibussowitsch {
7375c18a9d6SDave May   PetscInt f;
73852849c42SDave May   void    *data_p;
73952849c42SDave May   size_t   offset;
74052849c42SDave May 
741521f74f9SMatthew G. Knepley   PetscFunctionBegin;
74252849c42SDave May   offset = 0;
743521f74f9SMatthew G. Knepley   for (f = 0; f < db->nfields; ++f) {
74477048351SPatrick Sanan     DMSwarmDataField df = db->field[f];
74552849c42SDave May 
74652849c42SDave May     data_p = (void *)((char *)data + offset);
747835f2295SStefano Zampini     PetscCall(DMSwarmDataFieldInsertPoint(df, idx, data_p));
74852849c42SDave May     offset = offset + df->atomic_size;
74952849c42SDave May   }
7503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
75152849c42SDave May }
752