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