Actual source code: swarm.c

  1: #define PETSCDM_DLL
  2: #include <petsc/private/dmswarmimpl.h>
  3: #include <petsc/private/hashsetij.h>
  4: #include <petsc/private/petscfeimpl.h>
  5: #include <petscviewer.h>
  6: #include <petscdraw.h>
  7: #include <petscdmplex.h>
  8: #include <petscblaslapack.h>
  9: #include "../src/dm/impls/swarm/data_bucket.h"
 10: #include <petscdmlabel.h>
 11: #include <petscsection.h>

 13: PetscLogEvent DMSWARM_Migrate, DMSWARM_SetSizes, DMSWARM_AddPoints, DMSWARM_RemovePoints, DMSWARM_Sort;
 14: PetscLogEvent DMSWARM_DataExchangerTopologySetup, DMSWARM_DataExchangerBegin, DMSWARM_DataExchangerEnd;
 15: PetscLogEvent DMSWARM_DataExchangerSendCount, DMSWARM_DataExchangerPack;

 17: const char *DMSwarmTypeNames[]          = {"basic", "pic", NULL};
 18: const char *DMSwarmMigrateTypeNames[]   = {"basic", "dmcellnscatter", "dmcellexact", "user", NULL};
 19: const char *DMSwarmCollectTypeNames[]   = {"basic", "boundingbox", "general", "user", NULL};
 20: const char *DMSwarmPICLayoutTypeNames[] = {"regular", "gauss", "subdivision", NULL};

 22: const char DMSwarmField_pid[]       = "DMSwarm_pid";
 23: const char DMSwarmField_rank[]      = "DMSwarm_rank";
 24: const char DMSwarmPICField_coor[]   = "DMSwarmPIC_coor";
 25: const char DMSwarmPICField_cellid[] = "DMSwarm_cellid";

 27: PetscInt SwarmDataFieldId = -1;

 29: #if defined(PETSC_HAVE_HDF5)
 30: #include <petscviewerhdf5.h>

 32: static PetscErrorCode VecView_Swarm_HDF5_Internal(Vec v, PetscViewer viewer)
 33: {
 34:   DM        dm;
 35:   PetscReal seqval;
 36:   PetscInt  seqnum, bs;
 37:   PetscBool isseq, ists;

 39:   PetscFunctionBegin;
 40:   PetscCall(VecGetDM(v, &dm));
 41:   PetscCall(VecGetBlockSize(v, &bs));
 42:   PetscCall(PetscViewerHDF5PushGroup(viewer, "/particle_fields"));
 43:   PetscCall(PetscObjectTypeCompare((PetscObject)v, VECSEQ, &isseq));
 44:   PetscCall(DMGetOutputSequenceNumber(dm, &seqnum, &seqval));
 45:   PetscCall(PetscViewerHDF5IsTimestepping(viewer, &ists));
 46:   if (ists) PetscCall(PetscViewerHDF5SetTimestep(viewer, seqnum));
 47:   /* PetscCall(DMSequenceView_HDF5(dm, "time", seqnum, (PetscScalar) seqval, viewer)); */
 48:   PetscCall(VecViewNative(v, viewer));
 49:   PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject)v, "Nc", PETSC_INT, (void *)&bs));
 50:   PetscCall(PetscViewerHDF5PopGroup(viewer));
 51:   PetscFunctionReturn(PETSC_SUCCESS);
 52: }

 54: static PetscErrorCode DMSwarmView_HDF5(DM dm, PetscViewer viewer)
 55: {
 56:   Vec       coordinates;
 57:   PetscInt  Np;
 58:   PetscBool isseq;

 60:   PetscFunctionBegin;
 61:   PetscCall(DMSwarmGetSize(dm, &Np));
 62:   PetscCall(DMSwarmCreateGlobalVectorFromField(dm, DMSwarmPICField_coor, &coordinates));
 63:   PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
 64:   PetscCall(PetscViewerHDF5PushGroup(viewer, "/particles"));
 65:   PetscCall(PetscObjectTypeCompare((PetscObject)coordinates, VECSEQ, &isseq));
 66:   PetscCall(VecViewNative(coordinates, viewer));
 67:   PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject)coordinates, "Np", PETSC_INT, (void *)&Np));
 68:   PetscCall(PetscViewerHDF5PopGroup(viewer));
 69:   PetscCall(DMSwarmDestroyGlobalVectorFromField(dm, DMSwarmPICField_coor, &coordinates));
 70:   PetscFunctionReturn(PETSC_SUCCESS);
 71: }
 72: #endif

 74: static PetscErrorCode VecView_Swarm(Vec v, PetscViewer viewer)
 75: {
 76:   DM dm;
 77: #if defined(PETSC_HAVE_HDF5)
 78:   PetscBool ishdf5;
 79: #endif

 81:   PetscFunctionBegin;
 82:   PetscCall(VecGetDM(v, &dm));
 83:   PetscCheck(dm, PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
 84: #if defined(PETSC_HAVE_HDF5)
 85:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
 86:   if (ishdf5) {
 87:     PetscCall(VecView_Swarm_HDF5_Internal(v, viewer));
 88:     PetscFunctionReturn(PETSC_SUCCESS);
 89:   }
 90: #endif
 91:   PetscCall(VecViewNative(v, viewer));
 92:   PetscFunctionReturn(PETSC_SUCCESS);
 93: }

 95: /*@
 96:   DMSwarmVectorGetField - Gets the field from which to define a `Vec` object
 97:   when `DMCreateLocalVector()`, or `DMCreateGlobalVector()` is called

 99:   Not collective

101:   Input Parameter:
102: . dm - a `DMSWARM`

104:   Output Parameter:
105: . fieldname - the textual name given to a registered field, or the empty string if it has not been set

107:   Level: beginner

109: .seealso: `DM`, `DMSWARM`, `DMSwarmVectorDefineField()` `DMSwarmRegisterPetscDatatypeField()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`
110: @*/
111: PetscErrorCode DMSwarmVectorGetField(DM dm, const char *fieldname[])
112: {
113:   PetscFunctionBegin;
115:   PetscAssertPointer(fieldname, 2);
116:   *fieldname = ((DM_Swarm *)dm->data)->vec_field_name;
117:   PetscFunctionReturn(PETSC_SUCCESS);
118: }

120: /*@
121:   DMSwarmVectorDefineField - Sets the field from which to define a `Vec` object
122:   when `DMCreateLocalVector()`, or `DMCreateGlobalVector()` is called

124:   Collective

126:   Input Parameters:
127: + dm        - a `DMSWARM`
128: - fieldname - the textual name given to a registered field

130:   Level: beginner

132:   Notes:
133:   The field with name `fieldname` must be defined as having a data type of `PetscScalar`.

135:   This function must be called prior to calling `DMCreateLocalVector()`, `DMCreateGlobalVector()`.
136:   Multiple calls to `DMSwarmVectorDefineField()` are permitted.

138: .seealso: `DM`, `DMSWARM`, `DMSwarmVectorGetField()`, `DMSwarmRegisterPetscDatatypeField()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`
139: @*/
140: PetscErrorCode DMSwarmVectorDefineField(DM dm, const char fieldname[])
141: {
142:   DM_Swarm     *swarm = (DM_Swarm *)dm->data;
143:   PetscInt      bs, n;
144:   PetscScalar  *array;
145:   PetscDataType type;

147:   PetscFunctionBegin;
149:   if (fieldname) PetscAssertPointer(fieldname, 2);
150:   if (!swarm->issetup) PetscCall(DMSetUp(dm));
151:   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &n, NULL, NULL));
152:   PetscCall(DMSwarmGetField(dm, fieldname, &bs, &type, (void **)&array));

154:   /* Check all fields are of type PETSC_REAL or PETSC_SCALAR */
155:   PetscCheck(type == PETSC_REAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only valid for PETSC_REAL");
156:   PetscCall(PetscSNPrintf(swarm->vec_field_name, PETSC_MAX_PATH_LEN - 1, "%s", fieldname));
157:   swarm->vec_field_set    = PETSC_TRUE;
158:   swarm->vec_field_bs     = bs;
159:   swarm->vec_field_nlocal = n;
160:   PetscCall(DMSwarmRestoreField(dm, fieldname, &bs, &type, (void **)&array));
161:   PetscFunctionReturn(PETSC_SUCCESS);
162: }

164: /* requires DMSwarmDefineFieldVector has been called */
165: static PetscErrorCode DMCreateGlobalVector_Swarm(DM dm, Vec *vec)
166: {
167:   DM_Swarm *swarm = (DM_Swarm *)dm->data;
168:   Vec       x;
169:   char      name[PETSC_MAX_PATH_LEN];

171:   PetscFunctionBegin;
172:   if (!swarm->issetup) PetscCall(DMSetUp(dm));
173:   PetscCheck(swarm->vec_field_set, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMSwarmVectorDefineField first");
174:   PetscCheck(swarm->vec_field_nlocal == swarm->db->L, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSwarm sizes have changed since last call to VectorDefineField first"); /* Stale data */

176:   PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "DMSwarmField_%s", swarm->vec_field_name));
177:   PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), &x));
178:   PetscCall(PetscObjectSetName((PetscObject)x, name));
179:   PetscCall(VecSetSizes(x, swarm->db->L * swarm->vec_field_bs, PETSC_DETERMINE));
180:   PetscCall(VecSetBlockSize(x, swarm->vec_field_bs));
181:   PetscCall(VecSetDM(x, dm));
182:   PetscCall(VecSetFromOptions(x));
183:   PetscCall(VecSetDM(x, dm));
184:   PetscCall(VecSetOperation(x, VECOP_VIEW, (void (*)(void))VecView_Swarm));
185:   *vec = x;
186:   PetscFunctionReturn(PETSC_SUCCESS);
187: }

189: /* requires DMSwarmDefineFieldVector has been called */
190: static PetscErrorCode DMCreateLocalVector_Swarm(DM dm, Vec *vec)
191: {
192:   DM_Swarm *swarm = (DM_Swarm *)dm->data;
193:   Vec       x;
194:   char      name[PETSC_MAX_PATH_LEN];

196:   PetscFunctionBegin;
197:   if (!swarm->issetup) PetscCall(DMSetUp(dm));
198:   PetscCheck(swarm->vec_field_set, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMSwarmVectorDefineField first");
199:   PetscCheck(swarm->vec_field_nlocal == swarm->db->L, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSwarm sizes have changed since last call to VectorDefineField first");

201:   PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "DMSwarmField_%s", swarm->vec_field_name));
202:   PetscCall(VecCreate(PETSC_COMM_SELF, &x));
203:   PetscCall(PetscObjectSetName((PetscObject)x, name));
204:   PetscCall(VecSetSizes(x, swarm->db->L * swarm->vec_field_bs, PETSC_DETERMINE));
205:   PetscCall(VecSetBlockSize(x, swarm->vec_field_bs));
206:   PetscCall(VecSetDM(x, dm));
207:   PetscCall(VecSetFromOptions(x));
208:   *vec = x;
209:   PetscFunctionReturn(PETSC_SUCCESS);
210: }

212: static PetscErrorCode DMSwarmDestroyVectorFromField_Private(DM dm, const char fieldname[], Vec *vec)
213: {
214:   DM_Swarm        *swarm = (DM_Swarm *)dm->data;
215:   DMSwarmDataField gfield;
216:   PetscInt         bs, nlocal, fid = -1, cfid = -2;
217:   PetscBool        flg;

219:   PetscFunctionBegin;
220:   /* check vector is an inplace array */
221:   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldIdByName(swarm->db, fieldname, &fid));
222:   PetscCall(PetscObjectComposedDataGetInt((PetscObject)*vec, SwarmDataFieldId, cfid, flg));
223:   PetscCheck(cfid == fid, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Vector being destroyed was not created from DMSwarm field(%s)! %" PetscInt_FMT " != %" PetscInt_FMT, fieldname, cfid, fid);
224:   PetscCall(VecGetLocalSize(*vec, &nlocal));
225:   PetscCall(VecGetBlockSize(*vec, &bs));
226:   PetscCheck(nlocal / bs == swarm->db->L, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSwarm sizes have changed since vector was created - cannot ensure pointers are valid");
227:   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
228:   PetscCall(DMSwarmDataFieldRestoreAccess(gfield));
229:   PetscCall(VecResetArray(*vec));
230:   PetscCall(VecDestroy(vec));
231:   PetscFunctionReturn(PETSC_SUCCESS);
232: }

234: static PetscErrorCode DMSwarmCreateVectorFromField_Private(DM dm, const char fieldname[], MPI_Comm comm, Vec *vec)
235: {
236:   DM_Swarm     *swarm = (DM_Swarm *)dm->data;
237:   PetscDataType type;
238:   PetscScalar  *array;
239:   PetscInt      bs, n, fid;
240:   char          name[PETSC_MAX_PATH_LEN];
241:   PetscMPIInt   size;
242:   PetscBool     iscuda, iskokkos;

244:   PetscFunctionBegin;
245:   if (!swarm->issetup) PetscCall(DMSetUp(dm));
246:   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &n, NULL, NULL));
247:   PetscCall(DMSwarmGetField(dm, fieldname, &bs, &type, (void **)&array));
248:   PetscCheck(type == PETSC_REAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only valid for PETSC_REAL");

250:   PetscCallMPI(MPI_Comm_size(comm, &size));
251:   PetscCall(PetscStrcmp(dm->vectype, VECKOKKOS, &iskokkos));
252:   PetscCall(PetscStrcmp(dm->vectype, VECCUDA, &iscuda));
253:   PetscCall(VecCreate(comm, vec));
254:   PetscCall(VecSetSizes(*vec, n * bs, PETSC_DETERMINE));
255:   PetscCall(VecSetBlockSize(*vec, bs));
256:   if (iskokkos) PetscCall(VecSetType(*vec, VECKOKKOS));
257:   else if (iscuda) PetscCall(VecSetType(*vec, VECCUDA));
258:   else PetscCall(VecSetType(*vec, VECSTANDARD));
259:   PetscCall(VecPlaceArray(*vec, array));

261:   PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "DMSwarmSharedField_%s", fieldname));
262:   PetscCall(PetscObjectSetName((PetscObject)*vec, name));

264:   /* Set guard */
265:   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldIdByName(swarm->db, fieldname, &fid));
266:   PetscCall(PetscObjectComposedDataSetInt((PetscObject)*vec, SwarmDataFieldId, fid));

268:   PetscCall(VecSetDM(*vec, dm));
269:   PetscCall(VecSetOperation(*vec, VECOP_VIEW, (void (*)(void))VecView_Swarm));
270:   PetscFunctionReturn(PETSC_SUCCESS);
271: }

273: /* This creates a "mass matrix" between a finite element and particle space. If a finite element interpolant is given by

275:      \hat f = \sum_i f_i \phi_i

277:    and a particle function is given by

279:      f = \sum_i w_i \delta(x - x_i)

281:    then we want to require that

283:      M \hat f = M_p f

285:    where the particle mass matrix is given by

287:      (M_p)_{ij} = \int \phi_i \delta(x - x_j)

289:    The way Dave May does particles, they amount to quadratue weights rather than delta functions, so he has |J| is in
290:    his integral. We allow this with the boolean flag.
291: */
292: static PetscErrorCode DMSwarmComputeMassMatrix_Private(DM dmc, DM dmf, Mat mass, PetscBool useDeltaFunction, void *ctx)
293: {
294:   const char  *name = "Mass Matrix";
295:   MPI_Comm     comm;
296:   PetscDS      prob;
297:   PetscSection fsection, globalFSection;
298:   PetscHSetIJ  ht;
299:   PetscLayout  rLayout, colLayout;
300:   PetscInt    *dnz, *onz;
301:   PetscInt     locRows, locCols, rStart, colStart, colEnd, *rowIDXs;
302:   PetscReal   *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0};
303:   PetscScalar *elemMat;
304:   PetscInt     dim, Nf, field, cStart, cEnd, cell, totDim, maxC = 0, totNc = 0;

306:   PetscFunctionBegin;
307:   PetscCall(PetscObjectGetComm((PetscObject)mass, &comm));
308:   PetscCall(DMGetCoordinateDim(dmf, &dim));
309:   PetscCall(DMGetDS(dmf, &prob));
310:   PetscCall(PetscDSGetNumFields(prob, &Nf));
311:   PetscCall(PetscDSGetTotalDimension(prob, &totDim));
312:   PetscCall(PetscMalloc3(dim, &v0, dim * dim, &J, dim * dim, &invJ));
313:   PetscCall(DMGetLocalSection(dmf, &fsection));
314:   PetscCall(DMGetGlobalSection(dmf, &globalFSection));
315:   PetscCall(DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd));
316:   PetscCall(MatGetLocalSize(mass, &locRows, &locCols));

318:   PetscCall(PetscLayoutCreate(comm, &colLayout));
319:   PetscCall(PetscLayoutSetLocalSize(colLayout, locCols));
320:   PetscCall(PetscLayoutSetBlockSize(colLayout, 1));
321:   PetscCall(PetscLayoutSetUp(colLayout));
322:   PetscCall(PetscLayoutGetRange(colLayout, &colStart, &colEnd));
323:   PetscCall(PetscLayoutDestroy(&colLayout));

325:   PetscCall(PetscLayoutCreate(comm, &rLayout));
326:   PetscCall(PetscLayoutSetLocalSize(rLayout, locRows));
327:   PetscCall(PetscLayoutSetBlockSize(rLayout, 1));
328:   PetscCall(PetscLayoutSetUp(rLayout));
329:   PetscCall(PetscLayoutGetRange(rLayout, &rStart, NULL));
330:   PetscCall(PetscLayoutDestroy(&rLayout));

332:   PetscCall(PetscCalloc2(locRows, &dnz, locRows, &onz));
333:   PetscCall(PetscHSetIJCreate(&ht));

335:   PetscCall(PetscSynchronizedFlush(comm, NULL));
336:   for (field = 0; field < Nf; ++field) {
337:     PetscObject  obj;
338:     PetscClassId id;
339:     PetscInt     Nc;

341:     PetscCall(PetscDSGetDiscretization(prob, field, &obj));
342:     PetscCall(PetscObjectGetClassId(obj, &id));
343:     if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
344:     else PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc));
345:     totNc += Nc;
346:   }
347:   /* count non-zeros */
348:   PetscCall(DMSwarmSortGetAccess(dmc));
349:   for (field = 0; field < Nf; ++field) {
350:     PetscObject  obj;
351:     PetscClassId id;
352:     PetscInt     Nc;

354:     PetscCall(PetscDSGetDiscretization(prob, field, &obj));
355:     PetscCall(PetscObjectGetClassId(obj, &id));
356:     if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
357:     else PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc));

359:     for (cell = cStart; cell < cEnd; ++cell) {
360:       PetscInt *findices, *cindices; /* fine is vertices, coarse is particles */
361:       PetscInt  numFIndices, numCIndices;

363:       PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
364:       PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
365:       maxC = PetscMax(maxC, numCIndices);
366:       {
367:         PetscHashIJKey key;
368:         PetscBool      missing;
369:         for (PetscInt i = 0; i < numFIndices; ++i) {
370:           key.j = findices[i]; /* global column (from Plex) */
371:           if (key.j >= 0) {
372:             /* Get indices for coarse elements */
373:             for (PetscInt j = 0; j < numCIndices; ++j) {
374:               for (PetscInt c = 0; c < Nc; ++c) {
375:                 // TODO Need field offset on particle here
376:                 key.i = cindices[j] * totNc + c + rStart; /* global cols (from Swarm) */
377:                 if (key.i < 0) continue;
378:                 PetscCall(PetscHSetIJQueryAdd(ht, key, &missing));
379:                 if (missing) {
380:                   if ((key.j >= colStart) && (key.j < colEnd)) ++dnz[key.i - rStart];
381:                   else ++onz[key.i - rStart];
382:                 } else SETERRQ(PetscObjectComm((PetscObject)dmf), PETSC_ERR_SUP, "Set new value at %" PetscInt_FMT ",%" PetscInt_FMT, key.i, key.j);
383:               }
384:             }
385:           }
386:         }
387:         PetscCall(PetscFree(cindices));
388:       }
389:       PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
390:     }
391:   }
392:   PetscCall(PetscHSetIJDestroy(&ht));
393:   PetscCall(MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL));
394:   PetscCall(MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
395:   PetscCall(PetscFree2(dnz, onz));
396:   PetscCall(PetscMalloc3(maxC * totNc * totDim, &elemMat, maxC * totNc, &rowIDXs, maxC * dim, &xi));
397:   for (field = 0; field < Nf; ++field) {
398:     PetscTabulation Tcoarse;
399:     PetscObject     obj;
400:     PetscClassId    id;
401:     PetscReal      *fieldVals;
402:     PetscInt        Nc;

404:     PetscCall(PetscDSGetDiscretization(prob, field, &obj));
405:     PetscCall(PetscObjectGetClassId(obj, &id));
406:     if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
407:     else PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc));

409:     PetscCall(DMSwarmGetField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **)&fieldVals));
410:     for (cell = cStart; cell < cEnd; ++cell) {
411:       PetscInt *findices, *cindices;
412:       PetscInt  numFIndices, numCIndices;

414:       /* TODO: Use DMField instead of assuming affine */
415:       PetscCall(DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ));
416:       PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
417:       PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
418:       for (PetscInt j = 0; j < numCIndices; ++j) CoordinatesRealToRef(dim, dim, v0ref, v0, invJ, &fieldVals[cindices[j] * dim], &xi[j * dim]);
419:       if (id == PETSCFE_CLASSID) PetscCall(PetscFECreateTabulation((PetscFE)obj, 1, numCIndices, xi, 0, &Tcoarse));
420:       else PetscCall(PetscFVCreateTabulation((PetscFV)obj, 1, numCIndices, xi, 0, &Tcoarse));
421:       /* Get elemMat entries by multiplying by weight */
422:       PetscCall(PetscArrayzero(elemMat, numCIndices * Nc * totDim));
423:       for (PetscInt i = 0; i < numFIndices / Nc; ++i) {
424:         for (PetscInt j = 0; j < numCIndices; ++j) {
425:           for (PetscInt c = 0; c < Nc; ++c) {
426:             // TODO Need field offset on particle and field here
427:             /* B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c */
428:             elemMat[(j * totNc + c) * numFIndices + i * Nc + c] += Tcoarse->T[0][(j * numFIndices + i * Nc + c) * Nc + c] * (useDeltaFunction ? 1.0 : detJ);
429:           }
430:         }
431:       }
432:       for (PetscInt j = 0; j < numCIndices; ++j)
433:         // TODO Need field offset on particle here
434:         for (PetscInt c = 0; c < Nc; ++c) rowIDXs[j * Nc + c] = cindices[j] * totNc + c + rStart;
435:       if (0) PetscCall(DMPrintCellMatrix(cell, name, numCIndices * Nc, numFIndices, elemMat));
436:       PetscCall(MatSetValues(mass, numCIndices * Nc, rowIDXs, numFIndices, findices, elemMat, ADD_VALUES));
437:       PetscCall(PetscFree(cindices));
438:       PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
439:       PetscCall(PetscTabulationDestroy(&Tcoarse));
440:     }
441:     PetscCall(DMSwarmRestoreField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **)&fieldVals));
442:   }
443:   PetscCall(PetscFree3(elemMat, rowIDXs, xi));
444:   PetscCall(DMSwarmSortRestoreAccess(dmc));
445:   PetscCall(PetscFree3(v0, J, invJ));
446:   PetscCall(MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY));
447:   PetscCall(MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY));
448:   PetscFunctionReturn(PETSC_SUCCESS);
449: }

451: /* Returns empty matrix for use with SNES FD */
452: static PetscErrorCode DMCreateMatrix_Swarm(DM sw, Mat *m)
453: {
454:   Vec      field;
455:   PetscInt size;

457:   PetscFunctionBegin;
458:   PetscCall(DMGetGlobalVector(sw, &field));
459:   PetscCall(VecGetLocalSize(field, &size));
460:   PetscCall(DMRestoreGlobalVector(sw, &field));
461:   PetscCall(MatCreate(PETSC_COMM_WORLD, m));
462:   PetscCall(MatSetFromOptions(*m));
463:   PetscCall(MatSetSizes(*m, PETSC_DECIDE, PETSC_DECIDE, size, size));
464:   PetscCall(MatSeqAIJSetPreallocation(*m, 1, NULL));
465:   PetscCall(MatZeroEntries(*m));
466:   PetscCall(MatAssemblyBegin(*m, MAT_FINAL_ASSEMBLY));
467:   PetscCall(MatAssemblyEnd(*m, MAT_FINAL_ASSEMBLY));
468:   PetscCall(MatShift(*m, 1.0));
469:   PetscCall(MatSetDM(*m, sw));
470:   PetscFunctionReturn(PETSC_SUCCESS);
471: }

473: /* FEM cols, Particle rows */
474: static PetscErrorCode DMCreateMassMatrix_Swarm(DM dmCoarse, DM dmFine, Mat *mass)
475: {
476:   PetscSection gsf;
477:   PetscInt     m, n, Np, bs = 1;
478:   void        *ctx;
479:   PetscBool    set  = ((DM_Swarm *)dmCoarse->data)->vec_field_set;
480:   char        *name = ((DM_Swarm *)dmCoarse->data)->vec_field_name;

482:   PetscFunctionBegin;
483:   PetscCall(DMGetGlobalSection(dmFine, &gsf));
484:   PetscCall(PetscSectionGetConstrainedStorageSize(gsf, &m));
485:   PetscCall(DMSwarmGetLocalSize(dmCoarse, &Np));
486:   // TODO Include all fields
487:   if (set) PetscCall(DMSwarmGetFieldInfo(dmCoarse, name, &bs, NULL));
488:   n = Np * bs;
489:   PetscCall(MatCreate(PetscObjectComm((PetscObject)dmCoarse), mass));
490:   PetscCall(MatSetSizes(*mass, n, m, PETSC_DETERMINE, PETSC_DETERMINE));
491:   PetscCall(MatSetType(*mass, dmCoarse->mattype));
492:   PetscCall(DMGetApplicationContext(dmFine, &ctx));

494:   PetscCall(DMSwarmComputeMassMatrix_Private(dmCoarse, dmFine, *mass, PETSC_TRUE, ctx));
495:   PetscCall(MatViewFromOptions(*mass, NULL, "-mass_mat_view"));
496:   PetscFunctionReturn(PETSC_SUCCESS);
497: }

499: static PetscErrorCode DMSwarmComputeMassMatrixSquare_Private(DM dmc, DM dmf, Mat mass, PetscBool useDeltaFunction, void *ctx)
500: {
501:   const char  *name = "Mass Matrix Square";
502:   MPI_Comm     comm;
503:   PetscDS      prob;
504:   PetscSection fsection, globalFSection;
505:   PetscHSetIJ  ht;
506:   PetscLayout  rLayout, colLayout;
507:   PetscInt    *dnz, *onz, *adj, depth, maxConeSize, maxSupportSize, maxAdjSize;
508:   PetscInt     locRows, locCols, rStart, colStart, colEnd, *rowIDXs;
509:   PetscReal   *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0};
510:   PetscScalar *elemMat, *elemMatSq;
511:   PetscInt     cdim, Nf, field, cStart, cEnd, cell, totDim, maxC = 0;

513:   PetscFunctionBegin;
514:   PetscCall(PetscObjectGetComm((PetscObject)mass, &comm));
515:   PetscCall(DMGetCoordinateDim(dmf, &cdim));
516:   PetscCall(DMGetDS(dmf, &prob));
517:   PetscCall(PetscDSGetNumFields(prob, &Nf));
518:   PetscCall(PetscDSGetTotalDimension(prob, &totDim));
519:   PetscCall(PetscMalloc3(cdim, &v0, cdim * cdim, &J, cdim * cdim, &invJ));
520:   PetscCall(DMGetLocalSection(dmf, &fsection));
521:   PetscCall(DMGetGlobalSection(dmf, &globalFSection));
522:   PetscCall(DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd));
523:   PetscCall(MatGetLocalSize(mass, &locRows, &locCols));

525:   PetscCall(PetscLayoutCreate(comm, &colLayout));
526:   PetscCall(PetscLayoutSetLocalSize(colLayout, locCols));
527:   PetscCall(PetscLayoutSetBlockSize(colLayout, 1));
528:   PetscCall(PetscLayoutSetUp(colLayout));
529:   PetscCall(PetscLayoutGetRange(colLayout, &colStart, &colEnd));
530:   PetscCall(PetscLayoutDestroy(&colLayout));

532:   PetscCall(PetscLayoutCreate(comm, &rLayout));
533:   PetscCall(PetscLayoutSetLocalSize(rLayout, locRows));
534:   PetscCall(PetscLayoutSetBlockSize(rLayout, 1));
535:   PetscCall(PetscLayoutSetUp(rLayout));
536:   PetscCall(PetscLayoutGetRange(rLayout, &rStart, NULL));
537:   PetscCall(PetscLayoutDestroy(&rLayout));

539:   PetscCall(DMPlexGetDepth(dmf, &depth));
540:   PetscCall(DMPlexGetMaxSizes(dmf, &maxConeSize, &maxSupportSize));
541:   maxAdjSize = PetscPowInt(maxConeSize * maxSupportSize, depth);
542:   PetscCall(PetscMalloc1(maxAdjSize, &adj));

544:   PetscCall(PetscCalloc2(locRows, &dnz, locRows, &onz));
545:   PetscCall(PetscHSetIJCreate(&ht));
546:   /* Count nonzeros
547:        This is just FVM++, but we cannot use the Plex P0 allocation since unknowns in a cell will not be contiguous
548:   */
549:   PetscCall(DMSwarmSortGetAccess(dmc));
550:   for (cell = cStart; cell < cEnd; ++cell) {
551:     PetscInt  i;
552:     PetscInt *cindices;
553:     PetscInt  numCIndices;
554: #if 0
555:     PetscInt  adjSize = maxAdjSize, a, j;
556: #endif

558:     PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
559:     maxC = PetscMax(maxC, numCIndices);
560:     /* Diagonal block */
561:     for (i = 0; i < numCIndices; ++i) dnz[cindices[i]] += numCIndices;
562: #if 0
563:     /* Off-diagonal blocks */
564:     PetscCall(DMPlexGetAdjacency(dmf, cell, &adjSize, &adj));
565:     for (a = 0; a < adjSize; ++a) {
566:       if (adj[a] >= cStart && adj[a] < cEnd && adj[a] != cell) {
567:         const PetscInt ncell = adj[a];
568:         PetscInt      *ncindices;
569:         PetscInt       numNCIndices;

571:         PetscCall(DMSwarmSortGetPointsPerCell(dmc, ncell, &numNCIndices, &ncindices));
572:         {
573:           PetscHashIJKey key;
574:           PetscBool      missing;

576:           for (i = 0; i < numCIndices; ++i) {
577:             key.i = cindices[i] + rStart; /* global rows (from Swarm) */
578:             if (key.i < 0) continue;
579:             for (j = 0; j < numNCIndices; ++j) {
580:               key.j = ncindices[j] + rStart; /* global column (from Swarm) */
581:               if (key.j < 0) continue;
582:               PetscCall(PetscHSetIJQueryAdd(ht, key, &missing));
583:               if (missing) {
584:                 if ((key.j >= colStart) && (key.j < colEnd)) ++dnz[key.i - rStart];
585:                 else                                         ++onz[key.i - rStart];
586:               }
587:             }
588:           }
589:         }
590:         PetscCall(PetscFree(ncindices));
591:       }
592:     }
593: #endif
594:     PetscCall(PetscFree(cindices));
595:   }
596:   PetscCall(PetscHSetIJDestroy(&ht));
597:   PetscCall(MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL));
598:   PetscCall(MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
599:   PetscCall(PetscFree2(dnz, onz));
600:   PetscCall(PetscMalloc4(maxC * totDim, &elemMat, maxC * maxC, &elemMatSq, maxC, &rowIDXs, maxC * cdim, &xi));
601:   /* Fill in values
602:        Each entry is a sum of terms \phi_i(x_p) \phi_i(x_q)
603:        Start just by producing block diagonal
604:        Could loop over adjacent cells
605:          Produce neighboring element matrix
606:          TODO Determine which columns and rows correspond to shared dual vector
607:          Do MatMatMult with rectangular matrices
608:          Insert block
609:   */
610:   for (field = 0; field < Nf; ++field) {
611:     PetscTabulation Tcoarse;
612:     PetscObject     obj;
613:     PetscReal      *coords;
614:     PetscInt        Nc, i;

616:     PetscCall(PetscDSGetDiscretization(prob, field, &obj));
617:     PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
618:     PetscCheck(Nc == 1, PetscObjectComm((PetscObject)dmf), PETSC_ERR_SUP, "Can only interpolate a scalar field from particles, Nc = %" PetscInt_FMT, Nc);
619:     PetscCall(DMSwarmGetField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
620:     for (cell = cStart; cell < cEnd; ++cell) {
621:       PetscInt *findices, *cindices;
622:       PetscInt  numFIndices, numCIndices;
623:       PetscInt  p, c;

625:       /* TODO: Use DMField instead of assuming affine */
626:       PetscCall(DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ));
627:       PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
628:       PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
629:       for (p = 0; p < numCIndices; ++p) CoordinatesRealToRef(cdim, cdim, v0ref, v0, invJ, &coords[cindices[p] * cdim], &xi[p * cdim]);
630:       PetscCall(PetscFECreateTabulation((PetscFE)obj, 1, numCIndices, xi, 0, &Tcoarse));
631:       /* Get elemMat entries by multiplying by weight */
632:       PetscCall(PetscArrayzero(elemMat, numCIndices * totDim));
633:       for (i = 0; i < numFIndices; ++i) {
634:         for (p = 0; p < numCIndices; ++p) {
635:           for (c = 0; c < Nc; ++c) {
636:             /* B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c */
637:             elemMat[p * numFIndices + i] += Tcoarse->T[0][(p * numFIndices + i) * Nc + c] * (useDeltaFunction ? 1.0 : detJ);
638:           }
639:         }
640:       }
641:       PetscCall(PetscTabulationDestroy(&Tcoarse));
642:       for (p = 0; p < numCIndices; ++p) rowIDXs[p] = cindices[p] + rStart;
643:       if (0) PetscCall(DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat));
644:       /* Block diagonal */
645:       if (numCIndices) {
646:         PetscBLASInt blasn, blask;
647:         PetscScalar  one = 1.0, zero = 0.0;

649:         PetscCall(PetscBLASIntCast(numCIndices, &blasn));
650:         PetscCall(PetscBLASIntCast(numFIndices, &blask));
651:         PetscCallBLAS("BLASgemm", BLASgemm_("T", "N", &blasn, &blasn, &blask, &one, elemMat, &blask, elemMat, &blask, &zero, elemMatSq, &blasn));
652:       }
653:       PetscCall(MatSetValues(mass, numCIndices, rowIDXs, numCIndices, rowIDXs, elemMatSq, ADD_VALUES));
654:       /* TODO off-diagonal */
655:       PetscCall(PetscFree(cindices));
656:       PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
657:     }
658:     PetscCall(DMSwarmRestoreField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
659:   }
660:   PetscCall(PetscFree4(elemMat, elemMatSq, rowIDXs, xi));
661:   PetscCall(PetscFree(adj));
662:   PetscCall(DMSwarmSortRestoreAccess(dmc));
663:   PetscCall(PetscFree3(v0, J, invJ));
664:   PetscCall(MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY));
665:   PetscCall(MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY));
666:   PetscFunctionReturn(PETSC_SUCCESS);
667: }

669: /*@
670:   DMSwarmCreateMassMatrixSquare - Creates the block-diagonal of the square, M^T_p M_p, of the particle mass matrix M_p

672:   Collective

674:   Input Parameters:
675: + dmCoarse - a `DMSWARM`
676: - dmFine   - a `DMPLEX`

678:   Output Parameter:
679: . mass - the square of the particle mass matrix

681:   Level: advanced

683:   Note:
684:   We only compute the block diagonal since this provides a good preconditioner and is completely local. It would be possible in the
685:   future to compute the full normal equations.

687: .seealso: `DM`, `DMSWARM`, `DMCreateMassMatrix()`
688: @*/
689: PetscErrorCode DMSwarmCreateMassMatrixSquare(DM dmCoarse, DM dmFine, Mat *mass)
690: {
691:   PetscInt n;
692:   void    *ctx;

694:   PetscFunctionBegin;
695:   PetscCall(DMSwarmGetLocalSize(dmCoarse, &n));
696:   PetscCall(MatCreate(PetscObjectComm((PetscObject)dmCoarse), mass));
697:   PetscCall(MatSetSizes(*mass, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
698:   PetscCall(MatSetType(*mass, dmCoarse->mattype));
699:   PetscCall(DMGetApplicationContext(dmFine, &ctx));

701:   PetscCall(DMSwarmComputeMassMatrixSquare_Private(dmCoarse, dmFine, *mass, PETSC_TRUE, ctx));
702:   PetscCall(MatViewFromOptions(*mass, NULL, "-mass_sq_mat_view"));
703:   PetscFunctionReturn(PETSC_SUCCESS);
704: }

706: /*@
707:   DMSwarmCreateGlobalVectorFromField - Creates a `Vec` object sharing the array associated with a given field

709:   Collective

711:   Input Parameters:
712: + dm        - a `DMSWARM`
713: - fieldname - the textual name given to a registered field

715:   Output Parameter:
716: . vec - the vector

718:   Level: beginner

720:   Note:
721:   The vector must be returned using a matching call to `DMSwarmDestroyGlobalVectorFromField()`.

723: .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmDestroyGlobalVectorFromField()`
724: @*/
725: PetscErrorCode DMSwarmCreateGlobalVectorFromField(DM dm, const char fieldname[], Vec *vec)
726: {
727:   MPI_Comm comm = PetscObjectComm((PetscObject)dm);

729:   PetscFunctionBegin;
731:   PetscCall(DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec));
732:   PetscFunctionReturn(PETSC_SUCCESS);
733: }

735: /*@
736:   DMSwarmDestroyGlobalVectorFromField - Destroys the `Vec` object which share the array associated with a given field

738:   Collective

740:   Input Parameters:
741: + dm        - a `DMSWARM`
742: - fieldname - the textual name given to a registered field

744:   Output Parameter:
745: . vec - the vector

747:   Level: beginner

749: .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmCreateGlobalVectorFromField()`
750: @*/
751: PetscErrorCode DMSwarmDestroyGlobalVectorFromField(DM dm, const char fieldname[], Vec *vec)
752: {
753:   PetscFunctionBegin;
755:   PetscCall(DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec));
756:   PetscFunctionReturn(PETSC_SUCCESS);
757: }

759: /*@
760:   DMSwarmCreateLocalVectorFromField - Creates a `Vec` object sharing the array associated with a given field

762:   Collective

764:   Input Parameters:
765: + dm        - a `DMSWARM`
766: - fieldname - the textual name given to a registered field

768:   Output Parameter:
769: . vec - the vector

771:   Level: beginner

773:   Note:
774:   The vector must be returned using a matching call to DMSwarmDestroyLocalVectorFromField().

776: .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmDestroyLocalVectorFromField()`
777: @*/
778: PetscErrorCode DMSwarmCreateLocalVectorFromField(DM dm, const char fieldname[], Vec *vec)
779: {
780:   MPI_Comm comm = PETSC_COMM_SELF;

782:   PetscFunctionBegin;
783:   PetscCall(DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec));
784:   PetscFunctionReturn(PETSC_SUCCESS);
785: }

787: /*@
788:   DMSwarmDestroyLocalVectorFromField - Destroys the `Vec` object which share the array associated with a given field

790:   Collective

792:   Input Parameters:
793: + dm        - a `DMSWARM`
794: - fieldname - the textual name given to a registered field

796:   Output Parameter:
797: . vec - the vector

799:   Level: beginner

801: .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmCreateLocalVectorFromField()`
802: @*/
803: PetscErrorCode DMSwarmDestroyLocalVectorFromField(DM dm, const char fieldname[], Vec *vec)
804: {
805:   PetscFunctionBegin;
807:   PetscCall(DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec));
808:   PetscFunctionReturn(PETSC_SUCCESS);
809: }

811: /*@
812:   DMSwarmInitializeFieldRegister - Initiates the registration of fields to a `DMSWARM`

814:   Collective

816:   Input Parameter:
817: . dm - a `DMSWARM`

819:   Level: beginner

821:   Note:
822:   After all fields have been registered, you must call `DMSwarmFinalizeFieldRegister()`.

824: .seealso: `DM`, `DMSWARM`, `DMSwarmFinalizeFieldRegister()`, `DMSwarmRegisterPetscDatatypeField()`,
825:           `DMSwarmRegisterUserStructField()`, `DMSwarmRegisterUserDatatypeField()`
826: @*/
827: PetscErrorCode DMSwarmInitializeFieldRegister(DM dm)
828: {
829:   DM_Swarm *swarm = (DM_Swarm *)dm->data;

831:   PetscFunctionBegin;
832:   if (!swarm->field_registration_initialized) {
833:     swarm->field_registration_initialized = PETSC_TRUE;
834:     PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmField_pid, 1, PETSC_INT64)); /* unique identifier */
835:     PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmField_rank, 1, PETSC_INT));  /* used for communication */
836:   }
837:   PetscFunctionReturn(PETSC_SUCCESS);
838: }

840: /*@
841:   DMSwarmFinalizeFieldRegister - Finalizes the registration of fields to a `DMSWARM`

843:   Collective

845:   Input Parameter:
846: . dm - a `DMSWARM`

848:   Level: beginner

850:   Note:
851:   After `DMSwarmFinalizeFieldRegister()` has been called, no new fields can be defined on the `DMSWARM`.

853: .seealso: `DM`, `DMSWARM`, `DMSwarmInitializeFieldRegister()`, `DMSwarmRegisterPetscDatatypeField()`,
854:           `DMSwarmRegisterUserStructField()`, `DMSwarmRegisterUserDatatypeField()`
855: @*/
856: PetscErrorCode DMSwarmFinalizeFieldRegister(DM dm)
857: {
858:   DM_Swarm *swarm = (DM_Swarm *)dm->data;

860:   PetscFunctionBegin;
861:   if (!swarm->field_registration_finalized) PetscCall(DMSwarmDataBucketFinalize(swarm->db));
862:   swarm->field_registration_finalized = PETSC_TRUE;
863:   PetscFunctionReturn(PETSC_SUCCESS);
864: }

866: /*@
867:   DMSwarmSetLocalSizes - Sets the length of all registered fields on the `DMSWARM`

869:   Not Collective

871:   Input Parameters:
872: + dm     - a `DMSWARM`
873: . nlocal - the length of each registered field
874: - buffer - the length of the buffer used to efficient dynamic re-sizing

876:   Level: beginner

878: .seealso: `DM`, `DMSWARM`, `DMSwarmGetLocalSize()`
879: @*/
880: PetscErrorCode DMSwarmSetLocalSizes(DM dm, PetscInt nlocal, PetscInt buffer)
881: {
882:   DM_Swarm *swarm = (DM_Swarm *)dm->data;

884:   PetscFunctionBegin;
885:   PetscCall(PetscLogEventBegin(DMSWARM_SetSizes, 0, 0, 0, 0));
886:   PetscCall(DMSwarmDataBucketSetSizes(swarm->db, nlocal, buffer));
887:   PetscCall(PetscLogEventEnd(DMSWARM_SetSizes, 0, 0, 0, 0));
888:   PetscFunctionReturn(PETSC_SUCCESS);
889: }

891: /*@
892:   DMSwarmSetCellDM - Attaches a `DM` to a `DMSWARM`

894:   Collective

896:   Input Parameters:
897: + dm     - a `DMSWARM`
898: - dmcell - the `DM` to attach to the `DMSWARM`

900:   Level: beginner

902:   Note:
903:   The attached `DM` (dmcell) will be queried for point location and
904:   neighbor MPI-rank information if `DMSwarmMigrate()` is called.

906: .seealso: `DM`, `DMSWARM`, `DMSwarmSetType()`, `DMSwarmGetCellDM()`, `DMSwarmMigrate()`
907: @*/
908: PetscErrorCode DMSwarmSetCellDM(DM dm, DM dmcell)
909: {
910:   DM_Swarm *swarm = (DM_Swarm *)dm->data;

912:   PetscFunctionBegin;
913:   swarm->dmcell = dmcell;
914:   PetscFunctionReturn(PETSC_SUCCESS);
915: }

917: /*@
918:   DMSwarmGetCellDM - Fetches the attached cell `DM`

920:   Collective

922:   Input Parameter:
923: . dm - a `DMSWARM`

925:   Output Parameter:
926: . dmcell - the `DM` which was attached to the `DMSWARM`

928:   Level: beginner

930: .seealso: `DM`, `DMSWARM`, `DMSwarmSetCellDM()`
931: @*/
932: PetscErrorCode DMSwarmGetCellDM(DM dm, DM *dmcell)
933: {
934:   DM_Swarm *swarm = (DM_Swarm *)dm->data;

936:   PetscFunctionBegin;
937:   *dmcell = swarm->dmcell;
938:   PetscFunctionReturn(PETSC_SUCCESS);
939: }

941: /*@
942:   DMSwarmGetLocalSize - Retrieves the local length of fields registered

944:   Not Collective

946:   Input Parameter:
947: . dm - a `DMSWARM`

949:   Output Parameter:
950: . nlocal - the length of each registered field

952:   Level: beginner

954: .seealso: `DM`, `DMSWARM`, `DMSwarmGetSize()`, `DMSwarmSetLocalSizes()`
955: @*/
956: PetscErrorCode DMSwarmGetLocalSize(DM dm, PetscInt *nlocal)
957: {
958:   DM_Swarm *swarm = (DM_Swarm *)dm->data;

960:   PetscFunctionBegin;
961:   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, nlocal, NULL, NULL));
962:   PetscFunctionReturn(PETSC_SUCCESS);
963: }

965: /*@
966:   DMSwarmGetSize - Retrieves the total length of fields registered

968:   Collective

970:   Input Parameter:
971: . dm - a `DMSWARM`

973:   Output Parameter:
974: . n - the total length of each registered field

976:   Level: beginner

978:   Note:
979:   This calls `MPI_Allreduce()` upon each call (inefficient but safe)

981: .seealso: `DM`, `DMSWARM`, `DMSwarmGetLocalSize()`, `DMSwarmSetLocalSizes()`
982: @*/
983: PetscErrorCode DMSwarmGetSize(DM dm, PetscInt *n)
984: {
985:   DM_Swarm *swarm = (DM_Swarm *)dm->data;
986:   PetscInt  nlocal;

988:   PetscFunctionBegin;
989:   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &nlocal, NULL, NULL));
990:   PetscCall(MPIU_Allreduce(&nlocal, n, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm)));
991:   PetscFunctionReturn(PETSC_SUCCESS);
992: }

994: /*@
995:   DMSwarmRegisterPetscDatatypeField - Register a field to a `DMSWARM` with a native PETSc data type

997:   Collective

999:   Input Parameters:
1000: + dm        - a `DMSWARM`
1001: . fieldname - the textual name to identify this field
1002: . blocksize - the number of each data type
1003: - type      - a valid PETSc data type (`PETSC_CHAR`, `PETSC_SHORT`, `PETSC_INT`, `PETSC_FLOAT`, `PETSC_REAL`, `PETSC_LONG`)

1005:   Level: beginner

1007:   Notes:
1008:   The textual name for each registered field must be unique.

1010: .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterUserStructField()`, `DMSwarmRegisterUserDatatypeField()`
1011: @*/
1012: PetscErrorCode DMSwarmRegisterPetscDatatypeField(DM dm, const char fieldname[], PetscInt blocksize, PetscDataType type)
1013: {
1014:   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1015:   size_t    size;

1017:   PetscFunctionBegin;
1018:   PetscCheck(swarm->field_registration_initialized, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMSwarmInitializeFieldRegister() first");
1019:   PetscCheck(!swarm->field_registration_finalized, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Cannot register additional fields after calling DMSwarmFinalizeFieldRegister() first");

1021:   PetscCheck(type != PETSC_OBJECT, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
1022:   PetscCheck(type != PETSC_FUNCTION, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
1023:   PetscCheck(type != PETSC_STRING, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
1024:   PetscCheck(type != PETSC_STRUCT, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
1025:   PetscCheck(type != PETSC_DATATYPE_UNKNOWN, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");

1027:   PetscCall(PetscDataTypeGetSize(type, &size));
1028:   /* Load a specific data type into data bucket, specifying textual name and its size in bytes */
1029:   PetscCall(DMSwarmDataBucketRegisterField(swarm->db, "DMSwarmRegisterPetscDatatypeField", fieldname, blocksize * size, NULL));
1030:   {
1031:     DMSwarmDataField gfield;

1033:     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
1034:     PetscCall(DMSwarmDataFieldSetBlockSize(gfield, blocksize));
1035:   }
1036:   swarm->db->field[swarm->db->nfields - 1]->petsc_type = type;
1037:   PetscFunctionReturn(PETSC_SUCCESS);
1038: }

1040: /*@C
1041:   DMSwarmRegisterUserStructField - Register a user defined struct to a `DMSWARM`

1043:   Collective

1045:   Input Parameters:
1046: + dm        - a `DMSWARM`
1047: . fieldname - the textual name to identify this field
1048: - size      - the size in bytes of the user struct of each data type

1050:   Level: beginner

1052:   Note:
1053:   The textual name for each registered field must be unique.

1055: .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmRegisterUserDatatypeField()`
1056: @*/
1057: PetscErrorCode DMSwarmRegisterUserStructField(DM dm, const char fieldname[], size_t size)
1058: {
1059:   DM_Swarm *swarm = (DM_Swarm *)dm->data;

1061:   PetscFunctionBegin;
1062:   PetscCall(DMSwarmDataBucketRegisterField(swarm->db, "DMSwarmRegisterUserStructField", fieldname, size, NULL));
1063:   swarm->db->field[swarm->db->nfields - 1]->petsc_type = PETSC_STRUCT;
1064:   PetscFunctionReturn(PETSC_SUCCESS);
1065: }

1067: /*@C
1068:   DMSwarmRegisterUserDatatypeField - Register a user defined data type to a `DMSWARM`

1070:   Collective

1072:   Input Parameters:
1073: + dm        - a `DMSWARM`
1074: . fieldname - the textual name to identify this field
1075: . size      - the size in bytes of the user data type
1076: - blocksize - the number of each data type

1078:   Level: beginner

1080:   Note:
1081:   The textual name for each registered field must be unique.

1083: .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmRegisterUserStructField()`
1084: @*/
1085: PetscErrorCode DMSwarmRegisterUserDatatypeField(DM dm, const char fieldname[], size_t size, PetscInt blocksize)
1086: {
1087:   DM_Swarm *swarm = (DM_Swarm *)dm->data;

1089:   PetscFunctionBegin;
1090:   PetscCall(DMSwarmDataBucketRegisterField(swarm->db, "DMSwarmRegisterUserDatatypeField", fieldname, blocksize * size, NULL));
1091:   {
1092:     DMSwarmDataField gfield;

1094:     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
1095:     PetscCall(DMSwarmDataFieldSetBlockSize(gfield, blocksize));
1096:   }
1097:   swarm->db->field[swarm->db->nfields - 1]->petsc_type = PETSC_DATATYPE_UNKNOWN;
1098:   PetscFunctionReturn(PETSC_SUCCESS);
1099: }

1101: /*@C
1102:   DMSwarmGetField - Get access to the underlying array storing all entries associated with a registered field

1104:   Not Collective, No Fortran Support

1106:   Input Parameters:
1107: + dm        - a `DMSWARM`
1108: - fieldname - the textual name to identify this field

1110:   Output Parameters:
1111: + blocksize - the number of each data type
1112: . type      - the data type
1113: - data      - pointer to raw array

1115:   Level: beginner

1117:   Notes:
1118:   The array must be returned using a matching call to `DMSwarmRestoreField()`.

1120: .seealso: `DM`, `DMSWARM`, `DMSwarmRestoreField()`
1121: @*/
1122: PetscErrorCode DMSwarmGetField(DM dm, const char fieldname[], PetscInt *blocksize, PetscDataType *type, void **data)
1123: {
1124:   DM_Swarm        *swarm = (DM_Swarm *)dm->data;
1125:   DMSwarmDataField gfield;

1127:   PetscFunctionBegin;
1129:   if (!swarm->issetup) PetscCall(DMSetUp(dm));
1130:   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
1131:   PetscCall(DMSwarmDataFieldGetAccess(gfield));
1132:   PetscCall(DMSwarmDataFieldGetEntries(gfield, data));
1133:   if (blocksize) *blocksize = gfield->bs;
1134:   if (type) *type = gfield->petsc_type;
1135:   PetscFunctionReturn(PETSC_SUCCESS);
1136: }

1138: /*@C
1139:   DMSwarmRestoreField - Restore access to the underlying array storing all entries associated with a registered field

1141:   Not Collective, No Fortran Support

1143:   Input Parameters:
1144: + dm        - a `DMSWARM`
1145: - fieldname - the textual name to identify this field

1147:   Output Parameters:
1148: + blocksize - the number of each data type
1149: . type      - the data type
1150: - data      - pointer to raw array

1152:   Level: beginner

1154:   Notes:
1155:   The user must call `DMSwarmGetField()` prior to calling `DMSwarmRestoreField()`.

1157: .seealso: `DM`, `DMSWARM`, `DMSwarmGetField()`
1158: @*/
1159: PetscErrorCode DMSwarmRestoreField(DM dm, const char fieldname[], PetscInt *blocksize, PetscDataType *type, void **data)
1160: {
1161:   DM_Swarm        *swarm = (DM_Swarm *)dm->data;
1162:   DMSwarmDataField gfield;

1164:   PetscFunctionBegin;
1166:   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
1167:   PetscCall(DMSwarmDataFieldRestoreAccess(gfield));
1168:   if (data) *data = NULL;
1169:   PetscFunctionReturn(PETSC_SUCCESS);
1170: }

1172: PetscErrorCode DMSwarmGetFieldInfo(DM dm, const char fieldname[], PetscInt *blocksize, PetscDataType *type)
1173: {
1174:   DM_Swarm        *swarm = (DM_Swarm *)dm->data;
1175:   DMSwarmDataField gfield;

1177:   PetscFunctionBegin;
1179:   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
1180:   if (blocksize) *blocksize = gfield->bs;
1181:   if (type) *type = gfield->petsc_type;
1182:   PetscFunctionReturn(PETSC_SUCCESS);
1183: }

1185: /*@
1186:   DMSwarmAddPoint - Add space for one new point in the `DMSWARM`

1188:   Not Collective

1190:   Input Parameter:
1191: . dm - a `DMSWARM`

1193:   Level: beginner

1195:   Notes:
1196:   The new point will have all fields initialized to zero.

1198: .seealso: `DM`, `DMSWARM`, `DMSwarmAddNPoints()`
1199: @*/
1200: PetscErrorCode DMSwarmAddPoint(DM dm)
1201: {
1202:   DM_Swarm *swarm = (DM_Swarm *)dm->data;

1204:   PetscFunctionBegin;
1205:   if (!swarm->issetup) PetscCall(DMSetUp(dm));
1206:   PetscCall(PetscLogEventBegin(DMSWARM_AddPoints, 0, 0, 0, 0));
1207:   PetscCall(DMSwarmDataBucketAddPoint(swarm->db));
1208:   PetscCall(PetscLogEventEnd(DMSWARM_AddPoints, 0, 0, 0, 0));
1209:   PetscFunctionReturn(PETSC_SUCCESS);
1210: }

1212: /*@
1213:   DMSwarmAddNPoints - Add space for a number of new points in the `DMSWARM`

1215:   Not Collective

1217:   Input Parameters:
1218: + dm      - a `DMSWARM`
1219: - npoints - the number of new points to add

1221:   Level: beginner

1223:   Notes:
1224:   The new point will have all fields initialized to zero.

1226: .seealso: `DM`, `DMSWARM`, `DMSwarmAddPoint()`
1227: @*/
1228: PetscErrorCode DMSwarmAddNPoints(DM dm, PetscInt npoints)
1229: {
1230:   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1231:   PetscInt  nlocal;

1233:   PetscFunctionBegin;
1234:   PetscCall(PetscLogEventBegin(DMSWARM_AddPoints, 0, 0, 0, 0));
1235:   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &nlocal, NULL, NULL));
1236:   nlocal = nlocal + npoints;
1237:   PetscCall(DMSwarmDataBucketSetSizes(swarm->db, nlocal, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
1238:   PetscCall(PetscLogEventEnd(DMSWARM_AddPoints, 0, 0, 0, 0));
1239:   PetscFunctionReturn(PETSC_SUCCESS);
1240: }

1242: /*@
1243:   DMSwarmRemovePoint - Remove the last point from the `DMSWARM`

1245:   Not Collective

1247:   Input Parameter:
1248: . dm - a `DMSWARM`

1250:   Level: beginner

1252: .seealso: `DM`, `DMSWARM`, `DMSwarmRemovePointAtIndex()`
1253: @*/
1254: PetscErrorCode DMSwarmRemovePoint(DM dm)
1255: {
1256:   DM_Swarm *swarm = (DM_Swarm *)dm->data;

1258:   PetscFunctionBegin;
1259:   PetscCall(PetscLogEventBegin(DMSWARM_RemovePoints, 0, 0, 0, 0));
1260:   PetscCall(DMSwarmDataBucketRemovePoint(swarm->db));
1261:   PetscCall(PetscLogEventEnd(DMSWARM_RemovePoints, 0, 0, 0, 0));
1262:   PetscFunctionReturn(PETSC_SUCCESS);
1263: }

1265: /*@
1266:   DMSwarmRemovePointAtIndex - Removes a specific point from the `DMSWARM`

1268:   Not Collective

1270:   Input Parameters:
1271: + dm  - a `DMSWARM`
1272: - idx - index of point to remove

1274:   Level: beginner

1276: .seealso: `DM`, `DMSWARM`, `DMSwarmRemovePoint()`
1277: @*/
1278: PetscErrorCode DMSwarmRemovePointAtIndex(DM dm, PetscInt idx)
1279: {
1280:   DM_Swarm *swarm = (DM_Swarm *)dm->data;

1282:   PetscFunctionBegin;
1283:   PetscCall(PetscLogEventBegin(DMSWARM_RemovePoints, 0, 0, 0, 0));
1284:   PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db, idx));
1285:   PetscCall(PetscLogEventEnd(DMSWARM_RemovePoints, 0, 0, 0, 0));
1286:   PetscFunctionReturn(PETSC_SUCCESS);
1287: }

1289: /*@
1290:   DMSwarmCopyPoint - Copy point pj to point pi in the `DMSWARM`

1292:   Not Collective

1294:   Input Parameters:
1295: + dm - a `DMSWARM`
1296: . pi - the index of the point to copy
1297: - pj - the point index where the copy should be located

1299:   Level: beginner

1301: .seealso: `DM`, `DMSWARM`, `DMSwarmRemovePoint()`
1302: @*/
1303: PetscErrorCode DMSwarmCopyPoint(DM dm, PetscInt pi, PetscInt pj)
1304: {
1305:   DM_Swarm *swarm = (DM_Swarm *)dm->data;

1307:   PetscFunctionBegin;
1308:   if (!swarm->issetup) PetscCall(DMSetUp(dm));
1309:   PetscCall(DMSwarmDataBucketCopyPoint(swarm->db, pi, swarm->db, pj));
1310:   PetscFunctionReturn(PETSC_SUCCESS);
1311: }

1313: static PetscErrorCode DMSwarmMigrate_Basic(DM dm, PetscBool remove_sent_points)
1314: {
1315:   PetscFunctionBegin;
1316:   PetscCall(DMSwarmMigrate_Push_Basic(dm, remove_sent_points));
1317:   PetscFunctionReturn(PETSC_SUCCESS);
1318: }

1320: /*@
1321:   DMSwarmMigrate - Relocates points defined in the `DMSWARM` to other MPI-ranks

1323:   Collective

1325:   Input Parameters:
1326: + dm                 - the `DMSWARM`
1327: - remove_sent_points - flag indicating if sent points should be removed from the current MPI-rank

1329:   Level: advanced

1331:   Notes:
1332:   The `DM` will be modified to accommodate received points.
1333:   If `remove_sent_points` is `PETSC_TRUE`, any points that were sent will be removed from the `DM`.
1334:   Different styles of migration are supported. See `DMSwarmSetMigrateType()`.

1336: .seealso: `DM`, `DMSWARM`, `DMSwarmSetMigrateType()`
1337: @*/
1338: PetscErrorCode DMSwarmMigrate(DM dm, PetscBool remove_sent_points)
1339: {
1340:   DM_Swarm *swarm = (DM_Swarm *)dm->data;

1342:   PetscFunctionBegin;
1343:   PetscCall(PetscLogEventBegin(DMSWARM_Migrate, 0, 0, 0, 0));
1344:   switch (swarm->migrate_type) {
1345:   case DMSWARM_MIGRATE_BASIC:
1346:     PetscCall(DMSwarmMigrate_Basic(dm, remove_sent_points));
1347:     break;
1348:   case DMSWARM_MIGRATE_DMCELLNSCATTER:
1349:     PetscCall(DMSwarmMigrate_CellDMScatter(dm, remove_sent_points));
1350:     break;
1351:   case DMSWARM_MIGRATE_DMCELLEXACT:
1352:     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_MIGRATE_DMCELLEXACT not implemented");
1353:   case DMSWARM_MIGRATE_USER:
1354:     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_MIGRATE_USER not implemented");
1355:   default:
1356:     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_MIGRATE type unknown");
1357:   }
1358:   PetscCall(PetscLogEventEnd(DMSWARM_Migrate, 0, 0, 0, 0));
1359:   PetscCall(DMClearGlobalVectors(dm));
1360:   PetscFunctionReturn(PETSC_SUCCESS);
1361: }

1363: PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm, PetscInt *globalsize);

1365: /*
1366:  DMSwarmCollectViewCreate

1368:  * Applies a collection method and gathers point neighbour points into dm

1370:  Notes:
1371:  Users should call DMSwarmCollectViewDestroy() after
1372:  they have finished computations associated with the collected points
1373: */

1375: /*@
1376:   DMSwarmCollectViewCreate - Applies a collection method and gathers points
1377:   in neighbour ranks into the `DMSWARM`

1379:   Collective

1381:   Input Parameter:
1382: . dm - the `DMSWARM`

1384:   Level: advanced

1386:   Notes:
1387:   Users should call `DMSwarmCollectViewDestroy()` after
1388:   they have finished computations associated with the collected points

1390:   Different collect methods are supported. See `DMSwarmSetCollectType()`.

1392:   Developer Note:
1393:   Create and Destroy routines create new objects that can get destroyed, they do not change the state
1394:   of the current object.

1396: .seealso: `DM`, `DMSWARM`, `DMSwarmCollectViewDestroy()`, `DMSwarmSetCollectType()`
1397: @*/
1398: PetscErrorCode DMSwarmCollectViewCreate(DM dm)
1399: {
1400:   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1401:   PetscInt  ng;

1403:   PetscFunctionBegin;
1404:   PetscCheck(!swarm->collect_view_active, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "CollectView currently active");
1405:   PetscCall(DMSwarmGetLocalSize(dm, &ng));
1406:   switch (swarm->collect_type) {
1407:   case DMSWARM_COLLECT_BASIC:
1408:     PetscCall(DMSwarmMigrate_GlobalToLocal_Basic(dm, &ng));
1409:     break;
1410:   case DMSWARM_COLLECT_DMDABOUNDINGBOX:
1411:     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_COLLECT_DMDABOUNDINGBOX not implemented");
1412:   case DMSWARM_COLLECT_GENERAL:
1413:     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_COLLECT_GENERAL not implemented");
1414:   default:
1415:     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_COLLECT type unknown");
1416:   }
1417:   swarm->collect_view_active       = PETSC_TRUE;
1418:   swarm->collect_view_reset_nlocal = ng;
1419:   PetscFunctionReturn(PETSC_SUCCESS);
1420: }

1422: /*@
1423:   DMSwarmCollectViewDestroy - Resets the `DMSWARM` to the size prior to calling `DMSwarmCollectViewCreate()`

1425:   Collective

1427:   Input Parameters:
1428: . dm - the `DMSWARM`

1430:   Notes:
1431:   Users should call `DMSwarmCollectViewCreate()` before this function is called.

1433:   Level: advanced

1435:   Developer Note:
1436:   Create and Destroy routines create new objects that can get destroyed, they do not change the state
1437:   of the current object.

1439: .seealso: `DM`, `DMSWARM`, `DMSwarmCollectViewCreate()`, `DMSwarmSetCollectType()`
1440: @*/
1441: PetscErrorCode DMSwarmCollectViewDestroy(DM dm)
1442: {
1443:   DM_Swarm *swarm = (DM_Swarm *)dm->data;

1445:   PetscFunctionBegin;
1446:   PetscCheck(swarm->collect_view_active, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "CollectView is currently not active");
1447:   PetscCall(DMSwarmSetLocalSizes(dm, swarm->collect_view_reset_nlocal, -1));
1448:   swarm->collect_view_active = PETSC_FALSE;
1449:   PetscFunctionReturn(PETSC_SUCCESS);
1450: }

1452: static PetscErrorCode DMSwarmSetUpPIC(DM dm)
1453: {
1454:   PetscInt dim;

1456:   PetscFunctionBegin;
1457:   PetscCall(DMSwarmSetNumSpecies(dm, 1));
1458:   PetscCall(DMGetDimension(dm, &dim));
1459:   PetscCheck(dim >= 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Dimension must be 1,2,3 - found %" PetscInt_FMT, dim);
1460:   PetscCheck(dim <= 3, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Dimension must be 1,2,3 - found %" PetscInt_FMT, dim);
1461:   PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmPICField_coor, dim, PETSC_DOUBLE));
1462:   PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmPICField_cellid, 1, PETSC_INT));
1463:   PetscFunctionReturn(PETSC_SUCCESS);
1464: }

1466: /*@
1467:   DMSwarmSetPointCoordinatesRandom - Sets initial coordinates for particles in each cell

1469:   Collective

1471:   Input Parameters:
1472: + dm  - the `DMSWARM`
1473: - Npc - The number of particles per cell in the cell `DM`

1475:   Level: intermediate

1477:   Notes:
1478:   The user must use `DMSwarmSetCellDM()` to set the cell `DM` first. The particles are placed randomly inside each cell. If only
1479:   one particle is in each cell, it is placed at the centroid.

1481: .seealso: `DM`, `DMSWARM`, `DMSwarmSetCellDM()`
1482: @*/
1483: PetscErrorCode DMSwarmSetPointCoordinatesRandom(DM dm, PetscInt Npc)
1484: {
1485:   DM             cdm;
1486:   PetscRandom    rnd;
1487:   DMPolytopeType ct;
1488:   PetscBool      simplex;
1489:   PetscReal     *centroid, *coords, *xi0, *v0, *J, *invJ, detJ;
1490:   PetscInt       dim, d, cStart, cEnd, c, p;

1492:   PetscFunctionBeginUser;
1493:   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rnd));
1494:   PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0));
1495:   PetscCall(PetscRandomSetType(rnd, PETSCRAND48));

1497:   PetscCall(DMSwarmGetCellDM(dm, &cdm));
1498:   PetscCall(DMGetDimension(cdm, &dim));
1499:   PetscCall(DMPlexGetHeightStratum(cdm, 0, &cStart, &cEnd));
1500:   PetscCall(DMPlexGetCellType(cdm, cStart, &ct));
1501:   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;

1503:   PetscCall(PetscMalloc5(dim, &centroid, dim, &xi0, dim, &v0, dim * dim, &J, dim * dim, &invJ));
1504:   for (d = 0; d < dim; ++d) xi0[d] = -1.0;
1505:   PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
1506:   for (c = cStart; c < cEnd; ++c) {
1507:     if (Npc == 1) {
1508:       PetscCall(DMPlexComputeCellGeometryFVM(cdm, c, NULL, centroid, NULL));
1509:       for (d = 0; d < dim; ++d) coords[c * dim + d] = centroid[d];
1510:     } else {
1511:       PetscCall(DMPlexComputeCellGeometryFEM(cdm, c, NULL, v0, J, invJ, &detJ)); /* affine */
1512:       for (p = 0; p < Npc; ++p) {
1513:         const PetscInt n   = c * Npc + p;
1514:         PetscReal      sum = 0.0, refcoords[3];

1516:         for (d = 0; d < dim; ++d) {
1517:           PetscCall(PetscRandomGetValueReal(rnd, &refcoords[d]));
1518:           sum += refcoords[d];
1519:         }
1520:         if (simplex && sum > 0.0)
1521:           for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim) * sum;
1522:         CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n * dim]);
1523:       }
1524:     }
1525:   }
1526:   PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
1527:   PetscCall(PetscFree5(centroid, xi0, v0, J, invJ));
1528:   PetscCall(PetscRandomDestroy(&rnd));
1529:   PetscFunctionReturn(PETSC_SUCCESS);
1530: }

1532: /*@
1533:   DMSwarmSetType - Set particular flavor of `DMSWARM`

1535:   Collective

1537:   Input Parameters:
1538: + dm    - the `DMSWARM`
1539: - stype - the `DMSWARM` type (e.g. `DMSWARM_PIC`)

1541:   Level: advanced

1543: .seealso: `DM`, `DMSWARM`, `DMSwarmSetMigrateType()`, `DMSwarmSetCollectType()`, `DMSwarmType`, `DMSWARM_PIC`, `DMSWARM_BASIC`
1544: @*/
1545: PetscErrorCode DMSwarmSetType(DM dm, DMSwarmType stype)
1546: {
1547:   DM_Swarm *swarm = (DM_Swarm *)dm->data;

1549:   PetscFunctionBegin;
1550:   swarm->swarm_type = stype;
1551:   if (swarm->swarm_type == DMSWARM_PIC) PetscCall(DMSwarmSetUpPIC(dm));
1552:   PetscFunctionReturn(PETSC_SUCCESS);
1553: }

1555: static PetscErrorCode DMSetup_Swarm(DM dm)
1556: {
1557:   DM_Swarm   *swarm = (DM_Swarm *)dm->data;
1558:   PetscMPIInt rank;
1559:   PetscInt    p, npoints, *rankval;

1561:   PetscFunctionBegin;
1562:   if (swarm->issetup) PetscFunctionReturn(PETSC_SUCCESS);
1563:   swarm->issetup = PETSC_TRUE;

1565:   if (swarm->swarm_type == DMSWARM_PIC) {
1566:     /* check dmcell exists */
1567:     PetscCheck(swarm->dmcell, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSWARM_PIC requires you call DMSwarmSetCellDM");

1569:     if (swarm->dmcell->ops->locatepointssubdomain) {
1570:       /* check methods exists for exact ownership identificiation */
1571:       PetscCall(PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->ops->LocatePointsSubdomain\n"));
1572:       swarm->migrate_type = DMSWARM_MIGRATE_DMCELLEXACT;
1573:     } else {
1574:       /* check methods exist for point location AND rank neighbor identification */
1575:       if (swarm->dmcell->ops->locatepoints) {
1576:         PetscCall(PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->LocatePoints\n"));
1577:       } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSWARM_PIC requires the method CellDM->ops->locatepoints be defined");

1579:       if (swarm->dmcell->ops->getneighbors) {
1580:         PetscCall(PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->GetNeigbors\n"));
1581:       } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSWARM_PIC requires the method CellDM->ops->getneighbors be defined");

1583:       swarm->migrate_type = DMSWARM_MIGRATE_DMCELLNSCATTER;
1584:     }
1585:   }

1587:   PetscCall(DMSwarmFinalizeFieldRegister(dm));

1589:   /* check some fields were registered */
1590:   PetscCheck(swarm->db->nfields > 2, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "At least one field user must be registered via DMSwarmRegisterXXX()");

1592:   /* check local sizes were set */
1593:   PetscCheck(swarm->db->L != -1, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Local sizes must be set via DMSwarmSetLocalSizes()");

1595:   /* initialize values in pid and rank placeholders */
1596:   /* TODO: [pid - use MPI_Scan] */
1597:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
1598:   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
1599:   PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
1600:   for (p = 0; p < npoints; p++) rankval[p] = (PetscInt)rank;
1601:   PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
1602:   PetscFunctionReturn(PETSC_SUCCESS);
1603: }

1605: extern PetscErrorCode DMSwarmSortDestroy(DMSwarmSort *_ctx);

1607: static PetscErrorCode DMDestroy_Swarm(DM dm)
1608: {
1609:   DM_Swarm *swarm = (DM_Swarm *)dm->data;

1611:   PetscFunctionBegin;
1612:   if (--swarm->refct > 0) PetscFunctionReturn(PETSC_SUCCESS);
1613:   PetscCall(DMSwarmDataBucketDestroy(&swarm->db));
1614:   if (swarm->sort_context) PetscCall(DMSwarmSortDestroy(&swarm->sort_context));
1615:   PetscCall(PetscFree(swarm));
1616:   PetscFunctionReturn(PETSC_SUCCESS);
1617: }

1619: static PetscErrorCode DMSwarmView_Draw(DM dm, PetscViewer viewer)
1620: {
1621:   DM         cdm;
1622:   PetscDraw  draw;
1623:   PetscReal *coords, oldPause, radius = 0.01;
1624:   PetscInt   Np, p, bs;

1626:   PetscFunctionBegin;
1627:   PetscCall(PetscOptionsGetReal(NULL, ((PetscObject)dm)->prefix, "-dm_view_swarm_radius", &radius, NULL));
1628:   PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
1629:   PetscCall(DMSwarmGetCellDM(dm, &cdm));
1630:   PetscCall(PetscDrawGetPause(draw, &oldPause));
1631:   PetscCall(PetscDrawSetPause(draw, 0.0));
1632:   PetscCall(DMView(cdm, viewer));
1633:   PetscCall(PetscDrawSetPause(draw, oldPause));

1635:   PetscCall(DMSwarmGetLocalSize(dm, &Np));
1636:   PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, &bs, NULL, (void **)&coords));
1637:   for (p = 0; p < Np; ++p) {
1638:     const PetscInt i = p * bs;

1640:     PetscCall(PetscDrawEllipse(draw, coords[i], coords[i + 1], radius, radius, PETSC_DRAW_BLUE));
1641:   }
1642:   PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, &bs, NULL, (void **)&coords));
1643:   PetscCall(PetscDrawFlush(draw));
1644:   PetscCall(PetscDrawPause(draw));
1645:   PetscCall(PetscDrawSave(draw));
1646:   PetscFunctionReturn(PETSC_SUCCESS);
1647: }

1649: static PetscErrorCode DMView_Swarm_Ascii(DM dm, PetscViewer viewer)
1650: {
1651:   PetscViewerFormat format;
1652:   PetscInt         *sizes;
1653:   PetscInt          dim, Np, maxSize = 17;
1654:   MPI_Comm          comm;
1655:   PetscMPIInt       rank, size, p;
1656:   const char       *name;

1658:   PetscFunctionBegin;
1659:   PetscCall(PetscViewerGetFormat(viewer, &format));
1660:   PetscCall(DMGetDimension(dm, &dim));
1661:   PetscCall(DMSwarmGetLocalSize(dm, &Np));
1662:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1663:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1664:   PetscCallMPI(MPI_Comm_size(comm, &size));
1665:   PetscCall(PetscObjectGetName((PetscObject)dm, &name));
1666:   if (name) PetscCall(PetscViewerASCIIPrintf(viewer, "%s in %" PetscInt_FMT " dimension%s:\n", name, dim, dim == 1 ? "" : "s"));
1667:   else PetscCall(PetscViewerASCIIPrintf(viewer, "Swarm in %" PetscInt_FMT " dimension%s:\n", dim, dim == 1 ? "" : "s"));
1668:   if (size < maxSize) PetscCall(PetscCalloc1(size, &sizes));
1669:   else PetscCall(PetscCalloc1(3, &sizes));
1670:   if (size < maxSize) {
1671:     PetscCallMPI(MPI_Gather(&Np, 1, MPIU_INT, sizes, 1, MPIU_INT, 0, comm));
1672:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Number of particles per rank:"));
1673:     for (p = 0; p < size; ++p) {
1674:       if (rank == 0) PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, sizes[p]));
1675:     }
1676:   } else {
1677:     PetscInt locMinMax[2] = {Np, Np};

1679:     PetscCall(PetscGlobalMinMaxInt(comm, locMinMax, sizes));
1680:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Min/Max of particles per rank: %" PetscInt_FMT "/%" PetscInt_FMT, sizes[0], sizes[1]));
1681:   }
1682:   PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1683:   PetscCall(PetscFree(sizes));
1684:   if (format == PETSC_VIEWER_ASCII_INFO) {
1685:     PetscInt *cell;

1687:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Cells containing each particle:\n"));
1688:     PetscCall(PetscViewerASCIIPushSynchronized(viewer));
1689:     PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&cell));
1690:     for (p = 0; p < Np; ++p) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "  p%d: %" PetscInt_FMT "\n", p, cell[p]));
1691:     PetscCall(PetscViewerFlush(viewer));
1692:     PetscCall(PetscViewerASCIIPopSynchronized(viewer));
1693:     PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&cell));
1694:   }
1695:   PetscFunctionReturn(PETSC_SUCCESS);
1696: }

1698: static PetscErrorCode DMView_Swarm(DM dm, PetscViewer viewer)
1699: {
1700:   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1701:   PetscBool iascii, ibinary, isvtk, isdraw;
1702: #if defined(PETSC_HAVE_HDF5)
1703:   PetscBool ishdf5;
1704: #endif

1706:   PetscFunctionBegin;
1709:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1710:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &ibinary));
1711:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk));
1712: #if defined(PETSC_HAVE_HDF5)
1713:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
1714: #endif
1715:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
1716:   if (iascii) {
1717:     PetscViewerFormat format;

1719:     PetscCall(PetscViewerGetFormat(viewer, &format));
1720:     switch (format) {
1721:     case PETSC_VIEWER_ASCII_INFO_DETAIL:
1722:       PetscCall(DMSwarmDataBucketView(PetscObjectComm((PetscObject)dm), swarm->db, NULL, DATABUCKET_VIEW_STDOUT));
1723:       break;
1724:     default:
1725:       PetscCall(DMView_Swarm_Ascii(dm, viewer));
1726:     }
1727:   } else {
1728: #if defined(PETSC_HAVE_HDF5)
1729:     if (ishdf5) PetscCall(DMSwarmView_HDF5(dm, viewer));
1730: #endif
1731:     if (isdraw) PetscCall(DMSwarmView_Draw(dm, viewer));
1732:   }
1733:   PetscFunctionReturn(PETSC_SUCCESS);
1734: }

1736: /*@
1737:   DMSwarmGetCellSwarm - Extracts a single cell from the `DMSWARM` object, returns it as a single cell `DMSWARM`.
1738:   The cell `DM` is filtered for fields of that cell, and the filtered `DM` is used as the cell `DM` of the new swarm object.

1740:   Noncollective

1742:   Input Parameters:
1743: + sw        - the `DMSWARM`
1744: . cellID    - the integer id of the cell to be extracted and filtered
1745: - cellswarm - The `DMSWARM` to receive the cell

1747:   Level: beginner

1749:   Notes:
1750:   This presently only supports `DMSWARM_PIC` type

1752:   Should be restored with `DMSwarmRestoreCellSwarm()`

1754:   Changes to this cell of the swarm will be lost if they are made prior to restoring this cell.

1756: .seealso: `DM`, `DMSWARM`, `DMSwarmRestoreCellSwarm()`
1757: @*/
1758: PetscErrorCode DMSwarmGetCellSwarm(DM sw, PetscInt cellID, DM cellswarm)
1759: {
1760:   DM_Swarm *original = (DM_Swarm *)sw->data;
1761:   DMLabel   label;
1762:   DM        dmc, subdmc;
1763:   PetscInt *pids, particles, dim;

1765:   PetscFunctionBegin;
1766:   /* Configure new swarm */
1767:   PetscCall(DMSetType(cellswarm, DMSWARM));
1768:   PetscCall(DMGetDimension(sw, &dim));
1769:   PetscCall(DMSetDimension(cellswarm, dim));
1770:   PetscCall(DMSwarmSetType(cellswarm, DMSWARM_PIC));
1771:   /* Destroy the unused, unconfigured data bucket to prevent stragglers in memory */
1772:   PetscCall(DMSwarmDataBucketDestroy(&((DM_Swarm *)cellswarm->data)->db));
1773:   PetscCall(DMSwarmSortGetAccess(sw));
1774:   PetscCall(DMSwarmSortGetNumberOfPointsPerCell(sw, cellID, &particles));
1775:   PetscCall(DMSwarmSortGetPointsPerCell(sw, cellID, &particles, &pids));
1776:   PetscCall(DMSwarmDataBucketCreateFromSubset(original->db, particles, pids, &((DM_Swarm *)cellswarm->data)->db));
1777:   PetscCall(DMSwarmSortRestoreAccess(sw));
1778:   PetscCall(PetscFree(pids));
1779:   PetscCall(DMSwarmGetCellDM(sw, &dmc));
1780:   PetscCall(DMLabelCreate(PetscObjectComm((PetscObject)sw), "singlecell", &label));
1781:   PetscCall(DMAddLabel(dmc, label));
1782:   PetscCall(DMLabelSetValue(label, cellID, 1));
1783:   PetscCall(DMPlexFilter(dmc, label, 1, PETSC_FALSE, PETSC_FALSE, NULL, &subdmc));
1784:   PetscCall(DMSwarmSetCellDM(cellswarm, subdmc));
1785:   PetscCall(DMLabelDestroy(&label));
1786:   PetscFunctionReturn(PETSC_SUCCESS);
1787: }

1789: /*@
1790:   DMSwarmRestoreCellSwarm - Restores a `DMSWARM` object obtained with `DMSwarmGetCellSwarm()`. All fields are copied back into the parent swarm.

1792:   Noncollective

1794:   Input Parameters:
1795: + sw        - the parent `DMSWARM`
1796: . cellID    - the integer id of the cell to be copied back into the parent swarm
1797: - cellswarm - the cell swarm object

1799:   Level: beginner

1801:   Note:
1802:   This only supports `DMSWARM_PIC` types of `DMSWARM`s

1804: .seealso: `DM`, `DMSWARM`, `DMSwarmGetCellSwarm()`
1805: @*/
1806: PetscErrorCode DMSwarmRestoreCellSwarm(DM sw, PetscInt cellID, DM cellswarm)
1807: {
1808:   DM        dmc;
1809:   PetscInt *pids, particles, p;

1811:   PetscFunctionBegin;
1812:   PetscCall(DMSwarmSortGetAccess(sw));
1813:   PetscCall(DMSwarmSortGetPointsPerCell(sw, cellID, &particles, &pids));
1814:   PetscCall(DMSwarmSortRestoreAccess(sw));
1815:   /* Pointwise copy of each particle based on pid. The parent swarm may not be altered during this process. */
1816:   for (p = 0; p < particles; ++p) PetscCall(DMSwarmDataBucketCopyPoint(((DM_Swarm *)cellswarm->data)->db, pids[p], ((DM_Swarm *)sw->data)->db, pids[p]));
1817:   /* Free memory, destroy cell dm */
1818:   PetscCall(DMSwarmGetCellDM(cellswarm, &dmc));
1819:   PetscCall(DMDestroy(&dmc));
1820:   PetscCall(PetscFree(pids));
1821:   PetscFunctionReturn(PETSC_SUCCESS);
1822: }

1824: PETSC_INTERN PetscErrorCode DMClone_Swarm(DM, DM *);

1826: static PetscErrorCode DMInitialize_Swarm(DM sw)
1827: {
1828:   PetscFunctionBegin;
1829:   sw->dim                           = 0;
1830:   sw->ops->view                     = DMView_Swarm;
1831:   sw->ops->load                     = NULL;
1832:   sw->ops->setfromoptions           = NULL;
1833:   sw->ops->clone                    = DMClone_Swarm;
1834:   sw->ops->setup                    = DMSetup_Swarm;
1835:   sw->ops->createlocalsection       = NULL;
1836:   sw->ops->createsectionpermutation = NULL;
1837:   sw->ops->createdefaultconstraints = NULL;
1838:   sw->ops->createglobalvector       = DMCreateGlobalVector_Swarm;
1839:   sw->ops->createlocalvector        = DMCreateLocalVector_Swarm;
1840:   sw->ops->getlocaltoglobalmapping  = NULL;
1841:   sw->ops->createfieldis            = NULL;
1842:   sw->ops->createcoordinatedm       = NULL;
1843:   sw->ops->getcoloring              = NULL;
1844:   sw->ops->creatematrix             = DMCreateMatrix_Swarm;
1845:   sw->ops->createinterpolation      = NULL;
1846:   sw->ops->createinjection          = NULL;
1847:   sw->ops->createmassmatrix         = DMCreateMassMatrix_Swarm;
1848:   sw->ops->refine                   = NULL;
1849:   sw->ops->coarsen                  = NULL;
1850:   sw->ops->refinehierarchy          = NULL;
1851:   sw->ops->coarsenhierarchy         = NULL;
1852:   sw->ops->globaltolocalbegin       = NULL;
1853:   sw->ops->globaltolocalend         = NULL;
1854:   sw->ops->localtoglobalbegin       = NULL;
1855:   sw->ops->localtoglobalend         = NULL;
1856:   sw->ops->destroy                  = DMDestroy_Swarm;
1857:   sw->ops->createsubdm              = NULL;
1858:   sw->ops->getdimpoints             = NULL;
1859:   sw->ops->locatepoints             = NULL;
1860:   PetscFunctionReturn(PETSC_SUCCESS);
1861: }

1863: PETSC_INTERN PetscErrorCode DMClone_Swarm(DM dm, DM *newdm)
1864: {
1865:   DM_Swarm *swarm = (DM_Swarm *)dm->data;

1867:   PetscFunctionBegin;
1868:   swarm->refct++;
1869:   (*newdm)->data = swarm;
1870:   PetscCall(PetscObjectChangeTypeName((PetscObject)*newdm, DMSWARM));
1871:   PetscCall(DMInitialize_Swarm(*newdm));
1872:   (*newdm)->dim = dm->dim;
1873:   PetscFunctionReturn(PETSC_SUCCESS);
1874: }

1876: /*MC
1877:  DMSWARM = "swarm" - A `DM` object used to represent arrays of data (fields) of arbitrary data type.
1878:  This implementation was designed for particle methods in which the underlying
1879:  data required to be represented is both (i) dynamic in length, (ii) and of arbitrary data type.

1881:  Level: intermediate

1883:   Notes:
1884:  User data can be represented by `DMSWARM` through a registering "fields".
1885:  To register a field, the user must provide;
1886:  (a) a unique name;
1887:  (b) the data type (or size in bytes);
1888:  (c) the block size of the data.

1890:  For example, suppose the application requires a unique id, energy, momentum and density to be stored
1891:  on a set of particles. Then the following code could be used
1892: .vb
1893:     DMSwarmInitializeFieldRegister(dm)
1894:     DMSwarmRegisterPetscDatatypeField(dm,"uid",1,PETSC_LONG);
1895:     DMSwarmRegisterPetscDatatypeField(dm,"energy",1,PETSC_REAL);
1896:     DMSwarmRegisterPetscDatatypeField(dm,"momentum",3,PETSC_REAL);
1897:     DMSwarmRegisterPetscDatatypeField(dm,"density",1,PETSC_FLOAT);
1898:     DMSwarmFinalizeFieldRegister(dm)
1899: .ve

1901:  The fields represented by `DMSWARM` are dynamic and can be re-sized at any time.
1902:  The only restriction imposed by `DMSWARM` is that all fields contain the same number of points.

1904:  To support particle methods, "migration" techniques are provided. These methods migrate data
1905:  between ranks.

1907:  `DMSWARM` supports the methods `DMCreateGlobalVector()` and `DMCreateLocalVector()`.
1908:  As a `DMSWARM` may internally define and store values of different data types,
1909:  before calling `DMCreateGlobalVector()` or `DMCreateLocalVector()`, the user must inform `DMSWARM` which
1910:  fields should be used to define a `Vec` object via
1911:    `DMSwarmVectorDefineField()`
1912:  The specified field can be changed at any time - thereby permitting vectors
1913:  compatible with different fields to be created.

1915:  A dual representation of fields in the `DMSWARM` and a Vec object is permitted via
1916:    `DMSwarmCreateGlobalVectorFromField()`
1917:  Here the data defining the field in the `DMSWARM` is shared with a Vec.
1918:  This is inherently unsafe if you alter the size of the field at any time between
1919:  calls to `DMSwarmCreateGlobalVectorFromField()` and `DMSwarmDestroyGlobalVectorFromField()`.
1920:  If the local size of the `DMSWARM` does not match the local size of the global vector
1921:  when `DMSwarmDestroyGlobalVectorFromField()` is called, an error is thrown.

1923:  Additional high-level support is provided for Particle-In-Cell methods.
1924:  Please refer to `DMSwarmSetType()`.

1926: .seealso: `DM`, `DMSWARM`, `DMType`, `DMCreate()`, `DMSetType()`
1927: M*/

1929: PETSC_EXTERN PetscErrorCode DMCreate_Swarm(DM dm)
1930: {
1931:   DM_Swarm *swarm;

1933:   PetscFunctionBegin;
1935:   PetscCall(PetscNew(&swarm));
1936:   dm->data = swarm;
1937:   PetscCall(DMSwarmDataBucketCreate(&swarm->db));
1938:   PetscCall(DMSwarmInitializeFieldRegister(dm));
1939:   swarm->refct                          = 1;
1940:   swarm->vec_field_set                  = PETSC_FALSE;
1941:   swarm->issetup                        = PETSC_FALSE;
1942:   swarm->swarm_type                     = DMSWARM_BASIC;
1943:   swarm->migrate_type                   = DMSWARM_MIGRATE_BASIC;
1944:   swarm->collect_type                   = DMSWARM_COLLECT_BASIC;
1945:   swarm->migrate_error_on_missing_point = PETSC_FALSE;
1946:   swarm->dmcell                         = NULL;
1947:   swarm->collect_view_active            = PETSC_FALSE;
1948:   swarm->collect_view_reset_nlocal      = -1;
1949:   PetscCall(DMInitialize_Swarm(dm));
1950:   if (SwarmDataFieldId == -1) PetscCall(PetscObjectComposedDataRegister(&SwarmDataFieldId));
1951:   PetscFunctionReturn(PETSC_SUCCESS);
1952: }