xref: /petsc/src/dm/impls/swarm/swarm.c (revision d9acb416d05abeed0a33bde3a81aeb2ea0364f6a)
1 #define PETSCDM_DLL
2 #include <petsc/private/dmswarmimpl.h> /*I   "petscdmswarm.h"   I*/
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>
12 
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;
16 
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};
21 
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";
26 
27 PetscInt SwarmDataFieldId = -1;
28 
29 #if defined(PETSC_HAVE_HDF5)
30   #include <petscviewerhdf5.h>
31 
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;
38 
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(PetscViewerHDF5SetTimestep(viewer, seqnum));
46   /* PetscCall(DMSequenceView_HDF5(dm, "time", seqnum, (PetscScalar) seqval, viewer)); */
47   PetscCall(VecViewNative(v, viewer));
48   PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject)v, "Nc", PETSC_INT, (void *)&bs));
49   PetscCall(PetscViewerHDF5PopGroup(viewer));
50   PetscFunctionReturn(PETSC_SUCCESS);
51 }
52 
53 static PetscErrorCode DMSwarmView_HDF5(DM dm, PetscViewer viewer)
54 {
55   Vec       coordinates;
56   PetscInt  Np;
57   PetscBool isseq;
58 
59   PetscFunctionBegin;
60   PetscCall(DMSwarmGetSize(dm, &Np));
61   PetscCall(DMSwarmCreateGlobalVectorFromField(dm, DMSwarmPICField_coor, &coordinates));
62   PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
63   PetscCall(PetscViewerHDF5PushGroup(viewer, "/particles"));
64   PetscCall(PetscObjectTypeCompare((PetscObject)coordinates, VECSEQ, &isseq));
65   PetscCall(VecViewNative(coordinates, viewer));
66   PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject)coordinates, "Np", PETSC_INT, (void *)&Np));
67   PetscCall(PetscViewerHDF5PopGroup(viewer));
68   PetscCall(DMSwarmDestroyGlobalVectorFromField(dm, DMSwarmPICField_coor, &coordinates));
69   PetscFunctionReturn(PETSC_SUCCESS);
70 }
71 #endif
72 
73 static PetscErrorCode VecView_Swarm(Vec v, PetscViewer viewer)
74 {
75   DM dm;
76 #if defined(PETSC_HAVE_HDF5)
77   PetscBool ishdf5;
78 #endif
79 
80   PetscFunctionBegin;
81   PetscCall(VecGetDM(v, &dm));
82   PetscCheck(dm, PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
83 #if defined(PETSC_HAVE_HDF5)
84   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
85   if (ishdf5) {
86     PetscCall(VecView_Swarm_HDF5_Internal(v, viewer));
87     PetscFunctionReturn(PETSC_SUCCESS);
88   }
89 #endif
90   PetscCall(VecViewNative(v, viewer));
91   PetscFunctionReturn(PETSC_SUCCESS);
92 }
93 
94 /*@C
95   DMSwarmVectorDefineField - Sets the field from which to define a `Vec` object
96   when `DMCreateLocalVector()`, or `DMCreateGlobalVector()` is called
97 
98   Collective
99 
100   Input Parameters:
101 + dm        - a `DMSWARM`
102 - fieldname - the textual name given to a registered field
103 
104   Level: beginner
105 
106   Notes:
107   The field with name `fieldname` must be defined as having a data type of `PetscScalar`.
108 
109   This function must be called prior to calling `DMCreateLocalVector()`, `DMCreateGlobalVector()`.
110   Multiple calls to `DMSwarmVectorDefineField()` are permitted.
111 
112 .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`
113 @*/
114 PetscErrorCode DMSwarmVectorDefineField(DM dm, const char fieldname[])
115 {
116   DM_Swarm     *swarm = (DM_Swarm *)dm->data;
117   PetscInt      bs, n;
118   PetscScalar  *array;
119   PetscDataType type;
120 
121   PetscFunctionBegin;
122   if (!swarm->issetup) PetscCall(DMSetUp(dm));
123   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &n, NULL, NULL));
124   PetscCall(DMSwarmGetField(dm, fieldname, &bs, &type, (void **)&array));
125 
126   /* Check all fields are of type PETSC_REAL or PETSC_SCALAR */
127   PetscCheck(type == PETSC_REAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only valid for PETSC_REAL");
128   PetscCall(PetscSNPrintf(swarm->vec_field_name, PETSC_MAX_PATH_LEN - 1, "%s", fieldname));
129   swarm->vec_field_set    = PETSC_TRUE;
130   swarm->vec_field_bs     = bs;
131   swarm->vec_field_nlocal = n;
132   PetscCall(DMSwarmRestoreField(dm, fieldname, &bs, &type, (void **)&array));
133   PetscFunctionReturn(PETSC_SUCCESS);
134 }
135 
136 /* requires DMSwarmDefineFieldVector has been called */
137 static PetscErrorCode DMCreateGlobalVector_Swarm(DM dm, Vec *vec)
138 {
139   DM_Swarm *swarm = (DM_Swarm *)dm->data;
140   Vec       x;
141   char      name[PETSC_MAX_PATH_LEN];
142 
143   PetscFunctionBegin;
144   if (!swarm->issetup) PetscCall(DMSetUp(dm));
145   PetscCheck(swarm->vec_field_set, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMSwarmVectorDefineField first");
146   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 */
147 
148   PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "DMSwarmField_%s", swarm->vec_field_name));
149   PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), &x));
150   PetscCall(PetscObjectSetName((PetscObject)x, name));
151   PetscCall(VecSetSizes(x, swarm->db->L * swarm->vec_field_bs, PETSC_DETERMINE));
152   PetscCall(VecSetBlockSize(x, swarm->vec_field_bs));
153   PetscCall(VecSetDM(x, dm));
154   PetscCall(VecSetFromOptions(x));
155   PetscCall(VecSetDM(x, dm));
156   PetscCall(VecSetOperation(x, VECOP_VIEW, (void (*)(void))VecView_Swarm));
157   *vec = x;
158   PetscFunctionReturn(PETSC_SUCCESS);
159 }
160 
161 /* requires DMSwarmDefineFieldVector has been called */
162 static PetscErrorCode DMCreateLocalVector_Swarm(DM dm, Vec *vec)
163 {
164   DM_Swarm *swarm = (DM_Swarm *)dm->data;
165   Vec       x;
166   char      name[PETSC_MAX_PATH_LEN];
167 
168   PetscFunctionBegin;
169   if (!swarm->issetup) PetscCall(DMSetUp(dm));
170   PetscCheck(swarm->vec_field_set, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMSwarmVectorDefineField first");
171   PetscCheck(swarm->vec_field_nlocal == swarm->db->L, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSwarm sizes have changed since last call to VectorDefineField first");
172 
173   PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "DMSwarmField_%s", swarm->vec_field_name));
174   PetscCall(VecCreate(PETSC_COMM_SELF, &x));
175   PetscCall(PetscObjectSetName((PetscObject)x, name));
176   PetscCall(VecSetSizes(x, swarm->db->L * swarm->vec_field_bs, PETSC_DETERMINE));
177   PetscCall(VecSetBlockSize(x, swarm->vec_field_bs));
178   PetscCall(VecSetDM(x, dm));
179   PetscCall(VecSetFromOptions(x));
180   *vec = x;
181   PetscFunctionReturn(PETSC_SUCCESS);
182 }
183 
184 static PetscErrorCode DMSwarmDestroyVectorFromField_Private(DM dm, const char fieldname[], Vec *vec)
185 {
186   DM_Swarm        *swarm = (DM_Swarm *)dm->data;
187   DMSwarmDataField gfield;
188   PetscInt         bs, nlocal, fid = -1, cfid = -2;
189   PetscBool        flg;
190 
191   PetscFunctionBegin;
192   /* check vector is an inplace array */
193   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldIdByName(swarm->db, fieldname, &fid));
194   PetscCall(PetscObjectComposedDataGetInt((PetscObject)*vec, SwarmDataFieldId, cfid, flg));
195   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);
196   PetscCall(VecGetLocalSize(*vec, &nlocal));
197   PetscCall(VecGetBlockSize(*vec, &bs));
198   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");
199   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
200   PetscCall(DMSwarmDataFieldRestoreAccess(gfield));
201   PetscCall(VecResetArray(*vec));
202   PetscCall(VecDestroy(vec));
203   PetscFunctionReturn(PETSC_SUCCESS);
204 }
205 
206 static PetscErrorCode DMSwarmCreateVectorFromField_Private(DM dm, const char fieldname[], MPI_Comm comm, Vec *vec)
207 {
208   DM_Swarm     *swarm = (DM_Swarm *)dm->data;
209   PetscDataType type;
210   PetscScalar  *array;
211   PetscInt      bs, n, fid;
212   char          name[PETSC_MAX_PATH_LEN];
213   PetscMPIInt   size;
214   PetscBool     iscuda, iskokkos;
215 
216   PetscFunctionBegin;
217   if (!swarm->issetup) PetscCall(DMSetUp(dm));
218   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &n, NULL, NULL));
219   PetscCall(DMSwarmGetField(dm, fieldname, &bs, &type, (void **)&array));
220   PetscCheck(type == PETSC_REAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only valid for PETSC_REAL");
221 
222   PetscCallMPI(MPI_Comm_size(comm, &size));
223   PetscCall(PetscStrcmp(dm->vectype, VECKOKKOS, &iskokkos));
224   PetscCall(PetscStrcmp(dm->vectype, VECCUDA, &iscuda));
225   PetscCall(VecCreate(comm, vec));
226   PetscCall(VecSetSizes(*vec, n * bs, PETSC_DETERMINE));
227   PetscCall(VecSetBlockSize(*vec, bs));
228   if (iskokkos) PetscCall(VecSetType(*vec, VECKOKKOS));
229   else if (iscuda) PetscCall(VecSetType(*vec, VECCUDA));
230   else PetscCall(VecSetType(*vec, VECSTANDARD));
231   PetscCall(VecPlaceArray(*vec, array));
232 
233   PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "DMSwarmSharedField_%s", fieldname));
234   PetscCall(PetscObjectSetName((PetscObject)*vec, name));
235 
236   /* Set guard */
237   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldIdByName(swarm->db, fieldname, &fid));
238   PetscCall(PetscObjectComposedDataSetInt((PetscObject)*vec, SwarmDataFieldId, fid));
239 
240   PetscCall(VecSetDM(*vec, dm));
241   PetscCall(VecSetOperation(*vec, VECOP_VIEW, (void (*)(void))VecView_Swarm));
242   PetscFunctionReturn(PETSC_SUCCESS);
243 }
244 
245 /* This creates a "mass matrix" between a finite element and particle space. If a finite element interpolant is given by
246 
247      \hat f = \sum_i f_i \phi_i
248 
249    and a particle function is given by
250 
251      f = \sum_i w_i \delta(x - x_i)
252 
253    then we want to require that
254 
255      M \hat f = M_p f
256 
257    where the particle mass matrix is given by
258 
259      (M_p)_{ij} = \int \phi_i \delta(x - x_j)
260 
261    The way Dave May does particles, they amount to quadratue weights rather than delta functions, so he has |J| is in
262    his integral. We allow this with the boolean flag.
263 */
264 static PetscErrorCode DMSwarmComputeMassMatrix_Private(DM dmc, DM dmf, Mat mass, PetscBool useDeltaFunction, void *ctx)
265 {
266   const char  *name = "Mass Matrix";
267   MPI_Comm     comm;
268   PetscDS      prob;
269   PetscSection fsection, globalFSection;
270   PetscHSetIJ  ht;
271   PetscLayout  rLayout, colLayout;
272   PetscInt    *dnz, *onz;
273   PetscInt     locRows, locCols, rStart, colStart, colEnd, *rowIDXs;
274   PetscReal   *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0};
275   PetscScalar *elemMat;
276   PetscInt     dim, Nf, field, cStart, cEnd, cell, totDim, maxC = 0;
277 
278   PetscFunctionBegin;
279   PetscCall(PetscObjectGetComm((PetscObject)mass, &comm));
280   PetscCall(DMGetCoordinateDim(dmf, &dim));
281   PetscCall(DMGetDS(dmf, &prob));
282   PetscCall(PetscDSGetNumFields(prob, &Nf));
283   PetscCall(PetscDSGetTotalDimension(prob, &totDim));
284   PetscCall(PetscMalloc3(dim, &v0, dim * dim, &J, dim * dim, &invJ));
285   PetscCall(DMGetLocalSection(dmf, &fsection));
286   PetscCall(DMGetGlobalSection(dmf, &globalFSection));
287   PetscCall(DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd));
288   PetscCall(MatGetLocalSize(mass, &locRows, &locCols));
289 
290   PetscCall(PetscLayoutCreate(comm, &colLayout));
291   PetscCall(PetscLayoutSetLocalSize(colLayout, locCols));
292   PetscCall(PetscLayoutSetBlockSize(colLayout, 1));
293   PetscCall(PetscLayoutSetUp(colLayout));
294   PetscCall(PetscLayoutGetRange(colLayout, &colStart, &colEnd));
295   PetscCall(PetscLayoutDestroy(&colLayout));
296 
297   PetscCall(PetscLayoutCreate(comm, &rLayout));
298   PetscCall(PetscLayoutSetLocalSize(rLayout, locRows));
299   PetscCall(PetscLayoutSetBlockSize(rLayout, 1));
300   PetscCall(PetscLayoutSetUp(rLayout));
301   PetscCall(PetscLayoutGetRange(rLayout, &rStart, NULL));
302   PetscCall(PetscLayoutDestroy(&rLayout));
303 
304   PetscCall(PetscCalloc2(locRows, &dnz, locRows, &onz));
305   PetscCall(PetscHSetIJCreate(&ht));
306 
307   PetscCall(PetscSynchronizedFlush(comm, NULL));
308   /* count non-zeros */
309   PetscCall(DMSwarmSortGetAccess(dmc));
310   for (field = 0; field < Nf; ++field) {
311     for (cell = cStart; cell < cEnd; ++cell) {
312       PetscInt  c, i;
313       PetscInt *findices, *cindices; /* fine is vertices, coarse is particles */
314       PetscInt  numFIndices, numCIndices;
315 
316       PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
317       PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
318       maxC = PetscMax(maxC, numCIndices);
319       {
320         PetscHashIJKey key;
321         PetscBool      missing;
322         for (i = 0; i < numFIndices; ++i) {
323           key.j = findices[i]; /* global column (from Plex) */
324           if (key.j >= 0) {
325             /* Get indices for coarse elements */
326             for (c = 0; c < numCIndices; ++c) {
327               key.i = cindices[c] + rStart; /* global cols (from Swarm) */
328               if (key.i < 0) continue;
329               PetscCall(PetscHSetIJQueryAdd(ht, key, &missing));
330               if (missing) {
331                 if ((key.j >= colStart) && (key.j < colEnd)) ++dnz[key.i - rStart];
332                 else ++onz[key.i - rStart];
333               } else SETERRQ(PetscObjectComm((PetscObject)dmf), PETSC_ERR_SUP, "Set new value at %" PetscInt_FMT ",%" PetscInt_FMT, key.i, key.j);
334             }
335           }
336         }
337         PetscCall(PetscFree(cindices));
338       }
339       PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
340     }
341   }
342   PetscCall(PetscHSetIJDestroy(&ht));
343   PetscCall(MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL));
344   PetscCall(MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
345   PetscCall(PetscFree2(dnz, onz));
346   PetscCall(PetscMalloc3(maxC * totDim, &elemMat, maxC, &rowIDXs, maxC * dim, &xi));
347   for (field = 0; field < Nf; ++field) {
348     PetscTabulation Tcoarse;
349     PetscObject     obj;
350     PetscClassId    id;
351     PetscReal      *fieldVals;
352     PetscInt        Nc, i;
353 
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));
358     PetscCheck(Nc == 1, PetscObjectComm((PetscObject)dmf), PETSC_ERR_SUP, "Can only interpolate a scalar field from particles, Nc = %" PetscInt_FMT, Nc);
359 
360     PetscCall(DMSwarmGetField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **)&fieldVals));
361     for (cell = cStart; cell < cEnd; ++cell) {
362       PetscInt *findices, *cindices;
363       PetscInt  numFIndices, numCIndices;
364       PetscInt  p, c;
365 
366       /* TODO: Use DMField instead of assuming affine */
367       PetscCall(DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ));
368       PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
369       PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
370       for (p = 0; p < numCIndices; ++p) CoordinatesRealToRef(dim, dim, v0ref, v0, invJ, &fieldVals[cindices[p] * dim], &xi[p * dim]);
371       if (id == PETSCFE_CLASSID) PetscCall(PetscFECreateTabulation((PetscFE)obj, 1, numCIndices, xi, 0, &Tcoarse));
372       else PetscCall(PetscFVCreateTabulation((PetscFV)obj, 1, numCIndices, xi, 0, &Tcoarse));
373       /* Get elemMat entries by multiplying by weight */
374       PetscCall(PetscArrayzero(elemMat, numCIndices * totDim));
375       for (i = 0; i < numFIndices; ++i) {
376         for (p = 0; p < numCIndices; ++p) {
377           for (c = 0; c < Nc; ++c) {
378             /* B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c */
379             elemMat[p * numFIndices + i] += Tcoarse->T[0][(p * numFIndices + i) * Nc + c] * (useDeltaFunction ? 1.0 : detJ);
380           }
381         }
382       }
383       for (p = 0; p < numCIndices; ++p) rowIDXs[p] = cindices[p] + rStart;
384       if (0) PetscCall(DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat));
385       PetscCall(MatSetValues(mass, numCIndices, rowIDXs, numFIndices, findices, elemMat, ADD_VALUES));
386       PetscCall(PetscFree(cindices));
387       PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
388       PetscCall(PetscTabulationDestroy(&Tcoarse));
389     }
390     PetscCall(DMSwarmRestoreField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **)&fieldVals));
391   }
392   PetscCall(PetscFree3(elemMat, rowIDXs, xi));
393   PetscCall(DMSwarmSortRestoreAccess(dmc));
394   PetscCall(PetscFree3(v0, J, invJ));
395   PetscCall(MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY));
396   PetscCall(MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY));
397   PetscFunctionReturn(PETSC_SUCCESS);
398 }
399 
400 /* Returns empty matrix for use with SNES FD */
401 static PetscErrorCode DMCreateMatrix_Swarm(DM sw, Mat *m)
402 {
403   Vec      field;
404   PetscInt size;
405 
406   PetscFunctionBegin;
407   PetscCall(DMGetGlobalVector(sw, &field));
408   PetscCall(VecGetLocalSize(field, &size));
409   PetscCall(DMRestoreGlobalVector(sw, &field));
410   PetscCall(MatCreate(PETSC_COMM_WORLD, m));
411   PetscCall(MatSetFromOptions(*m));
412   PetscCall(MatSetSizes(*m, PETSC_DECIDE, PETSC_DECIDE, size, size));
413   PetscCall(MatSeqAIJSetPreallocation(*m, 1, NULL));
414   PetscCall(MatZeroEntries(*m));
415   PetscCall(MatAssemblyBegin(*m, MAT_FINAL_ASSEMBLY));
416   PetscCall(MatAssemblyEnd(*m, MAT_FINAL_ASSEMBLY));
417   PetscCall(MatShift(*m, 1.0));
418   PetscCall(MatSetDM(*m, sw));
419   PetscFunctionReturn(PETSC_SUCCESS);
420 }
421 
422 /* FEM cols, Particle rows */
423 static PetscErrorCode DMCreateMassMatrix_Swarm(DM dmCoarse, DM dmFine, Mat *mass)
424 {
425   PetscSection gsf;
426   PetscInt     m, n;
427   void        *ctx;
428 
429   PetscFunctionBegin;
430   PetscCall(DMGetGlobalSection(dmFine, &gsf));
431   PetscCall(PetscSectionGetConstrainedStorageSize(gsf, &m));
432   PetscCall(DMSwarmGetLocalSize(dmCoarse, &n));
433   PetscCall(MatCreate(PetscObjectComm((PetscObject)dmCoarse), mass));
434   PetscCall(MatSetSizes(*mass, n, m, PETSC_DETERMINE, PETSC_DETERMINE));
435   PetscCall(MatSetType(*mass, dmCoarse->mattype));
436   PetscCall(DMGetApplicationContext(dmFine, &ctx));
437 
438   PetscCall(DMSwarmComputeMassMatrix_Private(dmCoarse, dmFine, *mass, PETSC_TRUE, ctx));
439   PetscCall(MatViewFromOptions(*mass, NULL, "-mass_mat_view"));
440   PetscFunctionReturn(PETSC_SUCCESS);
441 }
442 
443 static PetscErrorCode DMSwarmComputeMassMatrixSquare_Private(DM dmc, DM dmf, Mat mass, PetscBool useDeltaFunction, void *ctx)
444 {
445   const char  *name = "Mass Matrix Square";
446   MPI_Comm     comm;
447   PetscDS      prob;
448   PetscSection fsection, globalFSection;
449   PetscHSetIJ  ht;
450   PetscLayout  rLayout, colLayout;
451   PetscInt    *dnz, *onz, *adj, depth, maxConeSize, maxSupportSize, maxAdjSize;
452   PetscInt     locRows, locCols, rStart, colStart, colEnd, *rowIDXs;
453   PetscReal   *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0};
454   PetscScalar *elemMat, *elemMatSq;
455   PetscInt     cdim, Nf, field, cStart, cEnd, cell, totDim, maxC = 0;
456 
457   PetscFunctionBegin;
458   PetscCall(PetscObjectGetComm((PetscObject)mass, &comm));
459   PetscCall(DMGetCoordinateDim(dmf, &cdim));
460   PetscCall(DMGetDS(dmf, &prob));
461   PetscCall(PetscDSGetNumFields(prob, &Nf));
462   PetscCall(PetscDSGetTotalDimension(prob, &totDim));
463   PetscCall(PetscMalloc3(cdim, &v0, cdim * cdim, &J, cdim * cdim, &invJ));
464   PetscCall(DMGetLocalSection(dmf, &fsection));
465   PetscCall(DMGetGlobalSection(dmf, &globalFSection));
466   PetscCall(DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd));
467   PetscCall(MatGetLocalSize(mass, &locRows, &locCols));
468 
469   PetscCall(PetscLayoutCreate(comm, &colLayout));
470   PetscCall(PetscLayoutSetLocalSize(colLayout, locCols));
471   PetscCall(PetscLayoutSetBlockSize(colLayout, 1));
472   PetscCall(PetscLayoutSetUp(colLayout));
473   PetscCall(PetscLayoutGetRange(colLayout, &colStart, &colEnd));
474   PetscCall(PetscLayoutDestroy(&colLayout));
475 
476   PetscCall(PetscLayoutCreate(comm, &rLayout));
477   PetscCall(PetscLayoutSetLocalSize(rLayout, locRows));
478   PetscCall(PetscLayoutSetBlockSize(rLayout, 1));
479   PetscCall(PetscLayoutSetUp(rLayout));
480   PetscCall(PetscLayoutGetRange(rLayout, &rStart, NULL));
481   PetscCall(PetscLayoutDestroy(&rLayout));
482 
483   PetscCall(DMPlexGetDepth(dmf, &depth));
484   PetscCall(DMPlexGetMaxSizes(dmf, &maxConeSize, &maxSupportSize));
485   maxAdjSize = PetscPowInt(maxConeSize * maxSupportSize, depth);
486   PetscCall(PetscMalloc1(maxAdjSize, &adj));
487 
488   PetscCall(PetscCalloc2(locRows, &dnz, locRows, &onz));
489   PetscCall(PetscHSetIJCreate(&ht));
490   /* Count nonzeros
491        This is just FVM++, but we cannot use the Plex P0 allocation since unknowns in a cell will not be contiguous
492   */
493   PetscCall(DMSwarmSortGetAccess(dmc));
494   for (cell = cStart; cell < cEnd; ++cell) {
495     PetscInt  i;
496     PetscInt *cindices;
497     PetscInt  numCIndices;
498 #if 0
499     PetscInt  adjSize = maxAdjSize, a, j;
500 #endif
501 
502     PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
503     maxC = PetscMax(maxC, numCIndices);
504     /* Diagonal block */
505     for (i = 0; i < numCIndices; ++i) dnz[cindices[i]] += numCIndices;
506 #if 0
507     /* Off-diagonal blocks */
508     PetscCall(DMPlexGetAdjacency(dmf, cell, &adjSize, &adj));
509     for (a = 0; a < adjSize; ++a) {
510       if (adj[a] >= cStart && adj[a] < cEnd && adj[a] != cell) {
511         const PetscInt ncell = adj[a];
512         PetscInt      *ncindices;
513         PetscInt       numNCIndices;
514 
515         PetscCall(DMSwarmSortGetPointsPerCell(dmc, ncell, &numNCIndices, &ncindices));
516         {
517           PetscHashIJKey key;
518           PetscBool      missing;
519 
520           for (i = 0; i < numCIndices; ++i) {
521             key.i = cindices[i] + rStart; /* global rows (from Swarm) */
522             if (key.i < 0) continue;
523             for (j = 0; j < numNCIndices; ++j) {
524               key.j = ncindices[j] + rStart; /* global column (from Swarm) */
525               if (key.j < 0) continue;
526               PetscCall(PetscHSetIJQueryAdd(ht, key, &missing));
527               if (missing) {
528                 if ((key.j >= colStart) && (key.j < colEnd)) ++dnz[key.i - rStart];
529                 else                                         ++onz[key.i - rStart];
530               }
531             }
532           }
533         }
534         PetscCall(PetscFree(ncindices));
535       }
536     }
537 #endif
538     PetscCall(PetscFree(cindices));
539   }
540   PetscCall(PetscHSetIJDestroy(&ht));
541   PetscCall(MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL));
542   PetscCall(MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
543   PetscCall(PetscFree2(dnz, onz));
544   PetscCall(PetscMalloc4(maxC * totDim, &elemMat, maxC * maxC, &elemMatSq, maxC, &rowIDXs, maxC * cdim, &xi));
545   /* Fill in values
546        Each entry is a sum of terms \phi_i(x_p) \phi_i(x_q)
547        Start just by producing block diagonal
548        Could loop over adjacent cells
549          Produce neighboring element matrix
550          TODO Determine which columns and rows correspond to shared dual vector
551          Do MatMatMult with rectangular matrices
552          Insert block
553   */
554   for (field = 0; field < Nf; ++field) {
555     PetscTabulation Tcoarse;
556     PetscObject     obj;
557     PetscReal      *coords;
558     PetscInt        Nc, i;
559 
560     PetscCall(PetscDSGetDiscretization(prob, field, &obj));
561     PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
562     PetscCheck(Nc == 1, PetscObjectComm((PetscObject)dmf), PETSC_ERR_SUP, "Can only interpolate a scalar field from particles, Nc = %" PetscInt_FMT, Nc);
563     PetscCall(DMSwarmGetField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
564     for (cell = cStart; cell < cEnd; ++cell) {
565       PetscInt *findices, *cindices;
566       PetscInt  numFIndices, numCIndices;
567       PetscInt  p, c;
568 
569       /* TODO: Use DMField instead of assuming affine */
570       PetscCall(DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ));
571       PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
572       PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
573       for (p = 0; p < numCIndices; ++p) CoordinatesRealToRef(cdim, cdim, v0ref, v0, invJ, &coords[cindices[p] * cdim], &xi[p * cdim]);
574       PetscCall(PetscFECreateTabulation((PetscFE)obj, 1, numCIndices, xi, 0, &Tcoarse));
575       /* Get elemMat entries by multiplying by weight */
576       PetscCall(PetscArrayzero(elemMat, numCIndices * totDim));
577       for (i = 0; i < numFIndices; ++i) {
578         for (p = 0; p < numCIndices; ++p) {
579           for (c = 0; c < Nc; ++c) {
580             /* B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c */
581             elemMat[p * numFIndices + i] += Tcoarse->T[0][(p * numFIndices + i) * Nc + c] * (useDeltaFunction ? 1.0 : detJ);
582           }
583         }
584       }
585       PetscCall(PetscTabulationDestroy(&Tcoarse));
586       for (p = 0; p < numCIndices; ++p) rowIDXs[p] = cindices[p] + rStart;
587       if (0) PetscCall(DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat));
588       /* Block diagonal */
589       if (numCIndices) {
590         PetscBLASInt blasn, blask;
591         PetscScalar  one = 1.0, zero = 0.0;
592 
593         PetscCall(PetscBLASIntCast(numCIndices, &blasn));
594         PetscCall(PetscBLASIntCast(numFIndices, &blask));
595         PetscCallBLAS("BLASgemm", BLASgemm_("T", "N", &blasn, &blasn, &blask, &one, elemMat, &blask, elemMat, &blask, &zero, elemMatSq, &blasn));
596       }
597       PetscCall(MatSetValues(mass, numCIndices, rowIDXs, numCIndices, rowIDXs, elemMatSq, ADD_VALUES));
598       /* TODO off-diagonal */
599       PetscCall(PetscFree(cindices));
600       PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
601     }
602     PetscCall(DMSwarmRestoreField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
603   }
604   PetscCall(PetscFree4(elemMat, elemMatSq, rowIDXs, xi));
605   PetscCall(PetscFree(adj));
606   PetscCall(DMSwarmSortRestoreAccess(dmc));
607   PetscCall(PetscFree3(v0, J, invJ));
608   PetscCall(MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY));
609   PetscCall(MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY));
610   PetscFunctionReturn(PETSC_SUCCESS);
611 }
612 
613 /*@C
614   DMSwarmCreateMassMatrixSquare - Creates the block-diagonal of the square, M^T_p M_p, of the particle mass matrix M_p
615 
616   Collective
617 
618   Input Parameters:
619 + dmCoarse - a `DMSWARM`
620 - dmFine   - a `DMPLEX`
621 
622   Output Parameter:
623 . mass - the square of the particle mass matrix
624 
625   Level: advanced
626 
627   Note:
628   We only compute the block diagonal since this provides a good preconditioner and is completely local. It would be possible in the
629   future to compute the full normal equations.
630 
631 .seealso: `DM`, `DMSWARM`, `DMCreateMassMatrix()`
632 @*/
633 PetscErrorCode DMSwarmCreateMassMatrixSquare(DM dmCoarse, DM dmFine, Mat *mass)
634 {
635   PetscInt n;
636   void    *ctx;
637 
638   PetscFunctionBegin;
639   PetscCall(DMSwarmGetLocalSize(dmCoarse, &n));
640   PetscCall(MatCreate(PetscObjectComm((PetscObject)dmCoarse), mass));
641   PetscCall(MatSetSizes(*mass, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
642   PetscCall(MatSetType(*mass, dmCoarse->mattype));
643   PetscCall(DMGetApplicationContext(dmFine, &ctx));
644 
645   PetscCall(DMSwarmComputeMassMatrixSquare_Private(dmCoarse, dmFine, *mass, PETSC_TRUE, ctx));
646   PetscCall(MatViewFromOptions(*mass, NULL, "-mass_sq_mat_view"));
647   PetscFunctionReturn(PETSC_SUCCESS);
648 }
649 
650 /*@C
651   DMSwarmCreateGlobalVectorFromField - Creates a `Vec` object sharing the array associated with a given field
652 
653   Collective
654 
655   Input Parameters:
656 + dm        - a `DMSWARM`
657 - fieldname - the textual name given to a registered field
658 
659   Output Parameter:
660 . vec - the vector
661 
662   Level: beginner
663 
664   Notes:
665   The vector must be returned using a matching call to `DMSwarmDestroyGlobalVectorFromField()`.
666 
667 .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmDestroyGlobalVectorFromField()`
668 @*/
669 PetscErrorCode DMSwarmCreateGlobalVectorFromField(DM dm, const char fieldname[], Vec *vec)
670 {
671   MPI_Comm comm = PetscObjectComm((PetscObject)dm);
672 
673   PetscFunctionBegin;
674   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
675   PetscCall(DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec));
676   PetscFunctionReturn(PETSC_SUCCESS);
677 }
678 
679 /*@C
680   DMSwarmDestroyGlobalVectorFromField - Destroys the `Vec` object which share the array associated with a given field
681 
682   Collective
683 
684   Input Parameters:
685 + dm        - a `DMSWARM`
686 - fieldname - the textual name given to a registered field
687 
688   Output Parameter:
689 . vec - the vector
690 
691   Level: beginner
692 
693 .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmCreateGlobalVectorFromField()`
694 @*/
695 PetscErrorCode DMSwarmDestroyGlobalVectorFromField(DM dm, const char fieldname[], Vec *vec)
696 {
697   PetscFunctionBegin;
698   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
699   PetscCall(DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec));
700   PetscFunctionReturn(PETSC_SUCCESS);
701 }
702 
703 /*@C
704   DMSwarmCreateLocalVectorFromField - Creates a `Vec` object sharing the array associated with a given field
705 
706   Collective
707 
708   Input Parameters:
709 + dm        - a `DMSWARM`
710 - fieldname - the textual name given to a registered field
711 
712   Output Parameter:
713 . vec - the vector
714 
715   Level: beginner
716 
717   Note:
718   The vector must be returned using a matching call to DMSwarmDestroyLocalVectorFromField().
719 
720 .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmDestroyLocalVectorFromField()`
721 @*/
722 PetscErrorCode DMSwarmCreateLocalVectorFromField(DM dm, const char fieldname[], Vec *vec)
723 {
724   MPI_Comm comm = PETSC_COMM_SELF;
725 
726   PetscFunctionBegin;
727   PetscCall(DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec));
728   PetscFunctionReturn(PETSC_SUCCESS);
729 }
730 
731 /*@C
732   DMSwarmDestroyLocalVectorFromField - Destroys the `Vec` object which share the array associated with a given field
733 
734   Collective
735 
736   Input Parameters:
737 + dm        - a `DMSWARM`
738 - fieldname - the textual name given to a registered field
739 
740   Output Parameter:
741 . vec - the vector
742 
743   Level: beginner
744 
745 .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmCreateLocalVectorFromField()`
746 @*/
747 PetscErrorCode DMSwarmDestroyLocalVectorFromField(DM dm, const char fieldname[], Vec *vec)
748 {
749   PetscFunctionBegin;
750   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
751   PetscCall(DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec));
752   PetscFunctionReturn(PETSC_SUCCESS);
753 }
754 
755 /*@C
756   DMSwarmInitializeFieldRegister - Initiates the registration of fields to a `DMSWARM`
757 
758   Collective
759 
760   Input Parameter:
761 . dm - a `DMSWARM`
762 
763   Level: beginner
764 
765   Note:
766   After all fields have been registered, you must call `DMSwarmFinalizeFieldRegister()`.
767 
768 .seealso: `DM`, `DMSWARM`, `DMSwarmFinalizeFieldRegister()`, `DMSwarmRegisterPetscDatatypeField()`,
769           `DMSwarmRegisterUserStructField()`, `DMSwarmRegisterUserDatatypeField()`
770 @*/
771 PetscErrorCode DMSwarmInitializeFieldRegister(DM dm)
772 {
773   DM_Swarm *swarm = (DM_Swarm *)dm->data;
774 
775   PetscFunctionBegin;
776   if (!swarm->field_registration_initialized) {
777     swarm->field_registration_initialized = PETSC_TRUE;
778     PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmField_pid, 1, PETSC_INT64)); /* unique identifier */
779     PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmField_rank, 1, PETSC_INT));  /* used for communication */
780   }
781   PetscFunctionReturn(PETSC_SUCCESS);
782 }
783 
784 /*@
785   DMSwarmFinalizeFieldRegister - Finalizes the registration of fields to a `DMSWARM`
786 
787   Collective
788 
789   Input Parameter:
790 . dm - a `DMSWARM`
791 
792   Level: beginner
793 
794   Note:
795   After `DMSwarmFinalizeFieldRegister()` has been called, no new fields can be defined on the `DMSWARM`.
796 
797 .seealso: `DM`, `DMSWARM`, `DMSwarmInitializeFieldRegister()`, `DMSwarmRegisterPetscDatatypeField()`,
798           `DMSwarmRegisterUserStructField()`, `DMSwarmRegisterUserDatatypeField()`
799 @*/
800 PetscErrorCode DMSwarmFinalizeFieldRegister(DM dm)
801 {
802   DM_Swarm *swarm = (DM_Swarm *)dm->data;
803 
804   PetscFunctionBegin;
805   if (!swarm->field_registration_finalized) PetscCall(DMSwarmDataBucketFinalize(swarm->db));
806   swarm->field_registration_finalized = PETSC_TRUE;
807   PetscFunctionReturn(PETSC_SUCCESS);
808 }
809 
810 /*@
811   DMSwarmSetLocalSizes - Sets the length of all registered fields on the `DMSWARM`
812 
813   Not Collective
814 
815   Input Parameters:
816 + dm     - a `DMSWARM`
817 . nlocal - the length of each registered field
818 - buffer - the length of the buffer used to efficient dynamic re-sizing
819 
820   Level: beginner
821 
822 .seealso: `DM`, `DMSWARM`, `DMSwarmGetLocalSize()`
823 @*/
824 PetscErrorCode DMSwarmSetLocalSizes(DM dm, PetscInt nlocal, PetscInt buffer)
825 {
826   DM_Swarm *swarm = (DM_Swarm *)dm->data;
827 
828   PetscFunctionBegin;
829   PetscCall(PetscLogEventBegin(DMSWARM_SetSizes, 0, 0, 0, 0));
830   PetscCall(DMSwarmDataBucketSetSizes(swarm->db, nlocal, buffer));
831   PetscCall(PetscLogEventEnd(DMSWARM_SetSizes, 0, 0, 0, 0));
832   PetscFunctionReturn(PETSC_SUCCESS);
833 }
834 
835 /*@
836   DMSwarmSetCellDM - Attaches a `DM` to a `DMSWARM`
837 
838   Collective
839 
840   Input Parameters:
841 + dm     - a `DMSWARM`
842 - dmcell - the `DM` to attach to the `DMSWARM`
843 
844   Level: beginner
845 
846   Note:
847   The attached `DM` (dmcell) will be queried for point location and
848   neighbor MPI-rank information if `DMSwarmMigrate()` is called.
849 
850 .seealso: `DM`, `DMSWARM`, `DMSwarmSetType()`, `DMSwarmGetCellDM()`, `DMSwarmMigrate()`
851 @*/
852 PetscErrorCode DMSwarmSetCellDM(DM dm, DM dmcell)
853 {
854   DM_Swarm *swarm = (DM_Swarm *)dm->data;
855 
856   PetscFunctionBegin;
857   swarm->dmcell = dmcell;
858   PetscFunctionReturn(PETSC_SUCCESS);
859 }
860 
861 /*@
862   DMSwarmGetCellDM - Fetches the attached cell `DM`
863 
864   Collective
865 
866   Input Parameter:
867 . dm - a `DMSWARM`
868 
869   Output Parameter:
870 . dmcell - the `DM` which was attached to the `DMSWARM`
871 
872   Level: beginner
873 
874 .seealso: `DM`, `DMSWARM`, `DMSwarmSetCellDM()`
875 @*/
876 PetscErrorCode DMSwarmGetCellDM(DM dm, DM *dmcell)
877 {
878   DM_Swarm *swarm = (DM_Swarm *)dm->data;
879 
880   PetscFunctionBegin;
881   *dmcell = swarm->dmcell;
882   PetscFunctionReturn(PETSC_SUCCESS);
883 }
884 
885 /*@
886   DMSwarmGetLocalSize - Retrieves the local length of fields registered
887 
888   Not Collective
889 
890   Input Parameter:
891 . dm - a `DMSWARM`
892 
893   Output Parameter:
894 . nlocal - the length of each registered field
895 
896   Level: beginner
897 
898 .seealso: `DM`, `DMSWARM`, `DMSwarmGetSize()`, `DMSwarmSetLocalSizes()`
899 @*/
900 PetscErrorCode DMSwarmGetLocalSize(DM dm, PetscInt *nlocal)
901 {
902   DM_Swarm *swarm = (DM_Swarm *)dm->data;
903 
904   PetscFunctionBegin;
905   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, nlocal, NULL, NULL));
906   PetscFunctionReturn(PETSC_SUCCESS);
907 }
908 
909 /*@
910   DMSwarmGetSize - Retrieves the total length of fields registered
911 
912   Collective
913 
914   Input Parameter:
915 . dm - a `DMSWARM`
916 
917   Output Parameter:
918 . n - the total length of each registered field
919 
920   Level: beginner
921 
922   Note:
923   This calls `MPI_Allreduce()` upon each call (inefficient but safe)
924 
925 .seealso: `DM`, `DMSWARM`, `DMSwarmGetLocalSize()`, `DMSwarmSetLocalSizes()`
926 @*/
927 PetscErrorCode DMSwarmGetSize(DM dm, PetscInt *n)
928 {
929   DM_Swarm *swarm = (DM_Swarm *)dm->data;
930   PetscInt  nlocal;
931 
932   PetscFunctionBegin;
933   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &nlocal, NULL, NULL));
934   PetscCall(MPIU_Allreduce(&nlocal, n, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm)));
935   PetscFunctionReturn(PETSC_SUCCESS);
936 }
937 
938 /*@C
939   DMSwarmRegisterPetscDatatypeField - Register a field to a `DMSWARM` with a native PETSc data type
940 
941   Collective
942 
943   Input Parameters:
944 + dm        - a `DMSWARM`
945 . fieldname - the textual name to identify this field
946 . blocksize - the number of each data type
947 - type      - a valid PETSc data type (`PETSC_CHAR`, `PETSC_SHORT`, `PETSC_INT`, `PETSC_FLOAT`, `PETSC_REAL`, `PETSC_LONG`)
948 
949   Level: beginner
950 
951   Notes:
952   The textual name for each registered field must be unique.
953 
954 .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterUserStructField()`, `DMSwarmRegisterUserDatatypeField()`
955 @*/
956 PetscErrorCode DMSwarmRegisterPetscDatatypeField(DM dm, const char fieldname[], PetscInt blocksize, PetscDataType type)
957 {
958   DM_Swarm *swarm = (DM_Swarm *)dm->data;
959   size_t    size;
960 
961   PetscFunctionBegin;
962   PetscCheck(swarm->field_registration_initialized, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMSwarmInitializeFieldRegister() first");
963   PetscCheck(!swarm->field_registration_finalized, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Cannot register additional fields after calling DMSwarmFinalizeFieldRegister() first");
964 
965   PetscCheck(type != PETSC_OBJECT, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
966   PetscCheck(type != PETSC_FUNCTION, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
967   PetscCheck(type != PETSC_STRING, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
968   PetscCheck(type != PETSC_STRUCT, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
969   PetscCheck(type != PETSC_DATATYPE_UNKNOWN, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
970 
971   PetscCall(PetscDataTypeGetSize(type, &size));
972   /* Load a specific data type into data bucket, specifying textual name and its size in bytes */
973   PetscCall(DMSwarmDataBucketRegisterField(swarm->db, "DMSwarmRegisterPetscDatatypeField", fieldname, blocksize * size, NULL));
974   {
975     DMSwarmDataField gfield;
976 
977     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
978     PetscCall(DMSwarmDataFieldSetBlockSize(gfield, blocksize));
979   }
980   swarm->db->field[swarm->db->nfields - 1]->petsc_type = type;
981   PetscFunctionReturn(PETSC_SUCCESS);
982 }
983 
984 /*@C
985   DMSwarmRegisterUserStructField - Register a user defined struct to a `DMSWARM`
986 
987   Collective
988 
989   Input Parameters:
990 + dm        - a `DMSWARM`
991 . fieldname - the textual name to identify this field
992 - size      - the size in bytes of the user struct of each data type
993 
994   Level: beginner
995 
996   Note:
997   The textual name for each registered field must be unique.
998 
999 .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmRegisterUserDatatypeField()`
1000 @*/
1001 PetscErrorCode DMSwarmRegisterUserStructField(DM dm, const char fieldname[], size_t size)
1002 {
1003   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1004 
1005   PetscFunctionBegin;
1006   PetscCall(DMSwarmDataBucketRegisterField(swarm->db, "DMSwarmRegisterUserStructField", fieldname, size, NULL));
1007   swarm->db->field[swarm->db->nfields - 1]->petsc_type = PETSC_STRUCT;
1008   PetscFunctionReturn(PETSC_SUCCESS);
1009 }
1010 
1011 /*@C
1012   DMSwarmRegisterUserDatatypeField - Register a user defined data type to a `DMSWARM`
1013 
1014   Collective
1015 
1016   Input Parameters:
1017 + dm        - a `DMSWARM`
1018 . fieldname - the textual name to identify this field
1019 . size      - the size in bytes of the user data type
1020 - blocksize - the number of each data type
1021 
1022   Level: beginner
1023 
1024   Note:
1025   The textual name for each registered field must be unique.
1026 
1027 .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmRegisterUserStructField()`
1028 @*/
1029 PetscErrorCode DMSwarmRegisterUserDatatypeField(DM dm, const char fieldname[], size_t size, PetscInt blocksize)
1030 {
1031   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1032 
1033   PetscFunctionBegin;
1034   PetscCall(DMSwarmDataBucketRegisterField(swarm->db, "DMSwarmRegisterUserDatatypeField", fieldname, blocksize * size, NULL));
1035   {
1036     DMSwarmDataField gfield;
1037 
1038     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
1039     PetscCall(DMSwarmDataFieldSetBlockSize(gfield, blocksize));
1040   }
1041   swarm->db->field[swarm->db->nfields - 1]->petsc_type = PETSC_DATATYPE_UNKNOWN;
1042   PetscFunctionReturn(PETSC_SUCCESS);
1043 }
1044 
1045 /*@C
1046   DMSwarmGetField - Get access to the underlying array storing all entries associated with a registered field
1047 
1048   Not Collective
1049 
1050   Input Parameters:
1051 + dm        - a `DMSWARM`
1052 - fieldname - the textual name to identify this field
1053 
1054   Output Parameters:
1055 + blocksize - the number of each data type
1056 . type      - the data type
1057 - data      - pointer to raw array
1058 
1059   Level: beginner
1060 
1061   Notes:
1062   The array must be returned using a matching call to `DMSwarmRestoreField()`.
1063 
1064 .seealso: `DM`, `DMSWARM`, `DMSwarmRestoreField()`
1065 @*/
1066 PetscErrorCode DMSwarmGetField(DM dm, const char fieldname[], PetscInt *blocksize, PetscDataType *type, void **data)
1067 {
1068   DM_Swarm        *swarm = (DM_Swarm *)dm->data;
1069   DMSwarmDataField gfield;
1070 
1071   PetscFunctionBegin;
1072   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
1073   if (!swarm->issetup) PetscCall(DMSetUp(dm));
1074   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
1075   PetscCall(DMSwarmDataFieldGetAccess(gfield));
1076   PetscCall(DMSwarmDataFieldGetEntries(gfield, data));
1077   if (blocksize) *blocksize = gfield->bs;
1078   if (type) *type = gfield->petsc_type;
1079   PetscFunctionReturn(PETSC_SUCCESS);
1080 }
1081 
1082 /*@C
1083   DMSwarmRestoreField - Restore access to the underlying array storing all entries associated with a registered field
1084 
1085   Not Collective
1086 
1087   Input Parameters:
1088 + dm        - a `DMSWARM`
1089 - fieldname - the textual name to identify this field
1090 
1091   Output Parameters:
1092 + blocksize - the number of each data type
1093 . type      - the data type
1094 - data      - pointer to raw array
1095 
1096   Level: beginner
1097 
1098   Notes:
1099   The user must call `DMSwarmGetField()` prior to calling `DMSwarmRestoreField()`.
1100 
1101 .seealso: `DM`, `DMSWARM`, `DMSwarmGetField()`
1102 @*/
1103 PetscErrorCode DMSwarmRestoreField(DM dm, const char fieldname[], PetscInt *blocksize, PetscDataType *type, void **data)
1104 {
1105   DM_Swarm        *swarm = (DM_Swarm *)dm->data;
1106   DMSwarmDataField gfield;
1107 
1108   PetscFunctionBegin;
1109   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
1110   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
1111   PetscCall(DMSwarmDataFieldRestoreAccess(gfield));
1112   if (data) *data = NULL;
1113   PetscFunctionReturn(PETSC_SUCCESS);
1114 }
1115 
1116 /*@
1117   DMSwarmAddPoint - Add space for one new point in the `DMSWARM`
1118 
1119   Not Collective
1120 
1121   Input Parameter:
1122 . dm - a `DMSWARM`
1123 
1124   Level: beginner
1125 
1126   Notes:
1127   The new point will have all fields initialized to zero.
1128 
1129 .seealso: `DM`, `DMSWARM`, `DMSwarmAddNPoints()`
1130 @*/
1131 PetscErrorCode DMSwarmAddPoint(DM dm)
1132 {
1133   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1134 
1135   PetscFunctionBegin;
1136   if (!swarm->issetup) PetscCall(DMSetUp(dm));
1137   PetscCall(PetscLogEventBegin(DMSWARM_AddPoints, 0, 0, 0, 0));
1138   PetscCall(DMSwarmDataBucketAddPoint(swarm->db));
1139   PetscCall(PetscLogEventEnd(DMSWARM_AddPoints, 0, 0, 0, 0));
1140   PetscFunctionReturn(PETSC_SUCCESS);
1141 }
1142 
1143 /*@
1144   DMSwarmAddNPoints - Add space for a number of new points in the `DMSWARM`
1145 
1146   Not Collective
1147 
1148   Input Parameters:
1149 + dm      - a `DMSWARM`
1150 - npoints - the number of new points to add
1151 
1152   Level: beginner
1153 
1154   Notes:
1155   The new point will have all fields initialized to zero.
1156 
1157 .seealso: `DM`, `DMSWARM`, `DMSwarmAddPoint()`
1158 @*/
1159 PetscErrorCode DMSwarmAddNPoints(DM dm, PetscInt npoints)
1160 {
1161   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1162   PetscInt  nlocal;
1163 
1164   PetscFunctionBegin;
1165   PetscCall(PetscLogEventBegin(DMSWARM_AddPoints, 0, 0, 0, 0));
1166   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &nlocal, NULL, NULL));
1167   nlocal = nlocal + npoints;
1168   PetscCall(DMSwarmDataBucketSetSizes(swarm->db, nlocal, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
1169   PetscCall(PetscLogEventEnd(DMSWARM_AddPoints, 0, 0, 0, 0));
1170   PetscFunctionReturn(PETSC_SUCCESS);
1171 }
1172 
1173 /*@
1174   DMSwarmRemovePoint - Remove the last point from the `DMSWARM`
1175 
1176   Not Collective
1177 
1178   Input Parameter:
1179 . dm - a `DMSWARM`
1180 
1181   Level: beginner
1182 
1183 .seealso: `DM`, `DMSWARM`, `DMSwarmRemovePointAtIndex()`
1184 @*/
1185 PetscErrorCode DMSwarmRemovePoint(DM dm)
1186 {
1187   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1188 
1189   PetscFunctionBegin;
1190   PetscCall(PetscLogEventBegin(DMSWARM_RemovePoints, 0, 0, 0, 0));
1191   PetscCall(DMSwarmDataBucketRemovePoint(swarm->db));
1192   PetscCall(PetscLogEventEnd(DMSWARM_RemovePoints, 0, 0, 0, 0));
1193   PetscFunctionReturn(PETSC_SUCCESS);
1194 }
1195 
1196 /*@
1197   DMSwarmRemovePointAtIndex - Removes a specific point from the `DMSWARM`
1198 
1199   Not Collective
1200 
1201   Input Parameters:
1202 + dm  - a `DMSWARM`
1203 - idx - index of point to remove
1204 
1205   Level: beginner
1206 
1207 .seealso: `DM`, `DMSWARM`, `DMSwarmRemovePoint()`
1208 @*/
1209 PetscErrorCode DMSwarmRemovePointAtIndex(DM dm, PetscInt idx)
1210 {
1211   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1212 
1213   PetscFunctionBegin;
1214   PetscCall(PetscLogEventBegin(DMSWARM_RemovePoints, 0, 0, 0, 0));
1215   PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db, idx));
1216   PetscCall(PetscLogEventEnd(DMSWARM_RemovePoints, 0, 0, 0, 0));
1217   PetscFunctionReturn(PETSC_SUCCESS);
1218 }
1219 
1220 /*@
1221   DMSwarmCopyPoint - Copy point pj to point pi in the `DMSWARM`
1222 
1223   Not Collective
1224 
1225   Input Parameters:
1226 + dm - a `DMSWARM`
1227 . pi - the index of the point to copy
1228 - pj - the point index where the copy should be located
1229 
1230   Level: beginner
1231 
1232 .seealso: `DM`, `DMSWARM`, `DMSwarmRemovePoint()`
1233 @*/
1234 PetscErrorCode DMSwarmCopyPoint(DM dm, PetscInt pi, PetscInt pj)
1235 {
1236   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1237 
1238   PetscFunctionBegin;
1239   if (!swarm->issetup) PetscCall(DMSetUp(dm));
1240   PetscCall(DMSwarmDataBucketCopyPoint(swarm->db, pi, swarm->db, pj));
1241   PetscFunctionReturn(PETSC_SUCCESS);
1242 }
1243 
1244 static PetscErrorCode DMSwarmMigrate_Basic(DM dm, PetscBool remove_sent_points)
1245 {
1246   PetscFunctionBegin;
1247   PetscCall(DMSwarmMigrate_Push_Basic(dm, remove_sent_points));
1248   PetscFunctionReturn(PETSC_SUCCESS);
1249 }
1250 
1251 /*@
1252   DMSwarmMigrate - Relocates points defined in the `DMSWARM` to other MPI-ranks
1253 
1254   Collective
1255 
1256   Input Parameters:
1257 + dm                 - the `DMSWARM`
1258 - remove_sent_points - flag indicating if sent points should be removed from the current MPI-rank
1259 
1260   Level: advanced
1261 
1262   Notes:
1263   The `DM` will be modified to accommodate received points.
1264   If `remove_sent_points` is `PETSC_TRUE`, any points that were sent will be removed from the `DM`.
1265   Different styles of migration are supported. See `DMSwarmSetMigrateType()`.
1266 
1267 .seealso: `DM`, `DMSWARM`, `DMSwarmSetMigrateType()`
1268 @*/
1269 PetscErrorCode DMSwarmMigrate(DM dm, PetscBool remove_sent_points)
1270 {
1271   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1272 
1273   PetscFunctionBegin;
1274   PetscCall(PetscLogEventBegin(DMSWARM_Migrate, 0, 0, 0, 0));
1275   switch (swarm->migrate_type) {
1276   case DMSWARM_MIGRATE_BASIC:
1277     PetscCall(DMSwarmMigrate_Basic(dm, remove_sent_points));
1278     break;
1279   case DMSWARM_MIGRATE_DMCELLNSCATTER:
1280     PetscCall(DMSwarmMigrate_CellDMScatter(dm, remove_sent_points));
1281     break;
1282   case DMSWARM_MIGRATE_DMCELLEXACT:
1283     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_MIGRATE_DMCELLEXACT not implemented");
1284   case DMSWARM_MIGRATE_USER:
1285     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_MIGRATE_USER not implemented");
1286   default:
1287     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_MIGRATE type unknown");
1288   }
1289   PetscCall(PetscLogEventEnd(DMSWARM_Migrate, 0, 0, 0, 0));
1290   PetscCall(DMClearGlobalVectors(dm));
1291   PetscFunctionReturn(PETSC_SUCCESS);
1292 }
1293 
1294 PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm, PetscInt *globalsize);
1295 
1296 /*
1297  DMSwarmCollectViewCreate
1298 
1299  * Applies a collection method and gathers point neighbour points into dm
1300 
1301  Notes:
1302  Users should call DMSwarmCollectViewDestroy() after
1303  they have finished computations associated with the collected points
1304 */
1305 
1306 /*@
1307   DMSwarmCollectViewCreate - Applies a collection method and gathers points
1308   in neighbour ranks into the `DMSWARM`
1309 
1310   Collective
1311 
1312   Input Parameter:
1313 . dm - the `DMSWARM`
1314 
1315   Level: advanced
1316 
1317   Notes:
1318   Users should call `DMSwarmCollectViewDestroy()` after
1319   they have finished computations associated with the collected points
1320   Different collect methods are supported. See `DMSwarmSetCollectType()`.
1321 
1322 .seealso: `DM`, `DMSWARM`, `DMSwarmCollectViewDestroy()`, `DMSwarmSetCollectType()`
1323 @*/
1324 PetscErrorCode DMSwarmCollectViewCreate(DM dm)
1325 {
1326   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1327   PetscInt  ng;
1328 
1329   PetscFunctionBegin;
1330   PetscCheck(!swarm->collect_view_active, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "CollectView currently active");
1331   PetscCall(DMSwarmGetLocalSize(dm, &ng));
1332   switch (swarm->collect_type) {
1333   case DMSWARM_COLLECT_BASIC:
1334     PetscCall(DMSwarmMigrate_GlobalToLocal_Basic(dm, &ng));
1335     break;
1336   case DMSWARM_COLLECT_DMDABOUNDINGBOX:
1337     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_COLLECT_DMDABOUNDINGBOX not implemented");
1338   case DMSWARM_COLLECT_GENERAL:
1339     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_COLLECT_GENERAL not implemented");
1340   default:
1341     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_COLLECT type unknown");
1342   }
1343   swarm->collect_view_active       = PETSC_TRUE;
1344   swarm->collect_view_reset_nlocal = ng;
1345   PetscFunctionReturn(PETSC_SUCCESS);
1346 }
1347 
1348 /*@
1349   DMSwarmCollectViewDestroy - Resets the `DMSWARM` to the size prior to calling `DMSwarmCollectViewCreate()`
1350 
1351   Collective
1352 
1353   Input Parameters:
1354 . dm - the `DMSWARM`
1355 
1356   Notes:
1357   Users should call `DMSwarmCollectViewCreate()` before this function is called.
1358 
1359   Level: advanced
1360 
1361 .seealso: `DM`, `DMSWARM`, `DMSwarmCollectViewCreate()`, `DMSwarmSetCollectType()`
1362 @*/
1363 PetscErrorCode DMSwarmCollectViewDestroy(DM dm)
1364 {
1365   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1366 
1367   PetscFunctionBegin;
1368   PetscCheck(swarm->collect_view_active, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "CollectView is currently not active");
1369   PetscCall(DMSwarmSetLocalSizes(dm, swarm->collect_view_reset_nlocal, -1));
1370   swarm->collect_view_active = PETSC_FALSE;
1371   PetscFunctionReturn(PETSC_SUCCESS);
1372 }
1373 
1374 static PetscErrorCode DMSwarmSetUpPIC(DM dm)
1375 {
1376   PetscInt dim;
1377 
1378   PetscFunctionBegin;
1379   PetscCall(DMSwarmSetNumSpecies(dm, 1));
1380   PetscCall(DMGetDimension(dm, &dim));
1381   PetscCheck(dim >= 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Dimension must be 1,2,3 - found %" PetscInt_FMT, dim);
1382   PetscCheck(dim <= 3, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Dimension must be 1,2,3 - found %" PetscInt_FMT, dim);
1383   PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmPICField_coor, dim, PETSC_DOUBLE));
1384   PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmPICField_cellid, 1, PETSC_INT));
1385   PetscFunctionReturn(PETSC_SUCCESS);
1386 }
1387 
1388 /*@
1389   DMSwarmSetPointCoordinatesRandom - Sets initial coordinates for particles in each cell
1390 
1391   Collective
1392 
1393   Input Parameters:
1394 + dm  - the `DMSWARM`
1395 - Npc - The number of particles per cell in the cell `DM`
1396 
1397   Level: intermediate
1398 
1399   Notes:
1400   The user must use `DMSwarmSetCellDM()` to set the cell `DM` first. The particles are placed randomly inside each cell. If only
1401   one particle is in each cell, it is placed at the centroid.
1402 
1403 .seealso: `DM`, `DMSWARM`, `DMSwarmSetCellDM()`
1404 @*/
1405 PetscErrorCode DMSwarmSetPointCoordinatesRandom(DM dm, PetscInt Npc)
1406 {
1407   DM             cdm;
1408   PetscRandom    rnd;
1409   DMPolytopeType ct;
1410   PetscBool      simplex;
1411   PetscReal     *centroid, *coords, *xi0, *v0, *J, *invJ, detJ;
1412   PetscInt       dim, d, cStart, cEnd, c, p;
1413 
1414   PetscFunctionBeginUser;
1415   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rnd));
1416   PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0));
1417   PetscCall(PetscRandomSetType(rnd, PETSCRAND48));
1418 
1419   PetscCall(DMSwarmGetCellDM(dm, &cdm));
1420   PetscCall(DMGetDimension(cdm, &dim));
1421   PetscCall(DMPlexGetHeightStratum(cdm, 0, &cStart, &cEnd));
1422   PetscCall(DMPlexGetCellType(cdm, cStart, &ct));
1423   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
1424 
1425   PetscCall(PetscMalloc5(dim, &centroid, dim, &xi0, dim, &v0, dim * dim, &J, dim * dim, &invJ));
1426   for (d = 0; d < dim; ++d) xi0[d] = -1.0;
1427   PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
1428   for (c = cStart; c < cEnd; ++c) {
1429     if (Npc == 1) {
1430       PetscCall(DMPlexComputeCellGeometryFVM(cdm, c, NULL, centroid, NULL));
1431       for (d = 0; d < dim; ++d) coords[c * dim + d] = centroid[d];
1432     } else {
1433       PetscCall(DMPlexComputeCellGeometryFEM(cdm, c, NULL, v0, J, invJ, &detJ)); /* affine */
1434       for (p = 0; p < Npc; ++p) {
1435         const PetscInt n   = c * Npc + p;
1436         PetscReal      sum = 0.0, refcoords[3];
1437 
1438         for (d = 0; d < dim; ++d) {
1439           PetscCall(PetscRandomGetValueReal(rnd, &refcoords[d]));
1440           sum += refcoords[d];
1441         }
1442         if (simplex && sum > 0.0)
1443           for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim) * sum;
1444         CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n * dim]);
1445       }
1446     }
1447   }
1448   PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
1449   PetscCall(PetscFree5(centroid, xi0, v0, J, invJ));
1450   PetscCall(PetscRandomDestroy(&rnd));
1451   PetscFunctionReturn(PETSC_SUCCESS);
1452 }
1453 
1454 /*@
1455   DMSwarmSetType - Set particular flavor of `DMSWARM`
1456 
1457   Collective
1458 
1459   Input Parameters:
1460 + dm    - the `DMSWARM`
1461 - stype - the `DMSWARM` type (e.g. `DMSWARM_PIC`)
1462 
1463   Level: advanced
1464 
1465 .seealso: `DM`, `DMSWARM`, `DMSwarmSetMigrateType()`, `DMSwarmSetCollectType()`, `DMSwarmType`, `DMSWARM_PIC`, `DMSWARM_BASIC`
1466 @*/
1467 PetscErrorCode DMSwarmSetType(DM dm, DMSwarmType stype)
1468 {
1469   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1470 
1471   PetscFunctionBegin;
1472   swarm->swarm_type = stype;
1473   if (swarm->swarm_type == DMSWARM_PIC) PetscCall(DMSwarmSetUpPIC(dm));
1474   PetscFunctionReturn(PETSC_SUCCESS);
1475 }
1476 
1477 static PetscErrorCode DMSetup_Swarm(DM dm)
1478 {
1479   DM_Swarm   *swarm = (DM_Swarm *)dm->data;
1480   PetscMPIInt rank;
1481   PetscInt    p, npoints, *rankval;
1482 
1483   PetscFunctionBegin;
1484   if (swarm->issetup) PetscFunctionReturn(PETSC_SUCCESS);
1485   swarm->issetup = PETSC_TRUE;
1486 
1487   if (swarm->swarm_type == DMSWARM_PIC) {
1488     /* check dmcell exists */
1489     PetscCheck(swarm->dmcell, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSWARM_PIC requires you call DMSwarmSetCellDM");
1490 
1491     if (swarm->dmcell->ops->locatepointssubdomain) {
1492       /* check methods exists for exact ownership identificiation */
1493       PetscCall(PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->ops->LocatePointsSubdomain\n"));
1494       swarm->migrate_type = DMSWARM_MIGRATE_DMCELLEXACT;
1495     } else {
1496       /* check methods exist for point location AND rank neighbor identification */
1497       if (swarm->dmcell->ops->locatepoints) {
1498         PetscCall(PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->LocatePoints\n"));
1499       } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSWARM_PIC requires the method CellDM->ops->locatepoints be defined");
1500 
1501       if (swarm->dmcell->ops->getneighbors) {
1502         PetscCall(PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->GetNeigbors\n"));
1503       } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSWARM_PIC requires the method CellDM->ops->getneighbors be defined");
1504 
1505       swarm->migrate_type = DMSWARM_MIGRATE_DMCELLNSCATTER;
1506     }
1507   }
1508 
1509   PetscCall(DMSwarmFinalizeFieldRegister(dm));
1510 
1511   /* check some fields were registered */
1512   PetscCheck(swarm->db->nfields > 2, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "At least one field user must be registered via DMSwarmRegisterXXX()");
1513 
1514   /* check local sizes were set */
1515   PetscCheck(swarm->db->L != -1, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Local sizes must be set via DMSwarmSetLocalSizes()");
1516 
1517   /* initialize values in pid and rank placeholders */
1518   /* TODO: [pid - use MPI_Scan] */
1519   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
1520   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL));
1521   PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
1522   for (p = 0; p < npoints; p++) rankval[p] = (PetscInt)rank;
1523   PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
1524   PetscFunctionReturn(PETSC_SUCCESS);
1525 }
1526 
1527 extern PetscErrorCode DMSwarmSortDestroy(DMSwarmSort *_ctx);
1528 
1529 static PetscErrorCode DMDestroy_Swarm(DM dm)
1530 {
1531   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1532 
1533   PetscFunctionBegin;
1534   if (--swarm->refct > 0) PetscFunctionReturn(PETSC_SUCCESS);
1535   PetscCall(DMSwarmDataBucketDestroy(&swarm->db));
1536   if (swarm->sort_context) PetscCall(DMSwarmSortDestroy(&swarm->sort_context));
1537   PetscCall(PetscFree(swarm));
1538   PetscFunctionReturn(PETSC_SUCCESS);
1539 }
1540 
1541 static PetscErrorCode DMSwarmView_Draw(DM dm, PetscViewer viewer)
1542 {
1543   DM         cdm;
1544   PetscDraw  draw;
1545   PetscReal *coords, oldPause, radius = 0.01;
1546   PetscInt   Np, p, bs;
1547 
1548   PetscFunctionBegin;
1549   PetscCall(PetscOptionsGetReal(NULL, ((PetscObject)dm)->prefix, "-dm_view_swarm_radius", &radius, NULL));
1550   PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
1551   PetscCall(DMSwarmGetCellDM(dm, &cdm));
1552   PetscCall(PetscDrawGetPause(draw, &oldPause));
1553   PetscCall(PetscDrawSetPause(draw, 0.0));
1554   PetscCall(DMView(cdm, viewer));
1555   PetscCall(PetscDrawSetPause(draw, oldPause));
1556 
1557   PetscCall(DMSwarmGetLocalSize(dm, &Np));
1558   PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, &bs, NULL, (void **)&coords));
1559   for (p = 0; p < Np; ++p) {
1560     const PetscInt i = p * bs;
1561 
1562     PetscCall(PetscDrawEllipse(draw, coords[i], coords[i + 1], radius, radius, PETSC_DRAW_BLUE));
1563   }
1564   PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, &bs, NULL, (void **)&coords));
1565   PetscCall(PetscDrawFlush(draw));
1566   PetscCall(PetscDrawPause(draw));
1567   PetscCall(PetscDrawSave(draw));
1568   PetscFunctionReturn(PETSC_SUCCESS);
1569 }
1570 
1571 static PetscErrorCode DMView_Swarm_Ascii(DM dm, PetscViewer viewer)
1572 {
1573   PetscViewerFormat format;
1574   PetscInt         *sizes;
1575   PetscInt          dim, Np, maxSize = 17;
1576   MPI_Comm          comm;
1577   PetscMPIInt       rank, size, p;
1578   const char       *name;
1579 
1580   PetscFunctionBegin;
1581   PetscCall(PetscViewerGetFormat(viewer, &format));
1582   PetscCall(DMGetDimension(dm, &dim));
1583   PetscCall(DMSwarmGetLocalSize(dm, &Np));
1584   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1585   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1586   PetscCallMPI(MPI_Comm_size(comm, &size));
1587   PetscCall(PetscObjectGetName((PetscObject)dm, &name));
1588   if (name) PetscCall(PetscViewerASCIIPrintf(viewer, "%s in %" PetscInt_FMT " dimension%s:\n", name, dim, dim == 1 ? "" : "s"));
1589   else PetscCall(PetscViewerASCIIPrintf(viewer, "Swarm in %" PetscInt_FMT " dimension%s:\n", dim, dim == 1 ? "" : "s"));
1590   if (size < maxSize) PetscCall(PetscCalloc1(size, &sizes));
1591   else PetscCall(PetscCalloc1(3, &sizes));
1592   if (size < maxSize) {
1593     PetscCallMPI(MPI_Gather(&Np, 1, MPIU_INT, sizes, 1, MPIU_INT, 0, comm));
1594     PetscCall(PetscViewerASCIIPrintf(viewer, "  Number of particles per rank:"));
1595     for (p = 0; p < size; ++p) {
1596       if (rank == 0) PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, sizes[p]));
1597     }
1598   } else {
1599     PetscInt locMinMax[2] = {Np, Np};
1600 
1601     PetscCall(PetscGlobalMinMaxInt(comm, locMinMax, sizes));
1602     PetscCall(PetscViewerASCIIPrintf(viewer, "  Min/Max of particles per rank: %" PetscInt_FMT "/%" PetscInt_FMT, sizes[0], sizes[1]));
1603   }
1604   PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1605   PetscCall(PetscFree(sizes));
1606   if (format == PETSC_VIEWER_ASCII_INFO) {
1607     PetscInt *cell;
1608 
1609     PetscCall(PetscViewerASCIIPrintf(viewer, "  Cells containing each particle:\n"));
1610     PetscCall(PetscViewerASCIIPushSynchronized(viewer));
1611     PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&cell));
1612     for (p = 0; p < Np; ++p) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "  p%d: %" PetscInt_FMT "\n", p, cell[p]));
1613     PetscCall(PetscViewerFlush(viewer));
1614     PetscCall(PetscViewerASCIIPopSynchronized(viewer));
1615     PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&cell));
1616   }
1617   PetscFunctionReturn(PETSC_SUCCESS);
1618 }
1619 
1620 static PetscErrorCode DMView_Swarm(DM dm, PetscViewer viewer)
1621 {
1622   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1623   PetscBool iascii, ibinary, isvtk, isdraw;
1624 #if defined(PETSC_HAVE_HDF5)
1625   PetscBool ishdf5;
1626 #endif
1627 
1628   PetscFunctionBegin;
1629   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1630   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1631   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1632   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &ibinary));
1633   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk));
1634 #if defined(PETSC_HAVE_HDF5)
1635   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
1636 #endif
1637   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
1638   if (iascii) {
1639     PetscViewerFormat format;
1640 
1641     PetscCall(PetscViewerGetFormat(viewer, &format));
1642     switch (format) {
1643     case PETSC_VIEWER_ASCII_INFO_DETAIL:
1644       PetscCall(DMSwarmDataBucketView(PetscObjectComm((PetscObject)dm), swarm->db, NULL, DATABUCKET_VIEW_STDOUT));
1645       break;
1646     default:
1647       PetscCall(DMView_Swarm_Ascii(dm, viewer));
1648     }
1649   } else {
1650 #if defined(PETSC_HAVE_HDF5)
1651     if (ishdf5) PetscCall(DMSwarmView_HDF5(dm, viewer));
1652 #endif
1653     if (isdraw) PetscCall(DMSwarmView_Draw(dm, viewer));
1654   }
1655   PetscFunctionReturn(PETSC_SUCCESS);
1656 }
1657 
1658 /*@C
1659   DMSwarmGetCellSwarm - Extracts a single cell from the `DMSWARM` object, returns it as a single cell `DMSWARM`.
1660   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.
1661 
1662   Noncollective
1663 
1664   Input Parameters:
1665 + sw        - the `DMSWARM`
1666 . cellID    - the integer id of the cell to be extracted and filtered
1667 - cellswarm - The `DMSWARM` to receive the cell
1668 
1669   Level: beginner
1670 
1671   Notes:
1672   This presently only supports `DMSWARM_PIC` type
1673 
1674   Should be restored with `DMSwarmRestoreCellSwarm()`
1675 
1676   Changes to this cell of the swarm will be lost if they are made prior to restoring this cell.
1677 
1678 .seealso: `DM`, `DMSWARM`, `DMSwarmRestoreCellSwarm()`
1679 @*/
1680 PETSC_EXTERN PetscErrorCode DMSwarmGetCellSwarm(DM sw, PetscInt cellID, DM cellswarm)
1681 {
1682   DM_Swarm *original = (DM_Swarm *)sw->data;
1683   DMLabel   label;
1684   DM        dmc, subdmc;
1685   PetscInt *pids, particles, dim;
1686 
1687   PetscFunctionBegin;
1688   /* Configure new swarm */
1689   PetscCall(DMSetType(cellswarm, DMSWARM));
1690   PetscCall(DMGetDimension(sw, &dim));
1691   PetscCall(DMSetDimension(cellswarm, dim));
1692   PetscCall(DMSwarmSetType(cellswarm, DMSWARM_PIC));
1693   /* Destroy the unused, unconfigured data bucket to prevent stragglers in memory */
1694   PetscCall(DMSwarmDataBucketDestroy(&((DM_Swarm *)cellswarm->data)->db));
1695   PetscCall(DMSwarmSortGetAccess(sw));
1696   PetscCall(DMSwarmSortGetNumberOfPointsPerCell(sw, cellID, &particles));
1697   PetscCall(DMSwarmSortGetPointsPerCell(sw, cellID, &particles, &pids));
1698   PetscCall(DMSwarmDataBucketCreateFromSubset(original->db, particles, pids, &((DM_Swarm *)cellswarm->data)->db));
1699   PetscCall(DMSwarmSortRestoreAccess(sw));
1700   PetscCall(PetscFree(pids));
1701   PetscCall(DMSwarmGetCellDM(sw, &dmc));
1702   PetscCall(DMLabelCreate(PetscObjectComm((PetscObject)sw), "singlecell", &label));
1703   PetscCall(DMAddLabel(dmc, label));
1704   PetscCall(DMLabelSetValue(label, cellID, 1));
1705   PetscCall(DMPlexFilter(dmc, label, 1, &subdmc));
1706   PetscCall(DMSwarmSetCellDM(cellswarm, subdmc));
1707   PetscCall(DMLabelDestroy(&label));
1708   PetscFunctionReturn(PETSC_SUCCESS);
1709 }
1710 
1711 /*@C
1712   DMSwarmRestoreCellSwarm - Restores a `DMSWARM` object obtained with `DMSwarmGetCellSwarm()`. All fields are copied back into the parent swarm.
1713 
1714   Noncollective
1715 
1716   Input Parameters:
1717 + sw        - the parent `DMSWARM`
1718 . cellID    - the integer id of the cell to be copied back into the parent swarm
1719 - cellswarm - the cell swarm object
1720 
1721   Level: beginner
1722 
1723   Note:
1724   This only supports `DMSWARM_PIC` types of `DMSWARM`s
1725 
1726 .seealso: `DM`, `DMSWARM`, `DMSwarmGetCellSwarm()`
1727 @*/
1728 PETSC_EXTERN PetscErrorCode DMSwarmRestoreCellSwarm(DM sw, PetscInt cellID, DM cellswarm)
1729 {
1730   DM        dmc;
1731   PetscInt *pids, particles, p;
1732 
1733   PetscFunctionBegin;
1734   PetscCall(DMSwarmSortGetAccess(sw));
1735   PetscCall(DMSwarmSortGetPointsPerCell(sw, cellID, &particles, &pids));
1736   PetscCall(DMSwarmSortRestoreAccess(sw));
1737   /* Pointwise copy of each particle based on pid. The parent swarm may not be altered during this process. */
1738   for (p = 0; p < particles; ++p) PetscCall(DMSwarmDataBucketCopyPoint(((DM_Swarm *)cellswarm->data)->db, pids[p], ((DM_Swarm *)sw->data)->db, pids[p]));
1739   /* Free memory, destroy cell dm */
1740   PetscCall(DMSwarmGetCellDM(cellswarm, &dmc));
1741   PetscCall(DMDestroy(&dmc));
1742   PetscCall(PetscFree(pids));
1743   PetscFunctionReturn(PETSC_SUCCESS);
1744 }
1745 
1746 PETSC_INTERN PetscErrorCode DMClone_Swarm(DM, DM *);
1747 
1748 static PetscErrorCode DMInitialize_Swarm(DM sw)
1749 {
1750   PetscFunctionBegin;
1751   sw->dim                           = 0;
1752   sw->ops->view                     = DMView_Swarm;
1753   sw->ops->load                     = NULL;
1754   sw->ops->setfromoptions           = NULL;
1755   sw->ops->clone                    = DMClone_Swarm;
1756   sw->ops->setup                    = DMSetup_Swarm;
1757   sw->ops->createlocalsection       = NULL;
1758   sw->ops->createsectionpermutation = NULL;
1759   sw->ops->createdefaultconstraints = NULL;
1760   sw->ops->createglobalvector       = DMCreateGlobalVector_Swarm;
1761   sw->ops->createlocalvector        = DMCreateLocalVector_Swarm;
1762   sw->ops->getlocaltoglobalmapping  = NULL;
1763   sw->ops->createfieldis            = NULL;
1764   sw->ops->createcoordinatedm       = NULL;
1765   sw->ops->getcoloring              = NULL;
1766   sw->ops->creatematrix             = DMCreateMatrix_Swarm;
1767   sw->ops->createinterpolation      = NULL;
1768   sw->ops->createinjection          = NULL;
1769   sw->ops->createmassmatrix         = DMCreateMassMatrix_Swarm;
1770   sw->ops->refine                   = NULL;
1771   sw->ops->coarsen                  = NULL;
1772   sw->ops->refinehierarchy          = NULL;
1773   sw->ops->coarsenhierarchy         = NULL;
1774   sw->ops->globaltolocalbegin       = NULL;
1775   sw->ops->globaltolocalend         = NULL;
1776   sw->ops->localtoglobalbegin       = NULL;
1777   sw->ops->localtoglobalend         = NULL;
1778   sw->ops->destroy                  = DMDestroy_Swarm;
1779   sw->ops->createsubdm              = NULL;
1780   sw->ops->getdimpoints             = NULL;
1781   sw->ops->locatepoints             = NULL;
1782   PetscFunctionReturn(PETSC_SUCCESS);
1783 }
1784 
1785 PETSC_INTERN PetscErrorCode DMClone_Swarm(DM dm, DM *newdm)
1786 {
1787   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1788 
1789   PetscFunctionBegin;
1790   swarm->refct++;
1791   (*newdm)->data = swarm;
1792   PetscCall(PetscObjectChangeTypeName((PetscObject)*newdm, DMSWARM));
1793   PetscCall(DMInitialize_Swarm(*newdm));
1794   (*newdm)->dim = dm->dim;
1795   PetscFunctionReturn(PETSC_SUCCESS);
1796 }
1797 
1798 /*MC
1799 
1800  DMSWARM = "swarm" - A `DM` object used to represent arrays of data (fields) of arbitrary data type.
1801  This implementation was designed for particle methods in which the underlying
1802  data required to be represented is both (i) dynamic in length, (ii) and of arbitrary data type.
1803 
1804  Level: intermediate
1805 
1806   Notes:
1807  User data can be represented by `DMSWARM` through a registering "fields".
1808  To register a field, the user must provide;
1809  (a) a unique name;
1810  (b) the data type (or size in bytes);
1811  (c) the block size of the data.
1812 
1813  For example, suppose the application requires a unique id, energy, momentum and density to be stored
1814  on a set of particles. Then the following code could be used
1815 .vb
1816     DMSwarmInitializeFieldRegister(dm)
1817     DMSwarmRegisterPetscDatatypeField(dm,"uid",1,PETSC_LONG);
1818     DMSwarmRegisterPetscDatatypeField(dm,"energy",1,PETSC_REAL);
1819     DMSwarmRegisterPetscDatatypeField(dm,"momentum",3,PETSC_REAL);
1820     DMSwarmRegisterPetscDatatypeField(dm,"density",1,PETSC_FLOAT);
1821     DMSwarmFinalizeFieldRegister(dm)
1822 .ve
1823 
1824  The fields represented by `DMSWARM` are dynamic and can be re-sized at any time.
1825  The only restriction imposed by `DMSWARM` is that all fields contain the same number of points.
1826 
1827  To support particle methods, "migration" techniques are provided. These methods migrate data
1828  between ranks.
1829 
1830  `DMSWARM` supports the methods `DMCreateGlobalVector()` and `DMCreateLocalVector()`.
1831  As a `DMSWARM` may internally define and store values of different data types,
1832  before calling `DMCreateGlobalVector()` or `DMCreateLocalVector()`, the user must inform `DMSWARM` which
1833  fields should be used to define a `Vec` object via
1834    `DMSwarmVectorDefineField()`
1835  The specified field can be changed at any time - thereby permitting vectors
1836  compatible with different fields to be created.
1837 
1838  A dual representation of fields in the `DMSWARM` and a Vec object is permitted via
1839    `DMSwarmCreateGlobalVectorFromField()`
1840  Here the data defining the field in the `DMSWARM` is shared with a Vec.
1841  This is inherently unsafe if you alter the size of the field at any time between
1842  calls to `DMSwarmCreateGlobalVectorFromField()` and `DMSwarmDestroyGlobalVectorFromField()`.
1843  If the local size of the `DMSWARM` does not match the local size of the global vector
1844  when `DMSwarmDestroyGlobalVectorFromField()` is called, an error is thrown.
1845 
1846  Additional high-level support is provided for Particle-In-Cell methods.
1847  Please refer to `DMSwarmSetType()`.
1848 
1849 .seealso: `DM`, `DMSWARM`, `DMType`, `DMCreate()`, `DMSetType()`
1850 M*/
1851 PETSC_EXTERN PetscErrorCode DMCreate_Swarm(DM dm)
1852 {
1853   DM_Swarm *swarm;
1854 
1855   PetscFunctionBegin;
1856   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1857   PetscCall(PetscNew(&swarm));
1858   dm->data = swarm;
1859   PetscCall(DMSwarmDataBucketCreate(&swarm->db));
1860   PetscCall(DMSwarmInitializeFieldRegister(dm));
1861   swarm->refct                          = 1;
1862   swarm->vec_field_set                  = PETSC_FALSE;
1863   swarm->issetup                        = PETSC_FALSE;
1864   swarm->swarm_type                     = DMSWARM_BASIC;
1865   swarm->migrate_type                   = DMSWARM_MIGRATE_BASIC;
1866   swarm->collect_type                   = DMSWARM_COLLECT_BASIC;
1867   swarm->migrate_error_on_missing_point = PETSC_FALSE;
1868   swarm->dmcell                         = NULL;
1869   swarm->collect_view_active            = PETSC_FALSE;
1870   swarm->collect_view_reset_nlocal      = -1;
1871   PetscCall(DMInitialize_Swarm(dm));
1872   if (SwarmDataFieldId == -1) PetscCall(PetscObjectComposedDataRegister(&SwarmDataFieldId));
1873   PetscFunctionReturn(PETSC_SUCCESS);
1874 }
1875