xref: /petsc/src/dm/impls/swarm/swarm.c (revision 9680b7cedbd3b2fa4b943552cd701a6633a2ad02)
1 #include "petscdmswarm.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 *DMSwarmRemapTypeNames[]     = {"none", "pfak", "colella", "DMSwarmRemapType", "DMSWARM_REMAP_", NULL};
22 const char *DMSwarmPICLayoutTypeNames[] = {"regular", "gauss", "subdivision", NULL};
23 
24 const char DMSwarmField_pid[]     = "DMSwarm_pid";
25 const char DMSwarmField_rank[]    = "DMSwarm_rank";
26 const char DMSwarmPICField_coor[] = "DMSwarmPIC_coor";
27 
28 PetscInt SwarmDataFieldId = -1;
29 
30 #if defined(PETSC_HAVE_HDF5)
31   #include <petscviewerhdf5.h>
32 
33 static PetscErrorCode VecView_Swarm_HDF5_Internal(Vec v, PetscViewer viewer)
34 {
35   DM        dm;
36   PetscReal seqval;
37   PetscInt  seqnum, bs;
38   PetscBool isseq, ists;
39 
40   PetscFunctionBegin;
41   PetscCall(VecGetDM(v, &dm));
42   PetscCall(VecGetBlockSize(v, &bs));
43   PetscCall(PetscViewerHDF5PushGroup(viewer, "/particle_fields"));
44   PetscCall(PetscObjectTypeCompare((PetscObject)v, VECSEQ, &isseq));
45   PetscCall(DMGetOutputSequenceNumber(dm, &seqnum, &seqval));
46   PetscCall(PetscViewerHDF5IsTimestepping(viewer, &ists));
47   if (ists) PetscCall(PetscViewerHDF5SetTimestep(viewer, seqnum));
48   /* PetscCall(DMSequenceView_HDF5(dm, "time", seqnum, (PetscScalar) seqval, viewer)); */
49   PetscCall(VecViewNative(v, viewer));
50   PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject)v, "Nc", PETSC_INT, (void *)&bs));
51   PetscCall(PetscViewerHDF5PopGroup(viewer));
52   PetscFunctionReturn(PETSC_SUCCESS);
53 }
54 
55 static PetscErrorCode DMSwarmView_HDF5(DM dm, PetscViewer viewer)
56 {
57   DMSwarmCellDM celldm;
58   Vec           coordinates;
59   PetscInt      Np, Nfc;
60   PetscBool     isseq;
61   const char  **coordFields;
62 
63   PetscFunctionBegin;
64   PetscCall(DMSwarmGetCellDMActive(dm, &celldm));
65   PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Nfc, &coordFields));
66   PetscCheck(Nfc == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "We only support a single coordinate field right now, not %" PetscInt_FMT, Nfc);
67   PetscCall(DMSwarmGetSize(dm, &Np));
68   PetscCall(DMSwarmCreateGlobalVectorFromField(dm, coordFields[0], &coordinates));
69   PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
70   PetscCall(PetscViewerHDF5PushGroup(viewer, "/particles"));
71   PetscCall(PetscObjectTypeCompare((PetscObject)coordinates, VECSEQ, &isseq));
72   PetscCall(VecViewNative(coordinates, viewer));
73   PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject)coordinates, "Np", PETSC_INT, (void *)&Np));
74   PetscCall(PetscViewerHDF5PopGroup(viewer));
75   PetscCall(DMSwarmDestroyGlobalVectorFromField(dm, coordFields[0], &coordinates));
76   PetscFunctionReturn(PETSC_SUCCESS);
77 }
78 #endif
79 
80 static PetscErrorCode VecView_Swarm(Vec v, PetscViewer viewer)
81 {
82   DM dm;
83 #if defined(PETSC_HAVE_HDF5)
84   PetscBool ishdf5;
85 #endif
86 
87   PetscFunctionBegin;
88   PetscCall(VecGetDM(v, &dm));
89   PetscCheck(dm, PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
90 #if defined(PETSC_HAVE_HDF5)
91   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
92   if (ishdf5) {
93     PetscCall(VecView_Swarm_HDF5_Internal(v, viewer));
94     PetscFunctionReturn(PETSC_SUCCESS);
95   }
96 #endif
97   PetscCall(VecViewNative(v, viewer));
98   PetscFunctionReturn(PETSC_SUCCESS);
99 }
100 
101 /*@C
102   DMSwarmVectorGetField - Gets the fields from which to define a `Vec` object
103   when `DMCreateLocalVector()`, or `DMCreateGlobalVector()` is called
104 
105   Not collective
106 
107   Input Parameter:
108 . sw - a `DMSWARM`
109 
110   Output Parameters:
111 + Nf         - the number of fields
112 - fieldnames - the textual name given to each registered field, or NULL if it has not been set
113 
114   Level: beginner
115 
116 .seealso: `DM`, `DMSWARM`, `DMSwarmVectorDefineField()`, `DMSwarmRegisterPetscDatatypeField()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`
117 @*/
118 PetscErrorCode DMSwarmVectorGetField(DM sw, PetscInt *Nf, const char **fieldnames[])
119 {
120   DMSwarmCellDM celldm;
121 
122   PetscFunctionBegin;
123   PetscValidHeaderSpecificType(sw, DM_CLASSID, 1, DMSWARM);
124   PetscCall(DMSwarmGetCellDMActive(sw, &celldm));
125   PetscCall(DMSwarmCellDMGetFields(celldm, Nf, fieldnames));
126   PetscFunctionReturn(PETSC_SUCCESS);
127 }
128 
129 /*@
130   DMSwarmVectorDefineField - Sets the field from which to define a `Vec` object
131   when `DMCreateLocalVector()`, or `DMCreateGlobalVector()` is called
132 
133   Collective
134 
135   Input Parameters:
136 + dm        - a `DMSWARM`
137 - fieldname - the textual name given to each registered field
138 
139   Level: beginner
140 
141   Notes:
142   The field with name `fieldname` must be defined as having a data type of `PetscScalar`.
143 
144   This function must be called prior to calling `DMCreateLocalVector()`, `DMCreateGlobalVector()`.
145   Multiple calls to `DMSwarmVectorDefineField()` are permitted.
146 
147 .seealso: `DM`, `DMSWARM`, `DMSwarmVectorDefineFields()`, `DMSwarmVectorGetField()`, `DMSwarmRegisterPetscDatatypeField()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`
148 @*/
149 PetscErrorCode DMSwarmVectorDefineField(DM dm, const char fieldname[])
150 {
151   PetscFunctionBegin;
152   PetscCall(DMSwarmVectorDefineFields(dm, 1, &fieldname));
153   PetscFunctionReturn(PETSC_SUCCESS);
154 }
155 
156 /*@C
157   DMSwarmVectorDefineFields - Sets the fields from which to define a `Vec` object
158   when `DMCreateLocalVector()`, or `DMCreateGlobalVector()` is called
159 
160   Collective, No Fortran support
161 
162   Input Parameters:
163 + sw         - a `DMSWARM`
164 . Nf         - the number of fields
165 - fieldnames - the textual name given to each registered field
166 
167   Level: beginner
168 
169   Notes:
170   Each field with name in `fieldnames` must be defined as having a data type of `PetscScalar`.
171 
172   This function must be called prior to calling `DMCreateLocalVector()`, `DMCreateGlobalVector()`.
173   Multiple calls to `DMSwarmVectorDefineField()` are permitted.
174 
175 .seealso: `DM`, `DMSWARM`, `DMSwarmVectorDefineField()`, `DMSwarmVectorGetField()`, `DMSwarmRegisterPetscDatatypeField()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`
176 @*/
177 PetscErrorCode DMSwarmVectorDefineFields(DM sw, PetscInt Nf, const char *fieldnames[])
178 {
179   DM_Swarm     *swarm = (DM_Swarm *)sw->data;
180   DMSwarmCellDM celldm;
181 
182   PetscFunctionBegin;
183   PetscValidHeaderSpecificType(sw, DM_CLASSID, 1, DMSWARM);
184   if (fieldnames) PetscAssertPointer(fieldnames, 3);
185   if (!swarm->issetup) PetscCall(DMSetUp(sw));
186   PetscCheck(Nf >= 0, PetscObjectComm((PetscObject)sw), PETSC_ERR_ARG_OUTOFRANGE, "Number of fields must be non-negative, not %" PetscInt_FMT, Nf);
187   // Create a dummy cell DM if none has been specified (I think we should not support this mode)
188   if (!swarm->activeCellDM) {
189     DM            dm;
190     DMSwarmCellDM celldm;
191 
192     PetscCall(DMCreate(PetscObjectComm((PetscObject)sw), &dm));
193     PetscCall(DMSetType(dm, DMSHELL));
194     PetscCall(PetscObjectSetName((PetscObject)dm, "dummy"));
195     PetscCall(DMSwarmCellDMCreate(dm, 0, NULL, 0, NULL, &celldm));
196     PetscCall(DMDestroy(&dm));
197     PetscCall(DMSwarmAddCellDM(sw, celldm));
198     PetscCall(DMSwarmCellDMDestroy(&celldm));
199     PetscCall(DMSwarmSetCellDMActive(sw, "dummy"));
200   }
201   PetscCall(DMSwarmGetCellDMActive(sw, &celldm));
202   for (PetscInt f = 0; f < celldm->Nf; ++f) PetscCall(PetscFree(celldm->dmFields[f]));
203   PetscCall(PetscFree(celldm->dmFields));
204 
205   celldm->Nf = Nf;
206   PetscCall(PetscMalloc1(Nf, &celldm->dmFields));
207   for (PetscInt f = 0; f < Nf; ++f) {
208     PetscDataType type;
209 
210     // Check all fields are of type PETSC_REAL or PETSC_SCALAR
211     PetscCall(DMSwarmGetFieldInfo(sw, fieldnames[f], NULL, &type));
212     PetscCheck(type == PETSC_REAL, PetscObjectComm((PetscObject)sw), PETSC_ERR_SUP, "Only valid for PETSC_REAL");
213     PetscCall(PetscStrallocpy(fieldnames[f], (char **)&celldm->dmFields[f]));
214   }
215   PetscFunctionReturn(PETSC_SUCCESS);
216 }
217 
218 /* requires DMSwarmDefineFieldVector has been called */
219 static PetscErrorCode DMCreateGlobalVector_Swarm(DM sw, Vec *vec)
220 {
221   DM_Swarm     *swarm = (DM_Swarm *)sw->data;
222   DMSwarmCellDM celldm;
223   Vec           x;
224   char          name[PETSC_MAX_PATH_LEN];
225   PetscInt      bs = 0, n;
226 
227   PetscFunctionBegin;
228   if (!swarm->issetup) PetscCall(DMSetUp(sw));
229   PetscCall(DMSwarmGetCellDMActive(sw, &celldm));
230   PetscCheck(celldm->Nf, PetscObjectComm((PetscObject)sw), PETSC_ERR_USER, "Active cell DM does not define any fields");
231   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &n, NULL, NULL));
232 
233   PetscCall(PetscStrncpy(name, "DMSwarmField", PETSC_MAX_PATH_LEN));
234   for (PetscInt f = 0; f < celldm->Nf; ++f) {
235     PetscInt fbs;
236     PetscCall(PetscStrlcat(name, "_", PETSC_MAX_PATH_LEN));
237     PetscCall(PetscStrlcat(name, celldm->dmFields[f], PETSC_MAX_PATH_LEN));
238     PetscCall(DMSwarmGetFieldInfo(sw, celldm->dmFields[f], &fbs, NULL));
239     bs += fbs;
240   }
241   PetscCall(VecCreate(PetscObjectComm((PetscObject)sw), &x));
242   PetscCall(PetscObjectSetName((PetscObject)x, name));
243   PetscCall(VecSetSizes(x, n * bs, PETSC_DETERMINE));
244   PetscCall(VecSetBlockSize(x, bs));
245   PetscCall(VecSetDM(x, sw));
246   PetscCall(VecSetFromOptions(x));
247   PetscCall(VecSetOperation(x, VECOP_VIEW, (void (*)(void))VecView_Swarm));
248   *vec = x;
249   PetscFunctionReturn(PETSC_SUCCESS);
250 }
251 
252 /* requires DMSwarmDefineFieldVector has been called */
253 static PetscErrorCode DMCreateLocalVector_Swarm(DM sw, Vec *vec)
254 {
255   DM_Swarm     *swarm = (DM_Swarm *)sw->data;
256   DMSwarmCellDM celldm;
257   Vec           x;
258   char          name[PETSC_MAX_PATH_LEN];
259   PetscInt      bs = 0, n;
260 
261   PetscFunctionBegin;
262   if (!swarm->issetup) PetscCall(DMSetUp(sw));
263   PetscCall(DMSwarmGetCellDMActive(sw, &celldm));
264   PetscCheck(celldm->Nf, PetscObjectComm((PetscObject)sw), PETSC_ERR_USER, "Active cell DM does not define any fields");
265   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &n, NULL, NULL));
266 
267   PetscCall(PetscStrncpy(name, "DMSwarmField", PETSC_MAX_PATH_LEN));
268   for (PetscInt f = 0; f < celldm->Nf; ++f) {
269     PetscInt fbs;
270     PetscCall(PetscStrlcat(name, "_", PETSC_MAX_PATH_LEN));
271     PetscCall(PetscStrlcat(name, celldm->dmFields[f], PETSC_MAX_PATH_LEN));
272     PetscCall(DMSwarmGetFieldInfo(sw, celldm->dmFields[f], &fbs, NULL));
273     bs += fbs;
274   }
275   PetscCall(VecCreate(PETSC_COMM_SELF, &x));
276   PetscCall(PetscObjectSetName((PetscObject)x, name));
277   PetscCall(VecSetSizes(x, n * bs, PETSC_DETERMINE));
278   PetscCall(VecSetBlockSize(x, bs));
279   PetscCall(VecSetDM(x, sw));
280   PetscCall(VecSetFromOptions(x));
281   *vec = x;
282   PetscFunctionReturn(PETSC_SUCCESS);
283 }
284 
285 static PetscErrorCode DMSwarmDestroyVectorFromField_Private(DM dm, const char fieldname[], Vec *vec)
286 {
287   DM_Swarm        *swarm = (DM_Swarm *)dm->data;
288   DMSwarmDataField gfield;
289   PetscInt         bs, nlocal, fid = -1, cfid = -2;
290   PetscBool        flg;
291 
292   PetscFunctionBegin;
293   /* check vector is an inplace array */
294   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldIdByName(swarm->db, fieldname, &fid));
295   PetscCall(PetscObjectComposedDataGetInt((PetscObject)*vec, SwarmDataFieldId, cfid, flg));
296   (void)flg; /* avoid compiler warning */
297   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);
298   PetscCall(VecGetLocalSize(*vec, &nlocal));
299   PetscCall(VecGetBlockSize(*vec, &bs));
300   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");
301   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
302   PetscCall(DMSwarmDataFieldRestoreAccess(gfield));
303   PetscCall(VecResetArray(*vec));
304   PetscCall(VecDestroy(vec));
305   PetscFunctionReturn(PETSC_SUCCESS);
306 }
307 
308 static PetscErrorCode DMSwarmCreateVectorFromField_Private(DM dm, const char fieldname[], MPI_Comm comm, Vec *vec)
309 {
310   DM_Swarm     *swarm = (DM_Swarm *)dm->data;
311   PetscDataType type;
312   PetscScalar  *array;
313   PetscInt      bs, n, fid;
314   char          name[PETSC_MAX_PATH_LEN];
315   PetscMPIInt   size;
316   PetscBool     iscuda, iskokkos, iship;
317 
318   PetscFunctionBegin;
319   if (!swarm->issetup) PetscCall(DMSetUp(dm));
320   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &n, NULL, NULL));
321   PetscCall(DMSwarmGetField(dm, fieldname, &bs, &type, (void **)&array));
322   PetscCheck(type == PETSC_REAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only valid for PETSC_REAL");
323 
324   PetscCallMPI(MPI_Comm_size(comm, &size));
325   PetscCall(PetscStrcmp(dm->vectype, VECKOKKOS, &iskokkos));
326   PetscCall(PetscStrcmp(dm->vectype, VECCUDA, &iscuda));
327   PetscCall(PetscStrcmp(dm->vectype, VECHIP, &iship));
328   PetscCall(VecCreate(comm, vec));
329   PetscCall(VecSetSizes(*vec, n * bs, PETSC_DETERMINE));
330   PetscCall(VecSetBlockSize(*vec, bs));
331   if (iskokkos) PetscCall(VecSetType(*vec, VECKOKKOS));
332   else if (iscuda) PetscCall(VecSetType(*vec, VECCUDA));
333   else if (iship) PetscCall(VecSetType(*vec, VECHIP));
334   else PetscCall(VecSetType(*vec, VECSTANDARD));
335   PetscCall(VecPlaceArray(*vec, array));
336 
337   PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "DMSwarmSharedField_%s", fieldname));
338   PetscCall(PetscObjectSetName((PetscObject)*vec, name));
339 
340   /* Set guard */
341   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldIdByName(swarm->db, fieldname, &fid));
342   PetscCall(PetscObjectComposedDataSetInt((PetscObject)*vec, SwarmDataFieldId, fid));
343 
344   PetscCall(VecSetDM(*vec, dm));
345   PetscCall(VecSetOperation(*vec, VECOP_VIEW, (void (*)(void))VecView_Swarm));
346   PetscFunctionReturn(PETSC_SUCCESS);
347 }
348 
349 static PetscErrorCode DMSwarmDestroyVectorFromFields_Private(DM sw, PetscInt Nf, const char *fieldnames[], Vec *vec)
350 {
351   DM_Swarm          *swarm = (DM_Swarm *)sw->data;
352   const PetscScalar *array;
353   PetscInt           bs, n, id = 0, cid = -2;
354   PetscBool          flg;
355 
356   PetscFunctionBegin;
357   for (PetscInt f = 0; f < Nf; ++f) {
358     PetscInt fid;
359 
360     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldIdByName(swarm->db, fieldnames[f], &fid));
361     id += fid;
362   }
363   PetscCall(PetscObjectComposedDataGetInt((PetscObject)*vec, SwarmDataFieldId, cid, flg));
364   (void)flg; /* avoid compiler warning */
365   PetscCheck(cid == id, PetscObjectComm((PetscObject)sw), PETSC_ERR_USER, "Vector being destroyed was not created from DMSwarm field(%s)! %" PetscInt_FMT " != %" PetscInt_FMT, fieldnames[0], cid, id);
366   PetscCall(VecGetLocalSize(*vec, &n));
367   PetscCall(VecGetBlockSize(*vec, &bs));
368   n /= bs;
369   PetscCheck(n == swarm->db->L, PetscObjectComm((PetscObject)sw), PETSC_ERR_USER, "DMSwarm sizes have changed since vector was created - cannot ensure pointers are valid");
370   PetscCall(VecGetArrayRead(*vec, &array));
371   for (PetscInt f = 0, off = 0; f < Nf; ++f) {
372     PetscScalar  *farray;
373     PetscDataType ftype;
374     PetscInt      fbs;
375 
376     PetscCall(DMSwarmGetField(sw, fieldnames[f], &fbs, &ftype, (void **)&farray));
377     PetscCheck(off + fbs <= bs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid blocksize %" PetscInt_FMT " < %" PetscInt_FMT, bs, off + fbs);
378     for (PetscInt i = 0; i < n; ++i) {
379       for (PetscInt b = 0; b < fbs; ++b) farray[i * fbs + b] = array[i * bs + off + b];
380     }
381     off += fbs;
382     PetscCall(DMSwarmRestoreField(sw, fieldnames[f], &fbs, &ftype, (void **)&farray));
383   }
384   PetscCall(VecRestoreArrayRead(*vec, &array));
385   PetscCall(VecDestroy(vec));
386   PetscFunctionReturn(PETSC_SUCCESS);
387 }
388 
389 static PetscErrorCode DMSwarmCreateVectorFromFields_Private(DM sw, PetscInt Nf, const char *fieldnames[], MPI_Comm comm, Vec *vec)
390 {
391   DM_Swarm    *swarm = (DM_Swarm *)sw->data;
392   PetscScalar *array;
393   PetscInt     n, bs = 0, id = 0;
394   char         name[PETSC_MAX_PATH_LEN];
395 
396   PetscFunctionBegin;
397   if (!swarm->issetup) PetscCall(DMSetUp(sw));
398   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &n, NULL, NULL));
399   for (PetscInt f = 0; f < Nf; ++f) {
400     PetscDataType ftype;
401     PetscInt      fbs;
402 
403     PetscCall(DMSwarmGetFieldInfo(sw, fieldnames[f], &fbs, &ftype));
404     PetscCheck(ftype == PETSC_REAL, PetscObjectComm((PetscObject)sw), PETSC_ERR_SUP, "Only valid for PETSC_REAL");
405     bs += fbs;
406   }
407 
408   PetscCall(VecCreate(comm, vec));
409   PetscCall(VecSetSizes(*vec, n * bs, PETSC_DETERMINE));
410   PetscCall(VecSetBlockSize(*vec, bs));
411   PetscCall(VecSetType(*vec, sw->vectype));
412 
413   PetscCall(VecGetArrayWrite(*vec, &array));
414   for (PetscInt f = 0, off = 0; f < Nf; ++f) {
415     PetscScalar  *farray;
416     PetscDataType ftype;
417     PetscInt      fbs;
418 
419     PetscCall(DMSwarmGetField(sw, fieldnames[f], &fbs, &ftype, (void **)&farray));
420     for (PetscInt i = 0; i < n; ++i) {
421       for (PetscInt b = 0; b < fbs; ++b) array[i * bs + off + b] = farray[i * fbs + b];
422     }
423     off += fbs;
424     PetscCall(DMSwarmRestoreField(sw, fieldnames[f], &fbs, &ftype, (void **)&farray));
425   }
426   PetscCall(VecRestoreArrayWrite(*vec, &array));
427 
428   PetscCall(PetscStrncpy(name, "DMSwarmField", PETSC_MAX_PATH_LEN));
429   for (PetscInt f = 0; f < Nf; ++f) {
430     PetscCall(PetscStrlcat(name, "_", PETSC_MAX_PATH_LEN));
431     PetscCall(PetscStrlcat(name, fieldnames[f], PETSC_MAX_PATH_LEN));
432   }
433   PetscCall(PetscObjectSetName((PetscObject)*vec, name));
434 
435   for (PetscInt f = 0; f < Nf; ++f) {
436     PetscInt fid;
437 
438     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldIdByName(swarm->db, fieldnames[f], &fid));
439     id += fid;
440   }
441   PetscCall(PetscObjectComposedDataSetInt((PetscObject)*vec, SwarmDataFieldId, id));
442 
443   PetscCall(VecSetDM(*vec, sw));
444   PetscCall(VecSetOperation(*vec, VECOP_VIEW, (void (*)(void))VecView_Swarm));
445   PetscFunctionReturn(PETSC_SUCCESS);
446 }
447 
448 /* This creates a "mass matrix" between a finite element and particle space. If a finite element interpolant is given by
449 
450      \hat f = \sum_i f_i \phi_i
451 
452    and a particle function is given by
453 
454      f = \sum_i w_i \delta(x - x_i)
455 
456    then we want to require that
457 
458      M \hat f = M_p f
459 
460    where the particle mass matrix is given by
461 
462      (M_p)_{ij} = \int \phi_i \delta(x - x_j)
463 
464    The way Dave May does particles, they amount to quadratue weights rather than delta functions, so he has |J| is in
465    his integral. We allow this with the boolean flag.
466 */
467 static PetscErrorCode DMSwarmComputeMassMatrix_Private(DM dmc, DM dmf, Mat mass, PetscBool useDeltaFunction, void *ctx)
468 {
469   const char   *name = "Mass Matrix";
470   MPI_Comm      comm;
471   DMSwarmCellDM celldm;
472   PetscDS       prob;
473   PetscSection  fsection, globalFSection;
474   PetscHSetIJ   ht;
475   PetscLayout   rLayout, colLayout;
476   PetscInt     *dnz, *onz;
477   PetscInt      locRows, locCols, rStart, colStart, colEnd, *rowIDXs;
478   PetscReal    *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0};
479   PetscScalar  *elemMat;
480   PetscInt      dim, Nf, Nfc, cStart, cEnd, totDim, maxC = 0, totNc = 0;
481   const char  **coordFields;
482   PetscReal   **coordVals;
483   PetscInt     *bs;
484 
485   PetscFunctionBegin;
486   PetscCall(PetscObjectGetComm((PetscObject)mass, &comm));
487   PetscCall(DMGetCoordinateDim(dmf, &dim));
488   PetscCall(DMGetDS(dmf, &prob));
489   PetscCall(PetscDSGetNumFields(prob, &Nf));
490   PetscCall(PetscDSGetTotalDimension(prob, &totDim));
491   PetscCall(PetscMalloc3(dim, &v0, dim * dim, &J, dim * dim, &invJ));
492   PetscCall(DMGetLocalSection(dmf, &fsection));
493   PetscCall(DMGetGlobalSection(dmf, &globalFSection));
494   PetscCall(DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd));
495   PetscCall(MatGetLocalSize(mass, &locRows, &locCols));
496 
497   PetscCall(DMSwarmGetCellDMActive(dmc, &celldm));
498   PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Nfc, &coordFields));
499   PetscCall(PetscMalloc2(Nfc, &coordVals, Nfc, &bs));
500 
501   PetscCall(PetscLayoutCreate(comm, &colLayout));
502   PetscCall(PetscLayoutSetLocalSize(colLayout, locCols));
503   PetscCall(PetscLayoutSetBlockSize(colLayout, 1));
504   PetscCall(PetscLayoutSetUp(colLayout));
505   PetscCall(PetscLayoutGetRange(colLayout, &colStart, &colEnd));
506   PetscCall(PetscLayoutDestroy(&colLayout));
507 
508   PetscCall(PetscLayoutCreate(comm, &rLayout));
509   PetscCall(PetscLayoutSetLocalSize(rLayout, locRows));
510   PetscCall(PetscLayoutSetBlockSize(rLayout, 1));
511   PetscCall(PetscLayoutSetUp(rLayout));
512   PetscCall(PetscLayoutGetRange(rLayout, &rStart, NULL));
513   PetscCall(PetscLayoutDestroy(&rLayout));
514 
515   PetscCall(PetscCalloc2(locRows, &dnz, locRows, &onz));
516   PetscCall(PetscHSetIJCreate(&ht));
517 
518   PetscCall(PetscSynchronizedFlush(comm, NULL));
519   for (PetscInt field = 0; field < Nf; ++field) {
520     PetscObject  obj;
521     PetscClassId id;
522     PetscInt     Nc;
523 
524     PetscCall(PetscDSGetDiscretization(prob, field, &obj));
525     PetscCall(PetscObjectGetClassId(obj, &id));
526     if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
527     else PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc));
528     totNc += Nc;
529   }
530   /* count non-zeros */
531   PetscCall(DMSwarmSortGetAccess(dmc));
532   for (PetscInt field = 0; field < Nf; ++field) {
533     PetscObject  obj;
534     PetscClassId id;
535     PetscInt     Nc;
536 
537     PetscCall(PetscDSGetDiscretization(prob, field, &obj));
538     PetscCall(PetscObjectGetClassId(obj, &id));
539     if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
540     else PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc));
541 
542     for (PetscInt cell = cStart; cell < cEnd; ++cell) {
543       PetscInt *findices, *cindices; /* fine is vertices, coarse is particles */
544       PetscInt  numFIndices, numCIndices;
545 
546       PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
547       PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
548       maxC = PetscMax(maxC, numCIndices);
549       {
550         PetscHashIJKey key;
551         PetscBool      missing;
552         for (PetscInt i = 0; i < numFIndices; ++i) {
553           key.j = findices[i]; /* global column (from Plex) */
554           if (key.j >= 0) {
555             /* Get indices for coarse elements */
556             for (PetscInt j = 0; j < numCIndices; ++j) {
557               for (PetscInt c = 0; c < Nc; ++c) {
558                 // TODO Need field offset on particle here
559                 key.i = cindices[j] * totNc + c + rStart; /* global cols (from Swarm) */
560                 if (key.i < 0) continue;
561                 PetscCall(PetscHSetIJQueryAdd(ht, key, &missing));
562                 if (missing) {
563                   if ((key.j >= colStart) && (key.j < colEnd)) ++dnz[key.i - rStart];
564                   else ++onz[key.i - rStart];
565                 } else SETERRQ(PetscObjectComm((PetscObject)dmf), PETSC_ERR_SUP, "Set new value at %" PetscInt_FMT ",%" PetscInt_FMT, key.i, key.j);
566               }
567             }
568           }
569         }
570         PetscCall(DMSwarmSortRestorePointsPerCell(dmc, cell, &numCIndices, &cindices));
571       }
572       PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
573     }
574   }
575   PetscCall(PetscHSetIJDestroy(&ht));
576   PetscCall(MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL));
577   PetscCall(MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
578   PetscCall(PetscFree2(dnz, onz));
579   PetscCall(PetscMalloc3(maxC * totNc * totDim, &elemMat, maxC * totNc, &rowIDXs, maxC * dim, &xi));
580   for (PetscInt field = 0; field < Nf; ++field) {
581     PetscTabulation Tcoarse;
582     PetscObject     obj;
583     PetscClassId    id;
584     PetscInt        Nc;
585 
586     PetscCall(PetscDSGetDiscretization(prob, field, &obj));
587     PetscCall(PetscObjectGetClassId(obj, &id));
588     if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
589     else PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc));
590 
591     for (PetscInt i = 0; i < Nfc; ++i) PetscCall(DMSwarmGetField(dmc, coordFields[i], &bs[i], NULL, (void **)&coordVals[i]));
592     for (PetscInt cell = cStart; cell < cEnd; ++cell) {
593       PetscInt *findices, *cindices;
594       PetscInt  numFIndices, numCIndices;
595 
596       /* TODO: Use DMField instead of assuming affine */
597       PetscCall(DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ));
598       PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
599       PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
600       for (PetscInt j = 0; j < numCIndices; ++j) {
601         PetscReal xr[8];
602         PetscInt  off = 0;
603 
604         for (PetscInt i = 0; i < Nfc; ++i) {
605           for (PetscInt b = 0; b < bs[i]; ++b, ++off) xr[off] = coordVals[i][cindices[j] * bs[i] + b];
606         }
607         PetscCheck(off == dim, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The total block size of coordinates is %" PetscInt_FMT " != %" PetscInt_FMT " the DM coordinate dimension", off, dim);
608         CoordinatesRealToRef(dim, dim, v0ref, v0, invJ, xr, &xi[j * dim]);
609       }
610       if (id == PETSCFE_CLASSID) PetscCall(PetscFECreateTabulation((PetscFE)obj, 1, numCIndices, xi, 0, &Tcoarse));
611       else PetscCall(PetscFVCreateTabulation((PetscFV)obj, 1, numCIndices, xi, 0, &Tcoarse));
612       /* Get elemMat entries by multiplying by weight */
613       PetscCall(PetscArrayzero(elemMat, numCIndices * Nc * totDim));
614       for (PetscInt i = 0; i < numFIndices / Nc; ++i) {
615         for (PetscInt j = 0; j < numCIndices; ++j) {
616           for (PetscInt c = 0; c < Nc; ++c) {
617             // TODO Need field offset on particle and field here
618             /* B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c */
619             elemMat[(j * totNc + c) * numFIndices + i * Nc + c] += Tcoarse->T[0][(j * numFIndices + i * Nc + c) * Nc + c] * (useDeltaFunction ? 1.0 : detJ);
620           }
621         }
622       }
623       for (PetscInt j = 0; j < numCIndices; ++j)
624         // TODO Need field offset on particle here
625         for (PetscInt c = 0; c < Nc; ++c) rowIDXs[j * Nc + c] = cindices[j] * totNc + c + rStart;
626       if (0) PetscCall(DMPrintCellMatrix(cell, name, numCIndices * Nc, numFIndices, elemMat));
627       PetscCall(MatSetValues(mass, numCIndices * Nc, rowIDXs, numFIndices, findices, elemMat, ADD_VALUES));
628       PetscCall(DMSwarmSortRestorePointsPerCell(dmc, cell, &numCIndices, &cindices));
629       PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
630       PetscCall(PetscTabulationDestroy(&Tcoarse));
631     }
632     for (PetscInt i = 0; i < Nfc; ++i) PetscCall(DMSwarmRestoreField(dmc, coordFields[i], &bs[i], NULL, (void **)&coordVals[i]));
633   }
634   PetscCall(PetscFree3(elemMat, rowIDXs, xi));
635   PetscCall(DMSwarmSortRestoreAccess(dmc));
636   PetscCall(PetscFree3(v0, J, invJ));
637   PetscCall(PetscFree2(coordVals, bs));
638   PetscCall(MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY));
639   PetscCall(MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY));
640   PetscFunctionReturn(PETSC_SUCCESS);
641 }
642 
643 /* Returns empty matrix for use with SNES FD */
644 static PetscErrorCode DMCreateMatrix_Swarm(DM sw, Mat *m)
645 {
646   Vec      field;
647   PetscInt size;
648 
649   PetscFunctionBegin;
650   PetscCall(DMGetGlobalVector(sw, &field));
651   PetscCall(VecGetLocalSize(field, &size));
652   PetscCall(DMRestoreGlobalVector(sw, &field));
653   PetscCall(MatCreate(PETSC_COMM_WORLD, m));
654   PetscCall(MatSetFromOptions(*m));
655   PetscCall(MatSetSizes(*m, PETSC_DECIDE, PETSC_DECIDE, size, size));
656   PetscCall(MatSeqAIJSetPreallocation(*m, 1, NULL));
657   PetscCall(MatZeroEntries(*m));
658   PetscCall(MatAssemblyBegin(*m, MAT_FINAL_ASSEMBLY));
659   PetscCall(MatAssemblyEnd(*m, MAT_FINAL_ASSEMBLY));
660   PetscCall(MatShift(*m, 1.0));
661   PetscCall(MatSetDM(*m, sw));
662   PetscFunctionReturn(PETSC_SUCCESS);
663 }
664 
665 /* FEM cols, Particle rows */
666 static PetscErrorCode DMCreateMassMatrix_Swarm(DM dmCoarse, DM dmFine, Mat *mass)
667 {
668   DMSwarmCellDM celldm;
669   PetscSection  gsf;
670   PetscInt      m, n, Np, bs;
671   void         *ctx;
672 
673   PetscFunctionBegin;
674   PetscCall(DMSwarmGetCellDMActive(dmCoarse, &celldm));
675   PetscCheck(celldm->Nf, PetscObjectComm((PetscObject)dmCoarse), PETSC_ERR_USER, "Active cell DM does not define any fields");
676   PetscCall(DMGetGlobalSection(dmFine, &gsf));
677   PetscCall(PetscSectionGetConstrainedStorageSize(gsf, &m));
678   PetscCall(DMSwarmGetLocalSize(dmCoarse, &Np));
679   PetscCall(DMSwarmCellDMGetBlockSize(celldm, dmCoarse, &bs));
680   n = Np * bs;
681   PetscCall(MatCreate(PetscObjectComm((PetscObject)dmCoarse), mass));
682   PetscCall(MatSetSizes(*mass, n, m, PETSC_DETERMINE, PETSC_DETERMINE));
683   PetscCall(MatSetType(*mass, dmCoarse->mattype));
684   PetscCall(DMGetApplicationContext(dmFine, &ctx));
685 
686   PetscCall(DMSwarmComputeMassMatrix_Private(dmCoarse, dmFine, *mass, PETSC_TRUE, ctx));
687   PetscCall(MatViewFromOptions(*mass, NULL, "-mass_mat_view"));
688   PetscFunctionReturn(PETSC_SUCCESS);
689 }
690 
691 static PetscErrorCode DMSwarmComputeMassMatrixSquare_Private(DM dmc, DM dmf, Mat mass, PetscBool useDeltaFunction, void *ctx)
692 {
693   const char   *name = "Mass Matrix Square";
694   MPI_Comm      comm;
695   DMSwarmCellDM celldm;
696   PetscDS       prob;
697   PetscSection  fsection, globalFSection;
698   PetscHSetIJ   ht;
699   PetscLayout   rLayout, colLayout;
700   PetscInt     *dnz, *onz, *adj, depth, maxConeSize, maxSupportSize, maxAdjSize;
701   PetscInt      locRows, locCols, rStart, colStart, colEnd, *rowIDXs;
702   PetscReal    *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0};
703   PetscScalar  *elemMat, *elemMatSq;
704   PetscInt      cdim, Nf, Nfc, cStart, cEnd, totDim, maxC = 0;
705   const char  **coordFields;
706   PetscReal   **coordVals;
707   PetscInt     *bs;
708 
709   PetscFunctionBegin;
710   PetscCall(PetscObjectGetComm((PetscObject)mass, &comm));
711   PetscCall(DMGetCoordinateDim(dmf, &cdim));
712   PetscCall(DMGetDS(dmf, &prob));
713   PetscCall(PetscDSGetNumFields(prob, &Nf));
714   PetscCall(PetscDSGetTotalDimension(prob, &totDim));
715   PetscCall(PetscMalloc3(cdim, &v0, cdim * cdim, &J, cdim * cdim, &invJ));
716   PetscCall(DMGetLocalSection(dmf, &fsection));
717   PetscCall(DMGetGlobalSection(dmf, &globalFSection));
718   PetscCall(DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd));
719   PetscCall(MatGetLocalSize(mass, &locRows, &locCols));
720 
721   PetscCall(DMSwarmGetCellDMActive(dmc, &celldm));
722   PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Nfc, &coordFields));
723   PetscCall(PetscMalloc2(Nfc, &coordVals, Nfc, &bs));
724 
725   PetscCall(PetscLayoutCreate(comm, &colLayout));
726   PetscCall(PetscLayoutSetLocalSize(colLayout, locCols));
727   PetscCall(PetscLayoutSetBlockSize(colLayout, 1));
728   PetscCall(PetscLayoutSetUp(colLayout));
729   PetscCall(PetscLayoutGetRange(colLayout, &colStart, &colEnd));
730   PetscCall(PetscLayoutDestroy(&colLayout));
731 
732   PetscCall(PetscLayoutCreate(comm, &rLayout));
733   PetscCall(PetscLayoutSetLocalSize(rLayout, locRows));
734   PetscCall(PetscLayoutSetBlockSize(rLayout, 1));
735   PetscCall(PetscLayoutSetUp(rLayout));
736   PetscCall(PetscLayoutGetRange(rLayout, &rStart, NULL));
737   PetscCall(PetscLayoutDestroy(&rLayout));
738 
739   PetscCall(DMPlexGetDepth(dmf, &depth));
740   PetscCall(DMPlexGetMaxSizes(dmf, &maxConeSize, &maxSupportSize));
741   maxAdjSize = PetscPowInt(maxConeSize * maxSupportSize, depth);
742   PetscCall(PetscMalloc1(maxAdjSize, &adj));
743 
744   PetscCall(PetscCalloc2(locRows, &dnz, locRows, &onz));
745   PetscCall(PetscHSetIJCreate(&ht));
746   /* Count nonzeros
747        This is just FVM++, but we cannot use the Plex P0 allocation since unknowns in a cell will not be contiguous
748   */
749   PetscCall(DMSwarmSortGetAccess(dmc));
750   for (PetscInt cell = cStart; cell < cEnd; ++cell) {
751     PetscInt *cindices;
752     PetscInt  numCIndices;
753 #if 0
754     PetscInt  adjSize = maxAdjSize, a, j;
755 #endif
756 
757     PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
758     maxC = PetscMax(maxC, numCIndices);
759     /* Diagonal block */
760     for (PetscInt i = 0; i < numCIndices; ++i) dnz[cindices[i]] += numCIndices;
761 #if 0
762     /* Off-diagonal blocks */
763     PetscCall(DMPlexGetAdjacency(dmf, cell, &adjSize, &adj));
764     for (a = 0; a < adjSize; ++a) {
765       if (adj[a] >= cStart && adj[a] < cEnd && adj[a] != cell) {
766         const PetscInt ncell = adj[a];
767         PetscInt      *ncindices;
768         PetscInt       numNCIndices;
769 
770         PetscCall(DMSwarmSortGetPointsPerCell(dmc, ncell, &numNCIndices, &ncindices));
771         {
772           PetscHashIJKey key;
773           PetscBool      missing;
774 
775           for (i = 0; i < numCIndices; ++i) {
776             key.i = cindices[i] + rStart; /* global rows (from Swarm) */
777             if (key.i < 0) continue;
778             for (j = 0; j < numNCIndices; ++j) {
779               key.j = ncindices[j] + rStart; /* global column (from Swarm) */
780               if (key.j < 0) continue;
781               PetscCall(PetscHSetIJQueryAdd(ht, key, &missing));
782               if (missing) {
783                 if ((key.j >= colStart) && (key.j < colEnd)) ++dnz[key.i - rStart];
784                 else                                         ++onz[key.i - rStart];
785               }
786             }
787           }
788         }
789         PetscCall(DMSwarmSortRestorePointsPerCell(dmc, ncell, &numNCIndices, &ncindices));
790       }
791     }
792 #endif
793     PetscCall(DMSwarmSortRestorePointsPerCell(dmc, cell, &numCIndices, &cindices));
794   }
795   PetscCall(PetscHSetIJDestroy(&ht));
796   PetscCall(MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL));
797   PetscCall(MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
798   PetscCall(PetscFree2(dnz, onz));
799   PetscCall(PetscMalloc4(maxC * totDim, &elemMat, maxC * maxC, &elemMatSq, maxC, &rowIDXs, maxC * cdim, &xi));
800   /* Fill in values
801        Each entry is a sum of terms \phi_i(x_p) \phi_i(x_q)
802        Start just by producing block diagonal
803        Could loop over adjacent cells
804          Produce neighboring element matrix
805          TODO Determine which columns and rows correspond to shared dual vector
806          Do MatMatMult with rectangular matrices
807          Insert block
808   */
809   for (PetscInt field = 0; field < Nf; ++field) {
810     PetscTabulation Tcoarse;
811     PetscObject     obj;
812     PetscInt        Nc;
813 
814     PetscCall(PetscDSGetDiscretization(prob, field, &obj));
815     PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
816     PetscCheck(Nc == 1, PetscObjectComm((PetscObject)dmf), PETSC_ERR_SUP, "Can only interpolate a scalar field from particles, Nc = %" PetscInt_FMT, Nc);
817     for (PetscInt i = 0; i < Nfc; ++i) PetscCall(DMSwarmGetField(dmc, coordFields[i], &bs[i], NULL, (void **)&coordVals[i]));
818     for (PetscInt cell = cStart; cell < cEnd; ++cell) {
819       PetscInt *findices, *cindices;
820       PetscInt  numFIndices, numCIndices;
821 
822       /* TODO: Use DMField instead of assuming affine */
823       PetscCall(DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ));
824       PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
825       PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices));
826       for (PetscInt p = 0; p < numCIndices; ++p) {
827         PetscReal xr[8];
828         PetscInt  off = 0;
829 
830         for (PetscInt i = 0; i < Nfc; ++i) {
831           for (PetscInt b = 0; b < bs[i]; ++b, ++off) xr[off] = coordVals[i][cindices[p] * bs[i] + b];
832         }
833         PetscCheck(off == cdim, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The total block size of coordinates is %" PetscInt_FMT " != %" PetscInt_FMT " the DM coordinate dimension", off, cdim);
834         CoordinatesRealToRef(cdim, cdim, v0ref, v0, invJ, xr, &xi[p * cdim]);
835       }
836       PetscCall(PetscFECreateTabulation((PetscFE)obj, 1, numCIndices, xi, 0, &Tcoarse));
837       /* Get elemMat entries by multiplying by weight */
838       PetscCall(PetscArrayzero(elemMat, numCIndices * totDim));
839       for (PetscInt i = 0; i < numFIndices; ++i) {
840         for (PetscInt p = 0; p < numCIndices; ++p) {
841           for (PetscInt c = 0; c < Nc; ++c) {
842             /* B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c */
843             elemMat[p * numFIndices + i] += Tcoarse->T[0][(p * numFIndices + i) * Nc + c] * (useDeltaFunction ? 1.0 : detJ);
844           }
845         }
846       }
847       PetscCall(PetscTabulationDestroy(&Tcoarse));
848       for (PetscInt p = 0; p < numCIndices; ++p) rowIDXs[p] = cindices[p] + rStart;
849       if (0) PetscCall(DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat));
850       /* Block diagonal */
851       if (numCIndices) {
852         PetscBLASInt blasn, blask;
853         PetscScalar  one = 1.0, zero = 0.0;
854 
855         PetscCall(PetscBLASIntCast(numCIndices, &blasn));
856         PetscCall(PetscBLASIntCast(numFIndices, &blask));
857         PetscCallBLAS("BLASgemm", BLASgemm_("T", "N", &blasn, &blasn, &blask, &one, elemMat, &blask, elemMat, &blask, &zero, elemMatSq, &blasn));
858       }
859       PetscCall(MatSetValues(mass, numCIndices, rowIDXs, numCIndices, rowIDXs, elemMatSq, ADD_VALUES));
860       /* TODO off-diagonal */
861       PetscCall(DMSwarmSortRestorePointsPerCell(dmc, cell, &numCIndices, &cindices));
862       PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
863     }
864     for (PetscInt i = 0; i < Nfc; ++i) PetscCall(DMSwarmRestoreField(dmc, coordFields[i], &bs[i], NULL, (void **)&coordVals[i]));
865   }
866   PetscCall(PetscFree4(elemMat, elemMatSq, rowIDXs, xi));
867   PetscCall(PetscFree(adj));
868   PetscCall(DMSwarmSortRestoreAccess(dmc));
869   PetscCall(PetscFree3(v0, J, invJ));
870   PetscCall(PetscFree2(coordVals, bs));
871   PetscCall(MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY));
872   PetscCall(MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY));
873   PetscFunctionReturn(PETSC_SUCCESS);
874 }
875 
876 /*@
877   DMSwarmCreateMassMatrixSquare - Creates the block-diagonal of the square, M^T_p M_p, of the particle mass matrix M_p
878 
879   Collective
880 
881   Input Parameters:
882 + dmCoarse - a `DMSWARM`
883 - dmFine   - a `DMPLEX`
884 
885   Output Parameter:
886 . mass - the square of the particle mass matrix
887 
888   Level: advanced
889 
890   Note:
891   We only compute the block diagonal since this provides a good preconditioner and is completely local. It would be possible in the
892   future to compute the full normal equations.
893 
894 .seealso: `DM`, `DMSWARM`, `DMCreateMassMatrix()`
895 @*/
896 PetscErrorCode DMSwarmCreateMassMatrixSquare(DM dmCoarse, DM dmFine, Mat *mass)
897 {
898   PetscInt n;
899   void    *ctx;
900 
901   PetscFunctionBegin;
902   PetscCall(DMSwarmGetLocalSize(dmCoarse, &n));
903   PetscCall(MatCreate(PetscObjectComm((PetscObject)dmCoarse), mass));
904   PetscCall(MatSetSizes(*mass, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
905   PetscCall(MatSetType(*mass, dmCoarse->mattype));
906   PetscCall(DMGetApplicationContext(dmFine, &ctx));
907 
908   PetscCall(DMSwarmComputeMassMatrixSquare_Private(dmCoarse, dmFine, *mass, PETSC_TRUE, ctx));
909   PetscCall(MatViewFromOptions(*mass, NULL, "-mass_sq_mat_view"));
910   PetscFunctionReturn(PETSC_SUCCESS);
911 }
912 
913 /* This creates a "gradient matrix" between a finite element and particle space, which is meant to enforce a weak divergence condition on the particle space. We are looking for a finite element field that has the same divergence as our particle field, so that
914 
915      \int_X \psi_i \nabla \cdot \hat f = \int_X \psi_i \nabla \cdot f
916 
917    and then integrate by parts
918 
919      \int_X \nabla \psi_i \cdot \hat f = \int_X \nabla \psi_i \cdot f
920 
921    where \psi is from a scalar FE space. If a finite element interpolant is given by
922 
923      \hat f^c = \sum_i f_i \phi^c_i
924 
925    and a particle function is given by
926 
927      f^c = \sum_p f^c_p \delta(x - x_p)
928 
929    then we want to require that
930 
931      D_f \hat f = D_p f
932 
933    where the gradient matrices are given by
934 
935      (D_f)_{i(jc)} = \int \partial_c \psi_i \phi_j
936      (D_p)_{i(jc)} = \int \partial_c \psi_i \delta(x - x_j)
937 
938    Thus we need two finite element spaces, a scalar and a vector. The vector space holds the representer for the
939    vector particle field. The scalar space holds the output of D_p or D_f, which is the weak divergence of the field.
940 
941    The way Dave May does particles, they amount to quadratue weights rather than delta functions, so he has |J| is in
942    his integral. We allow this with the boolean flag.
943 */
944 static PetscErrorCode DMSwarmComputeGradientMatrix_Private(DM sw, DM dm, Mat derv, PetscBool useDeltaFunction, void *ctx)
945 {
946   const char   *name = "Derivative Matrix";
947   MPI_Comm      comm;
948   DMSwarmCellDM celldm;
949   PetscDS       ds;
950   PetscSection  fsection, globalFSection;
951   PetscLayout   rLayout;
952   PetscInt      locRows, rStart, *rowIDXs;
953   PetscReal    *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0};
954   PetscScalar  *elemMat;
955   PetscInt      cdim, Nf, Nfc, cStart, cEnd, totDim, maxNpc = 0, totNc = 0;
956   const char  **coordFields;
957   PetscReal   **coordVals;
958   PetscInt     *bs;
959 
960   PetscFunctionBegin;
961   PetscCall(PetscObjectGetComm((PetscObject)derv, &comm));
962   PetscCall(DMGetCoordinateDim(dm, &cdim));
963   PetscCall(DMGetDS(dm, &ds));
964   PetscCall(PetscDSGetNumFields(ds, &Nf));
965   PetscCheck(Nf == 1, comm, PETSC_ERR_SUP, "Currently, we only support a single field");
966   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
967   PetscCall(PetscMalloc3(cdim, &v0, cdim * cdim, &J, cdim * cdim, &invJ));
968   PetscCall(DMGetLocalSection(dm, &fsection));
969   PetscCall(DMGetGlobalSection(dm, &globalFSection));
970   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
971   PetscCall(MatGetLocalSize(derv, &locRows, NULL));
972 
973   PetscCall(DMSwarmGetCellDMActive(sw, &celldm));
974   PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Nfc, &coordFields));
975   PetscCheck(Nfc == 1, comm, PETSC_ERR_SUP, "Currently, we only support a single field");
976   PetscCall(PetscMalloc2(Nfc, &coordVals, Nfc, &bs));
977 
978   PetscCall(PetscLayoutCreate(comm, &rLayout));
979   PetscCall(PetscLayoutSetLocalSize(rLayout, locRows));
980   PetscCall(PetscLayoutSetBlockSize(rLayout, cdim));
981   PetscCall(PetscLayoutSetUp(rLayout));
982   PetscCall(PetscLayoutGetRange(rLayout, &rStart, NULL));
983   PetscCall(PetscLayoutDestroy(&rLayout));
984 
985   for (PetscInt field = 0; field < Nf; ++field) {
986     PetscObject obj;
987     PetscInt    Nc;
988 
989     PetscCall(PetscDSGetDiscretization(ds, field, &obj));
990     PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
991     totNc += Nc;
992   }
993   PetscCheck(totNc == 1, comm, PETSC_ERR_ARG_WRONG, "The number of field components %" PetscInt_FMT " != 1", totNc);
994   /* count non-zeros */
995   PetscCall(DMSwarmSortGetAccess(sw));
996   for (PetscInt field = 0; field < Nf; ++field) {
997     PetscObject obj;
998     PetscInt    Nc;
999 
1000     PetscCall(PetscDSGetDiscretization(ds, field, &obj));
1001     PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
1002     for (PetscInt cell = cStart; cell < cEnd; ++cell) {
1003       PetscInt *pind;
1004       PetscInt  Npc;
1005 
1006       PetscCall(DMSwarmSortGetPointsPerCell(sw, cell, &Npc, &pind));
1007       maxNpc = PetscMax(maxNpc, Npc);
1008       PetscCall(DMSwarmSortRestorePointsPerCell(sw, cell, &Npc, &pind));
1009     }
1010   }
1011   PetscCall(PetscMalloc3(maxNpc * cdim * totDim, &elemMat, maxNpc * cdim, &rowIDXs, maxNpc * cdim, &xi));
1012   for (PetscInt field = 0; field < Nf; ++field) {
1013     PetscTabulation Tcoarse;
1014     PetscFE         fe;
1015 
1016     PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
1017     for (PetscInt i = 0; i < Nfc; ++i) PetscCall(DMSwarmGetField(sw, coordFields[i], &bs[i], NULL, (void **)&coordVals[i]));
1018     for (PetscInt cell = cStart; cell < cEnd; ++cell) {
1019       PetscInt *findices, *pind;
1020       PetscInt  numFIndices, Npc;
1021 
1022       /* TODO: Use DMField instead of assuming affine */
1023       PetscCall(DMPlexComputeCellGeometryFEM(dm, cell, NULL, v0, J, invJ, &detJ));
1024       PetscCall(DMPlexGetClosureIndices(dm, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
1025       PetscCall(DMSwarmSortGetPointsPerCell(sw, cell, &Npc, &pind));
1026       for (PetscInt j = 0; j < Npc; ++j) {
1027         PetscReal xr[8];
1028         PetscInt  off = 0;
1029 
1030         for (PetscInt i = 0; i < Nfc; ++i) {
1031           for (PetscInt b = 0; b < bs[i]; ++b, ++off) xr[off] = coordVals[i][pind[j] * bs[i] + b];
1032         }
1033         PetscCheck(off == cdim, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The total block size of coordinates is %" PetscInt_FMT " != %" PetscInt_FMT " the DM coordinate dimension", off, cdim);
1034         CoordinatesRealToRef(cdim, cdim, v0ref, v0, invJ, xr, &xi[j * cdim]);
1035       }
1036       PetscCall(PetscFECreateTabulation(fe, 1, Npc, xi, 1, &Tcoarse));
1037       /* Get elemMat entries by multiplying by weight */
1038       PetscCall(PetscArrayzero(elemMat, Npc * cdim * totDim));
1039       for (PetscInt i = 0; i < numFIndices; ++i) {
1040         for (PetscInt j = 0; j < Npc; ++j) {
1041           /* D[((p*pdim + i)*Nc + c)*cdim + d] is the value at point p for basis function i, component c, derviative d */
1042           for (PetscInt d = 0; d < cdim; ++d) {
1043             xi[d] = 0.;
1044             for (PetscInt e = 0; e < cdim; ++e) xi[d] += invJ[e * cdim + d] * Tcoarse->T[1][(j * numFIndices + i) * cdim + e];
1045             elemMat[(j * cdim + d) * numFIndices + i] += xi[d] * (useDeltaFunction ? 1.0 : detJ);
1046           }
1047         }
1048       }
1049       for (PetscInt j = 0; j < Npc; ++j)
1050         for (PetscInt d = 0; d < cdim; ++d) rowIDXs[j * cdim + d] = pind[j] * cdim + d + rStart;
1051       if (0) PetscCall(DMPrintCellMatrix(cell, name, Npc * cdim, numFIndices, elemMat));
1052       PetscCall(MatSetValues(derv, Npc * cdim, rowIDXs, numFIndices, findices, elemMat, ADD_VALUES));
1053       PetscCall(DMSwarmSortRestorePointsPerCell(sw, cell, &Npc, &pind));
1054       PetscCall(DMPlexRestoreClosureIndices(dm, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
1055       PetscCall(PetscTabulationDestroy(&Tcoarse));
1056     }
1057     for (PetscInt i = 0; i < Nfc; ++i) PetscCall(DMSwarmRestoreField(sw, coordFields[i], &bs[i], NULL, (void **)&coordVals[i]));
1058   }
1059   PetscCall(PetscFree3(elemMat, rowIDXs, xi));
1060   PetscCall(DMSwarmSortRestoreAccess(sw));
1061   PetscCall(PetscFree3(v0, J, invJ));
1062   PetscCall(PetscFree2(coordVals, bs));
1063   PetscCall(MatAssemblyBegin(derv, MAT_FINAL_ASSEMBLY));
1064   PetscCall(MatAssemblyEnd(derv, MAT_FINAL_ASSEMBLY));
1065   PetscFunctionReturn(PETSC_SUCCESS);
1066 }
1067 
1068 /* FEM cols:      this is a scalar space
1069    Particle rows: this is a vector space that contracts with the derivative
1070 */
1071 static PetscErrorCode DMCreateGradientMatrix_Swarm(DM sw, DM dm, Mat *derv)
1072 {
1073   DMSwarmCellDM celldm;
1074   PetscSection  gs;
1075   PetscInt      cdim, m, n, Np, bs;
1076   void         *ctx;
1077   MPI_Comm      comm;
1078 
1079   PetscFunctionBegin;
1080   PetscCall(PetscObjectGetComm((PetscObject)sw, &comm));
1081   PetscCall(DMGetCoordinateDim(dm, &cdim));
1082   PetscCall(DMSwarmGetCellDMActive(sw, &celldm));
1083   PetscCheck(celldm->Nf, comm, PETSC_ERR_USER, "Active cell DM does not define any fields");
1084   PetscCall(DMGetGlobalSection(dm, &gs));
1085   PetscCall(PetscSectionGetConstrainedStorageSize(gs, &n));
1086   PetscCall(DMSwarmGetLocalSize(sw, &Np));
1087   PetscCall(DMSwarmCellDMGetBlockSize(celldm, sw, &bs));
1088   PetscCheck(cdim == bs, comm, PETSC_ERR_ARG_WRONG, "Coordinate dimension %" PetscInt_FMT " != %" PetscInt_FMT " swarm field block size", cdim, bs);
1089   m = Np * bs;
1090   PetscCall(MatCreate(PetscObjectComm((PetscObject)sw), derv));
1091   PetscCall(PetscObjectSetName((PetscObject)*derv, "Swarm Derivative Matrix"));
1092   PetscCall(MatSetSizes(*derv, m, n, PETSC_DETERMINE, PETSC_DETERMINE));
1093   PetscCall(MatSetType(*derv, sw->mattype));
1094   PetscCall(DMGetApplicationContext(dm, &ctx));
1095 
1096   PetscCall(DMSwarmComputeGradientMatrix_Private(sw, dm, *derv, PETSC_TRUE, ctx));
1097   PetscCall(MatViewFromOptions(*derv, NULL, "-gradient_mat_view"));
1098   PetscFunctionReturn(PETSC_SUCCESS);
1099 }
1100 
1101 /*@
1102   DMSwarmCreateGlobalVectorFromField - Creates a `Vec` object sharing the array associated with a given field
1103 
1104   Collective
1105 
1106   Input Parameters:
1107 + dm        - a `DMSWARM`
1108 - fieldname - the textual name given to a registered field
1109 
1110   Output Parameter:
1111 . vec - the vector
1112 
1113   Level: beginner
1114 
1115   Note:
1116   The vector must be returned using a matching call to `DMSwarmDestroyGlobalVectorFromField()`.
1117 
1118 .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmDestroyGlobalVectorFromField()`
1119 @*/
1120 PetscErrorCode DMSwarmCreateGlobalVectorFromField(DM dm, const char fieldname[], Vec *vec)
1121 {
1122   MPI_Comm comm = PetscObjectComm((PetscObject)dm);
1123 
1124   PetscFunctionBegin;
1125   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
1126   PetscCall(DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec));
1127   PetscFunctionReturn(PETSC_SUCCESS);
1128 }
1129 
1130 /*@
1131   DMSwarmDestroyGlobalVectorFromField - Destroys the `Vec` object which share the array associated with a given field
1132 
1133   Collective
1134 
1135   Input Parameters:
1136 + dm        - a `DMSWARM`
1137 - fieldname - the textual name given to a registered field
1138 
1139   Output Parameter:
1140 . vec - the vector
1141 
1142   Level: beginner
1143 
1144 .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmCreateGlobalVectorFromField()`
1145 @*/
1146 PetscErrorCode DMSwarmDestroyGlobalVectorFromField(DM dm, const char fieldname[], Vec *vec)
1147 {
1148   PetscFunctionBegin;
1149   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
1150   PetscCall(DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec));
1151   PetscFunctionReturn(PETSC_SUCCESS);
1152 }
1153 
1154 /*@
1155   DMSwarmCreateLocalVectorFromField - Creates a `Vec` object sharing the array associated with a given field
1156 
1157   Collective
1158 
1159   Input Parameters:
1160 + dm        - a `DMSWARM`
1161 - fieldname - the textual name given to a registered field
1162 
1163   Output Parameter:
1164 . vec - the vector
1165 
1166   Level: beginner
1167 
1168   Note:
1169   The vector must be returned using a matching call to DMSwarmDestroyLocalVectorFromField().
1170 
1171 .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmDestroyLocalVectorFromField()`
1172 @*/
1173 PetscErrorCode DMSwarmCreateLocalVectorFromField(DM dm, const char fieldname[], Vec *vec)
1174 {
1175   MPI_Comm comm = PETSC_COMM_SELF;
1176 
1177   PetscFunctionBegin;
1178   PetscCall(DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec));
1179   PetscFunctionReturn(PETSC_SUCCESS);
1180 }
1181 
1182 /*@
1183   DMSwarmDestroyLocalVectorFromField - Destroys the `Vec` object which share the array associated with a given field
1184 
1185   Collective
1186 
1187   Input Parameters:
1188 + dm        - a `DMSWARM`
1189 - fieldname - the textual name given to a registered field
1190 
1191   Output Parameter:
1192 . vec - the vector
1193 
1194   Level: beginner
1195 
1196 .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmCreateLocalVectorFromField()`
1197 @*/
1198 PetscErrorCode DMSwarmDestroyLocalVectorFromField(DM dm, const char fieldname[], Vec *vec)
1199 {
1200   PetscFunctionBegin;
1201   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
1202   PetscCall(DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec));
1203   PetscFunctionReturn(PETSC_SUCCESS);
1204 }
1205 
1206 /*@
1207   DMSwarmCreateGlobalVectorFromFields - Creates a `Vec` object sharing the array associated with a given field set
1208 
1209   Collective
1210 
1211   Input Parameters:
1212 + dm         - a `DMSWARM`
1213 . Nf         - the number of fields
1214 - fieldnames - the textual names given to the registered fields
1215 
1216   Output Parameter:
1217 . vec - the vector
1218 
1219   Level: beginner
1220 
1221   Notes:
1222   The vector must be returned using a matching call to `DMSwarmDestroyGlobalVectorFromFields()`.
1223 
1224   This vector is copyin-copyout, rather than a direct pointer like `DMSwarmCreateGlobalVectorFromField()`
1225 
1226 .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmDestroyGlobalVectorFromFields()`
1227 @*/
1228 PetscErrorCode DMSwarmCreateGlobalVectorFromFields(DM dm, PetscInt Nf, const char *fieldnames[], Vec *vec)
1229 {
1230   MPI_Comm comm = PetscObjectComm((PetscObject)dm);
1231 
1232   PetscFunctionBegin;
1233   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
1234   PetscCall(DMSwarmCreateVectorFromFields_Private(dm, Nf, fieldnames, comm, vec));
1235   PetscFunctionReturn(PETSC_SUCCESS);
1236 }
1237 
1238 /*@
1239   DMSwarmDestroyGlobalVectorFromFields - Destroys the `Vec` object which share the array associated with a given field set
1240 
1241   Collective
1242 
1243   Input Parameters:
1244 + dm         - a `DMSWARM`
1245 . Nf         - the number of fields
1246 - fieldnames - the textual names given to the registered fields
1247 
1248   Output Parameter:
1249 . vec - the vector
1250 
1251   Level: beginner
1252 
1253 .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmCreateGlobalVectorFromField()`
1254 @*/
1255 PetscErrorCode DMSwarmDestroyGlobalVectorFromFields(DM dm, PetscInt Nf, const char *fieldnames[], Vec *vec)
1256 {
1257   PetscFunctionBegin;
1258   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
1259   PetscCall(DMSwarmDestroyVectorFromFields_Private(dm, Nf, fieldnames, vec));
1260   PetscFunctionReturn(PETSC_SUCCESS);
1261 }
1262 
1263 /*@
1264   DMSwarmCreateLocalVectorFromFields - Creates a `Vec` object sharing the array associated with a given field set
1265 
1266   Collective
1267 
1268   Input Parameters:
1269 + dm         - a `DMSWARM`
1270 . Nf         - the number of fields
1271 - fieldnames - the textual names given to the registered fields
1272 
1273   Output Parameter:
1274 . vec - the vector
1275 
1276   Level: beginner
1277 
1278   Notes:
1279   The vector must be returned using a matching call to DMSwarmDestroyLocalVectorFromField().
1280 
1281   This vector is copyin-copyout, rather than a direct pointer like `DMSwarmCreateLocalVectorFromField()`
1282 
1283 .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmDestroyLocalVectorFromField()`
1284 @*/
1285 PetscErrorCode DMSwarmCreateLocalVectorFromFields(DM dm, PetscInt Nf, const char *fieldnames[], Vec *vec)
1286 {
1287   MPI_Comm comm = PETSC_COMM_SELF;
1288 
1289   PetscFunctionBegin;
1290   PetscCall(DMSwarmCreateVectorFromFields_Private(dm, Nf, fieldnames, comm, vec));
1291   PetscFunctionReturn(PETSC_SUCCESS);
1292 }
1293 
1294 /*@
1295   DMSwarmDestroyLocalVectorFromFields - Destroys the `Vec` object which share the array associated with a given field set
1296 
1297   Collective
1298 
1299   Input Parameters:
1300 + dm         - a `DMSWARM`
1301 . Nf         - the number of fields
1302 - fieldnames - the textual names given to the registered fields
1303 
1304   Output Parameter:
1305 . vec - the vector
1306 
1307   Level: beginner
1308 
1309 .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmCreateLocalVectorFromFields()`
1310 @*/
1311 PetscErrorCode DMSwarmDestroyLocalVectorFromFields(DM dm, PetscInt Nf, const char *fieldnames[], Vec *vec)
1312 {
1313   PetscFunctionBegin;
1314   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
1315   PetscCall(DMSwarmDestroyVectorFromFields_Private(dm, Nf, fieldnames, vec));
1316   PetscFunctionReturn(PETSC_SUCCESS);
1317 }
1318 
1319 /*@
1320   DMSwarmInitializeFieldRegister - Initiates the registration of fields to a `DMSWARM`
1321 
1322   Collective
1323 
1324   Input Parameter:
1325 . dm - a `DMSWARM`
1326 
1327   Level: beginner
1328 
1329   Note:
1330   After all fields have been registered, you must call `DMSwarmFinalizeFieldRegister()`.
1331 
1332 .seealso: `DM`, `DMSWARM`, `DMSwarmFinalizeFieldRegister()`, `DMSwarmRegisterPetscDatatypeField()`,
1333           `DMSwarmRegisterUserStructField()`, `DMSwarmRegisterUserDatatypeField()`
1334 @*/
1335 PetscErrorCode DMSwarmInitializeFieldRegister(DM dm)
1336 {
1337   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1338 
1339   PetscFunctionBegin;
1340   if (!swarm->field_registration_initialized) {
1341     swarm->field_registration_initialized = PETSC_TRUE;
1342     PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmField_pid, 1, PETSC_INT64)); /* unique identifier */
1343     PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmField_rank, 1, PETSC_INT));  /* used for communication */
1344   }
1345   PetscFunctionReturn(PETSC_SUCCESS);
1346 }
1347 
1348 /*@
1349   DMSwarmFinalizeFieldRegister - Finalizes the registration of fields to a `DMSWARM`
1350 
1351   Collective
1352 
1353   Input Parameter:
1354 . dm - a `DMSWARM`
1355 
1356   Level: beginner
1357 
1358   Note:
1359   After `DMSwarmFinalizeFieldRegister()` has been called, no new fields can be defined on the `DMSWARM`.
1360 
1361 .seealso: `DM`, `DMSWARM`, `DMSwarmInitializeFieldRegister()`, `DMSwarmRegisterPetscDatatypeField()`,
1362           `DMSwarmRegisterUserStructField()`, `DMSwarmRegisterUserDatatypeField()`
1363 @*/
1364 PetscErrorCode DMSwarmFinalizeFieldRegister(DM dm)
1365 {
1366   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1367 
1368   PetscFunctionBegin;
1369   if (!swarm->field_registration_finalized) PetscCall(DMSwarmDataBucketFinalize(swarm->db));
1370   swarm->field_registration_finalized = PETSC_TRUE;
1371   PetscFunctionReturn(PETSC_SUCCESS);
1372 }
1373 
1374 /*@
1375   DMSwarmSetLocalSizes - Sets the length of all registered fields on the `DMSWARM`
1376 
1377   Not Collective
1378 
1379   Input Parameters:
1380 + sw     - a `DMSWARM`
1381 . nlocal - the length of each registered field
1382 - buffer - the length of the buffer used to efficient dynamic re-sizing
1383 
1384   Level: beginner
1385 
1386 .seealso: `DM`, `DMSWARM`, `DMSwarmGetLocalSize()`
1387 @*/
1388 PetscErrorCode DMSwarmSetLocalSizes(DM sw, PetscInt nlocal, PetscInt buffer)
1389 {
1390   DM_Swarm   *swarm = (DM_Swarm *)sw->data;
1391   PetscMPIInt rank;
1392   PetscInt   *rankval;
1393 
1394   PetscFunctionBegin;
1395   PetscCall(PetscLogEventBegin(DMSWARM_SetSizes, 0, 0, 0, 0));
1396   PetscCall(DMSwarmDataBucketSetSizes(swarm->db, nlocal, buffer));
1397   PetscCall(PetscLogEventEnd(DMSWARM_SetSizes, 0, 0, 0, 0));
1398 
1399   // Initialize values in pid and rank placeholders
1400   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sw), &rank));
1401   PetscCall(DMSwarmGetField(sw, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
1402   for (PetscInt p = 0; p < nlocal; p++) rankval[p] = rank;
1403   PetscCall(DMSwarmRestoreField(sw, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
1404   /* TODO: [pid - use MPI_Scan] */
1405   PetscFunctionReturn(PETSC_SUCCESS);
1406 }
1407 
1408 /*@
1409   DMSwarmSetCellDM - Attaches a `DM` to a `DMSWARM`
1410 
1411   Collective
1412 
1413   Input Parameters:
1414 + sw - a `DMSWARM`
1415 - dm - the `DM` to attach to the `DMSWARM`
1416 
1417   Level: beginner
1418 
1419   Note:
1420   The attached `DM` (dm) will be queried for point location and
1421   neighbor MPI-rank information if `DMSwarmMigrate()` is called.
1422 
1423 .seealso: `DM`, `DMSWARM`, `DMSwarmSetType()`, `DMSwarmGetCellDM()`, `DMSwarmMigrate()`
1424 @*/
1425 PetscErrorCode DMSwarmSetCellDM(DM sw, DM dm)
1426 {
1427   DMSwarmCellDM celldm;
1428   const char   *name;
1429   char         *coordName;
1430 
1431   PetscFunctionBegin;
1432   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
1433   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1434   PetscCall(PetscStrallocpy(DMSwarmPICField_coor, &coordName));
1435   PetscCall(DMSwarmCellDMCreate(dm, 0, NULL, 1, (const char **)&coordName, &celldm));
1436   PetscCall(PetscFree(coordName));
1437   PetscCall(PetscObjectGetName((PetscObject)celldm, &name));
1438   PetscCall(DMSwarmAddCellDM(sw, celldm));
1439   PetscCall(DMSwarmCellDMDestroy(&celldm));
1440   PetscCall(DMSwarmSetCellDMActive(sw, name));
1441   PetscFunctionReturn(PETSC_SUCCESS);
1442 }
1443 
1444 /*@
1445   DMSwarmGetCellDM - Fetches the active cell `DM`
1446 
1447   Collective
1448 
1449   Input Parameter:
1450 . sw - a `DMSWARM`
1451 
1452   Output Parameter:
1453 . dm - the active `DM` for the `DMSWARM`
1454 
1455   Level: beginner
1456 
1457 .seealso: `DM`, `DMSWARM`, `DMSwarmSetCellDM()`
1458 @*/
1459 PetscErrorCode DMSwarmGetCellDM(DM sw, DM *dm)
1460 {
1461   DM_Swarm     *swarm = (DM_Swarm *)sw->data;
1462   DMSwarmCellDM celldm;
1463 
1464   PetscFunctionBegin;
1465   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
1466   PetscCall(PetscObjectListFind(swarm->cellDMs, swarm->activeCellDM, (PetscObject *)&celldm));
1467   PetscCheck(celldm, PetscObjectComm((PetscObject)sw), PETSC_ERR_ARG_WRONG, "There is no cell DM named %s in this Swarm", swarm->activeCellDM);
1468   PetscCall(DMSwarmCellDMGetDM(celldm, dm));
1469   PetscFunctionReturn(PETSC_SUCCESS);
1470 }
1471 
1472 /*@C
1473   DMSwarmGetCellDMNames - Get the list of cell `DM` names
1474 
1475   Not collective
1476 
1477   Input Parameter:
1478 . sw - a `DMSWARM`
1479 
1480   Output Parameters:
1481 + Ndm     - the number of `DMSwarmCellDM` in the `DMSWARM`
1482 - celldms - the name of each `DMSwarmCellDM`
1483 
1484   Level: beginner
1485 
1486 .seealso: `DM`, `DMSWARM`, `DMSwarmSetCellDM()`, `DMSwarmGetCellDMByName()`
1487 @*/
1488 PetscErrorCode DMSwarmGetCellDMNames(DM sw, PetscInt *Ndm, const char **celldms[])
1489 {
1490   DM_Swarm       *swarm = (DM_Swarm *)sw->data;
1491   PetscObjectList next  = swarm->cellDMs;
1492   PetscInt        n     = 0;
1493 
1494   PetscFunctionBegin;
1495   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
1496   PetscAssertPointer(Ndm, 2);
1497   PetscAssertPointer(celldms, 3);
1498   while (next) {
1499     next = next->next;
1500     ++n;
1501   }
1502   PetscCall(PetscMalloc1(n, celldms));
1503   next = swarm->cellDMs;
1504   n    = 0;
1505   while (next) {
1506     (*celldms)[n] = (const char *)next->obj->name;
1507     next          = next->next;
1508     ++n;
1509   }
1510   *Ndm = n;
1511   PetscFunctionReturn(PETSC_SUCCESS);
1512 }
1513 
1514 /*@
1515   DMSwarmSetCellDMActive - Activates a cell `DM` for a `DMSWARM`
1516 
1517   Collective
1518 
1519   Input Parameters:
1520 + sw   - a `DMSWARM`
1521 - name - name of the cell `DM` to active for the `DMSWARM`
1522 
1523   Level: beginner
1524 
1525   Note:
1526   The attached `DM` (dmcell) will be queried for point location and
1527   neighbor MPI-rank information if `DMSwarmMigrate()` is called.
1528 
1529 .seealso: `DM`, `DMSWARM`, `DMSwarmCellDM`, `DMSwarmSetType()`, `DMSwarmAddCellDM()`, `DMSwarmSetCellDM()`, `DMSwarmMigrate()`
1530 @*/
1531 PetscErrorCode DMSwarmSetCellDMActive(DM sw, const char name[])
1532 {
1533   DM_Swarm     *swarm = (DM_Swarm *)sw->data;
1534   DMSwarmCellDM celldm;
1535 
1536   PetscFunctionBegin;
1537   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
1538   PetscCall(PetscInfo(sw, "Setting cell DM to %s\n", name));
1539   PetscCall(PetscFree(swarm->activeCellDM));
1540   PetscCall(PetscStrallocpy(name, (char **)&swarm->activeCellDM));
1541   PetscCall(DMSwarmGetCellDMActive(sw, &celldm));
1542   PetscFunctionReturn(PETSC_SUCCESS);
1543 }
1544 
1545 /*@
1546   DMSwarmGetCellDMActive - Returns the active cell `DM` for a `DMSWARM`
1547 
1548   Collective
1549 
1550   Input Parameter:
1551 . sw - a `DMSWARM`
1552 
1553   Output Parameter:
1554 . celldm - the active `DMSwarmCellDM`
1555 
1556   Level: beginner
1557 
1558 .seealso: `DM`, `DMSWARM`, `DMSwarmCellDM`, `DMSwarmSetType()`, `DMSwarmAddCellDM()`, `DMSwarmSetCellDM()`, `DMSwarmMigrate()`
1559 @*/
1560 PetscErrorCode DMSwarmGetCellDMActive(DM sw, DMSwarmCellDM *celldm)
1561 {
1562   DM_Swarm *swarm = (DM_Swarm *)sw->data;
1563 
1564   PetscFunctionBegin;
1565   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
1566   PetscAssertPointer(celldm, 2);
1567   PetscCheck(swarm->activeCellDM, PetscObjectComm((PetscObject)sw), PETSC_ERR_ARG_WRONGSTATE, "Swarm has no active cell DM");
1568   PetscCall(PetscObjectListFind(swarm->cellDMs, swarm->activeCellDM, (PetscObject *)celldm));
1569   PetscCheck(*celldm, PetscObjectComm((PetscObject)sw), PETSC_ERR_ARG_WRONGSTATE, "Swarm has no valid cell DM for %s", swarm->activeCellDM);
1570   PetscFunctionReturn(PETSC_SUCCESS);
1571 }
1572 
1573 /*@C
1574   DMSwarmGetCellDMByName - Get a `DMSwarmCellDM` from its name
1575 
1576   Not collective
1577 
1578   Input Parameters:
1579 + sw   - a `DMSWARM`
1580 - name - the name
1581 
1582   Output Parameter:
1583 . celldm - the `DMSwarmCellDM`
1584 
1585   Level: beginner
1586 
1587 .seealso: `DM`, `DMSWARM`, `DMSwarmSetCellDM()`, `DMSwarmGetCellDMNames()`
1588 @*/
1589 PetscErrorCode DMSwarmGetCellDMByName(DM sw, const char name[], DMSwarmCellDM *celldm)
1590 {
1591   DM_Swarm *swarm = (DM_Swarm *)sw->data;
1592 
1593   PetscFunctionBegin;
1594   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
1595   PetscAssertPointer(name, 2);
1596   PetscAssertPointer(celldm, 3);
1597   PetscCall(PetscObjectListFind(swarm->cellDMs, name, (PetscObject *)celldm));
1598   PetscCheck(*celldm, PetscObjectComm((PetscObject)sw), PETSC_ERR_ARG_WRONGSTATE, "Swarm has no valid cell DM for %s", name);
1599   PetscFunctionReturn(PETSC_SUCCESS);
1600 }
1601 
1602 /*@
1603   DMSwarmAddCellDM - Adds a cell `DM` to the `DMSWARM`
1604 
1605   Collective
1606 
1607   Input Parameters:
1608 + sw     - a `DMSWARM`
1609 - celldm - the `DMSwarmCellDM`
1610 
1611   Level: beginner
1612 
1613   Note:
1614   Cell DMs with the same name will share the cellid field
1615 
1616 .seealso: `DM`, `DMSWARM`, `DMSwarmSetType()`, `DMSwarmPushCellDM()`, `DMSwarmSetCellDM()`, `DMSwarmMigrate()`
1617 @*/
1618 PetscErrorCode DMSwarmAddCellDM(DM sw, DMSwarmCellDM celldm)
1619 {
1620   DM_Swarm   *swarm = (DM_Swarm *)sw->data;
1621   const char *name;
1622   PetscInt    dim;
1623   PetscBool   flg;
1624   MPI_Comm    comm;
1625 
1626   PetscFunctionBegin;
1627   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
1628   PetscCall(PetscObjectGetComm((PetscObject)sw, &comm));
1629   PetscValidHeaderSpecific(celldm, DMSWARMCELLDM_CLASSID, 2);
1630   PetscCall(PetscObjectGetName((PetscObject)celldm, &name));
1631   PetscCall(PetscObjectListAdd(&swarm->cellDMs, name, (PetscObject)celldm));
1632   PetscCall(DMGetDimension(sw, &dim));
1633   for (PetscInt f = 0; f < celldm->Nfc; ++f) {
1634     PetscCall(DMSwarmDataFieldStringInList(celldm->coordFields[f], swarm->db->nfields, (const DMSwarmDataField *)swarm->db->field, &flg));
1635     if (!flg) {
1636       PetscCall(DMSwarmRegisterPetscDatatypeField(sw, celldm->coordFields[f], dim, PETSC_DOUBLE));
1637     } else {
1638       PetscDataType dt;
1639       PetscInt      bs;
1640 
1641       PetscCall(DMSwarmGetFieldInfo(sw, celldm->coordFields[f], &bs, &dt));
1642       PetscCheck(bs == dim, comm, PETSC_ERR_ARG_WRONG, "Coordinate field %s has blocksize %" PetscInt_FMT " != %" PetscInt_FMT " spatial dimension", celldm->coordFields[f], bs, dim);
1643       PetscCheck(dt == PETSC_DOUBLE, comm, PETSC_ERR_ARG_WRONG, "Coordinate field %s has datatype %s != PETSC_DOUBLE", celldm->coordFields[f], PetscDataTypes[dt]);
1644     }
1645   }
1646   // Assume that DMs with the same name share the cellid field
1647   PetscCall(DMSwarmDataFieldStringInList(celldm->cellid, swarm->db->nfields, (const DMSwarmDataField *)swarm->db->field, &flg));
1648   if (!flg) {
1649     PetscBool   isShell, isDummy;
1650     const char *name;
1651 
1652     // Allow dummy DMSHELL (I don't think we should support this mode)
1653     PetscCall(PetscObjectTypeCompare((PetscObject)celldm->dm, DMSHELL, &isShell));
1654     PetscCall(PetscObjectGetName((PetscObject)celldm->dm, &name));
1655     PetscCall(PetscStrcmp(name, "dummy", &isDummy));
1656     if (!isShell || !isDummy) PetscCall(DMSwarmRegisterPetscDatatypeField(sw, celldm->cellid, 1, PETSC_INT));
1657   }
1658   PetscCall(DMSwarmSetCellDMActive(sw, name));
1659   PetscFunctionReturn(PETSC_SUCCESS);
1660 }
1661 
1662 /*@
1663   DMSwarmGetLocalSize - Retrieves the local length of fields registered
1664 
1665   Not Collective
1666 
1667   Input Parameter:
1668 . dm - a `DMSWARM`
1669 
1670   Output Parameter:
1671 . nlocal - the length of each registered field
1672 
1673   Level: beginner
1674 
1675 .seealso: `DM`, `DMSWARM`, `DMSwarmGetSize()`, `DMSwarmSetLocalSizes()`
1676 @*/
1677 PetscErrorCode DMSwarmGetLocalSize(DM dm, PetscInt *nlocal)
1678 {
1679   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1680 
1681   PetscFunctionBegin;
1682   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, nlocal, NULL, NULL));
1683   PetscFunctionReturn(PETSC_SUCCESS);
1684 }
1685 
1686 /*@
1687   DMSwarmGetSize - Retrieves the total length of fields registered
1688 
1689   Collective
1690 
1691   Input Parameter:
1692 . dm - a `DMSWARM`
1693 
1694   Output Parameter:
1695 . n - the total length of each registered field
1696 
1697   Level: beginner
1698 
1699   Note:
1700   This calls `MPI_Allreduce()` upon each call (inefficient but safe)
1701 
1702 .seealso: `DM`, `DMSWARM`, `DMSwarmGetLocalSize()`, `DMSwarmSetLocalSizes()`
1703 @*/
1704 PetscErrorCode DMSwarmGetSize(DM dm, PetscInt *n)
1705 {
1706   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1707   PetscInt  nlocal;
1708 
1709   PetscFunctionBegin;
1710   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &nlocal, NULL, NULL));
1711   PetscCallMPI(MPIU_Allreduce(&nlocal, n, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm)));
1712   PetscFunctionReturn(PETSC_SUCCESS);
1713 }
1714 
1715 /*@C
1716   DMSwarmRegisterPetscDatatypeField - Register a field to a `DMSWARM` with a native PETSc data type
1717 
1718   Collective
1719 
1720   Input Parameters:
1721 + dm        - a `DMSWARM`
1722 . fieldname - the textual name to identify this field
1723 . blocksize - the number of each data type
1724 - type      - a valid PETSc data type (`PETSC_CHAR`, `PETSC_SHORT`, `PETSC_INT`, `PETSC_FLOAT`, `PETSC_REAL`, `PETSC_LONG`)
1725 
1726   Level: beginner
1727 
1728   Notes:
1729   The textual name for each registered field must be unique.
1730 
1731 .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterUserStructField()`, `DMSwarmRegisterUserDatatypeField()`
1732 @*/
1733 PetscErrorCode DMSwarmRegisterPetscDatatypeField(DM dm, const char fieldname[], PetscInt blocksize, PetscDataType type)
1734 {
1735   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1736   size_t    size;
1737 
1738   PetscFunctionBegin;
1739   PetscCheck(swarm->field_registration_initialized, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMSwarmInitializeFieldRegister() first");
1740   PetscCheck(!swarm->field_registration_finalized, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Cannot register additional fields after calling DMSwarmFinalizeFieldRegister() first");
1741 
1742   PetscCheck(type != PETSC_OBJECT, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
1743   PetscCheck(type != PETSC_FUNCTION, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
1744   PetscCheck(type != PETSC_STRING, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
1745   PetscCheck(type != PETSC_STRUCT, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
1746   PetscCheck(type != PETSC_DATATYPE_UNKNOWN, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}");
1747 
1748   PetscCall(PetscDataTypeGetSize(type, &size));
1749   /* Load a specific data type into data bucket, specifying textual name and its size in bytes */
1750   PetscCall(DMSwarmDataBucketRegisterField(swarm->db, "DMSwarmRegisterPetscDatatypeField", fieldname, blocksize * size, NULL));
1751   {
1752     DMSwarmDataField gfield;
1753 
1754     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
1755     PetscCall(DMSwarmDataFieldSetBlockSize(gfield, blocksize));
1756   }
1757   swarm->db->field[swarm->db->nfields - 1]->petsc_type = type;
1758   PetscFunctionReturn(PETSC_SUCCESS);
1759 }
1760 
1761 /*@C
1762   DMSwarmRegisterUserStructField - Register a user defined struct to a `DMSWARM`
1763 
1764   Collective
1765 
1766   Input Parameters:
1767 + dm        - a `DMSWARM`
1768 . fieldname - the textual name to identify this field
1769 - size      - the size in bytes of the user struct of each data type
1770 
1771   Level: beginner
1772 
1773   Note:
1774   The textual name for each registered field must be unique.
1775 
1776 .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmRegisterUserDatatypeField()`
1777 @*/
1778 PetscErrorCode DMSwarmRegisterUserStructField(DM dm, const char fieldname[], size_t size)
1779 {
1780   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1781 
1782   PetscFunctionBegin;
1783   PetscCall(DMSwarmDataBucketRegisterField(swarm->db, "DMSwarmRegisterUserStructField", fieldname, size, NULL));
1784   swarm->db->field[swarm->db->nfields - 1]->petsc_type = PETSC_STRUCT;
1785   PetscFunctionReturn(PETSC_SUCCESS);
1786 }
1787 
1788 /*@C
1789   DMSwarmRegisterUserDatatypeField - Register a user defined data type to a `DMSWARM`
1790 
1791   Collective
1792 
1793   Input Parameters:
1794 + dm        - a `DMSWARM`
1795 . fieldname - the textual name to identify this field
1796 . size      - the size in bytes of the user data type
1797 - blocksize - the number of each data type
1798 
1799   Level: beginner
1800 
1801   Note:
1802   The textual name for each registered field must be unique.
1803 
1804 .seealso: `DM`, `DMSWARM`, `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmRegisterUserStructField()`
1805 @*/
1806 PetscErrorCode DMSwarmRegisterUserDatatypeField(DM dm, const char fieldname[], size_t size, PetscInt blocksize)
1807 {
1808   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1809 
1810   PetscFunctionBegin;
1811   PetscCall(DMSwarmDataBucketRegisterField(swarm->db, "DMSwarmRegisterUserDatatypeField", fieldname, blocksize * size, NULL));
1812   {
1813     DMSwarmDataField gfield;
1814 
1815     PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
1816     PetscCall(DMSwarmDataFieldSetBlockSize(gfield, blocksize));
1817   }
1818   swarm->db->field[swarm->db->nfields - 1]->petsc_type = PETSC_DATATYPE_UNKNOWN;
1819   PetscFunctionReturn(PETSC_SUCCESS);
1820 }
1821 
1822 /*@C
1823   DMSwarmGetField - Get access to the underlying array storing all entries associated with a registered field
1824 
1825   Not Collective, No Fortran Support
1826 
1827   Input Parameters:
1828 + dm        - a `DMSWARM`
1829 - fieldname - the textual name to identify this field
1830 
1831   Output Parameters:
1832 + blocksize - the number of each data type
1833 . type      - the data type
1834 - data      - pointer to raw array
1835 
1836   Level: beginner
1837 
1838   Notes:
1839   The array must be returned using a matching call to `DMSwarmRestoreField()`.
1840 
1841   Fortran Note:
1842   Only works for `type` of `PETSC_SCALAR`
1843 
1844 .seealso: `DM`, `DMSWARM`, `DMSwarmRestoreField()`
1845 @*/
1846 PetscErrorCode DMSwarmGetField(DM dm, const char fieldname[], PetscInt *blocksize, PetscDataType *type, void **data) PeNS
1847 {
1848   DM_Swarm        *swarm = (DM_Swarm *)dm->data;
1849   DMSwarmDataField gfield;
1850 
1851   PetscFunctionBegin;
1852   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
1853   if (!swarm->issetup) PetscCall(DMSetUp(dm));
1854   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
1855   PetscCall(DMSwarmDataFieldGetAccess(gfield));
1856   PetscCall(DMSwarmDataFieldGetEntries(gfield, data));
1857   if (blocksize) *blocksize = gfield->bs;
1858   if (type) *type = gfield->petsc_type;
1859   PetscFunctionReturn(PETSC_SUCCESS);
1860 }
1861 
1862 /*@C
1863   DMSwarmRestoreField - Restore access to the underlying array storing all entries associated with a registered field
1864 
1865   Not Collective
1866 
1867   Input Parameters:
1868 + dm        - a `DMSWARM`
1869 - fieldname - the textual name to identify this field
1870 
1871   Output Parameters:
1872 + blocksize - the number of each data type
1873 . type      - the data type
1874 - data      - pointer to raw array
1875 
1876   Level: beginner
1877 
1878   Notes:
1879   The user must call `DMSwarmGetField()` prior to calling `DMSwarmRestoreField()`.
1880 
1881   Fortran Note:
1882   Only works for `type` of `PETSC_SCALAR`
1883 
1884 .seealso: `DM`, `DMSWARM`, `DMSwarmGetField()`
1885 @*/
1886 PetscErrorCode DMSwarmRestoreField(DM dm, const char fieldname[], PetscInt *blocksize, PetscDataType *type, void **data) PeNS
1887 {
1888   DM_Swarm        *swarm = (DM_Swarm *)dm->data;
1889   DMSwarmDataField gfield;
1890 
1891   PetscFunctionBegin;
1892   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
1893   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
1894   PetscCall(DMSwarmDataFieldRestoreAccess(gfield));
1895   if (data) *data = NULL;
1896   PetscFunctionReturn(PETSC_SUCCESS);
1897 }
1898 
1899 PetscErrorCode DMSwarmGetFieldInfo(DM dm, const char fieldname[], PetscInt *blocksize, PetscDataType *type)
1900 {
1901   DM_Swarm        *swarm = (DM_Swarm *)dm->data;
1902   DMSwarmDataField gfield;
1903 
1904   PetscFunctionBegin;
1905   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM);
1906   PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield));
1907   if (blocksize) *blocksize = gfield->bs;
1908   if (type) *type = gfield->petsc_type;
1909   PetscFunctionReturn(PETSC_SUCCESS);
1910 }
1911 
1912 /*@
1913   DMSwarmAddPoint - Add space for one new point in the `DMSWARM`
1914 
1915   Not Collective
1916 
1917   Input Parameter:
1918 . dm - a `DMSWARM`
1919 
1920   Level: beginner
1921 
1922   Notes:
1923   The new point will have all fields initialized to zero.
1924 
1925 .seealso: `DM`, `DMSWARM`, `DMSwarmAddNPoints()`
1926 @*/
1927 PetscErrorCode DMSwarmAddPoint(DM dm)
1928 {
1929   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1930 
1931   PetscFunctionBegin;
1932   if (!swarm->issetup) PetscCall(DMSetUp(dm));
1933   PetscCall(PetscLogEventBegin(DMSWARM_AddPoints, 0, 0, 0, 0));
1934   PetscCall(DMSwarmDataBucketAddPoint(swarm->db));
1935   PetscCall(PetscLogEventEnd(DMSWARM_AddPoints, 0, 0, 0, 0));
1936   PetscFunctionReturn(PETSC_SUCCESS);
1937 }
1938 
1939 /*@
1940   DMSwarmAddNPoints - Add space for a number of new points in the `DMSWARM`
1941 
1942   Not Collective
1943 
1944   Input Parameters:
1945 + dm      - a `DMSWARM`
1946 - npoints - the number of new points to add
1947 
1948   Level: beginner
1949 
1950   Notes:
1951   The new point will have all fields initialized to zero.
1952 
1953 .seealso: `DM`, `DMSWARM`, `DMSwarmAddPoint()`
1954 @*/
1955 PetscErrorCode DMSwarmAddNPoints(DM dm, PetscInt npoints)
1956 {
1957   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1958   PetscInt  nlocal;
1959 
1960   PetscFunctionBegin;
1961   PetscCall(PetscLogEventBegin(DMSWARM_AddPoints, 0, 0, 0, 0));
1962   PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &nlocal, NULL, NULL));
1963   nlocal = PetscMax(nlocal, 0) + npoints;
1964   PetscCall(DMSwarmDataBucketSetSizes(swarm->db, nlocal, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT));
1965   PetscCall(PetscLogEventEnd(DMSWARM_AddPoints, 0, 0, 0, 0));
1966   PetscFunctionReturn(PETSC_SUCCESS);
1967 }
1968 
1969 /*@
1970   DMSwarmRemovePoint - Remove the last point from the `DMSWARM`
1971 
1972   Not Collective
1973 
1974   Input Parameter:
1975 . dm - a `DMSWARM`
1976 
1977   Level: beginner
1978 
1979 .seealso: `DM`, `DMSWARM`, `DMSwarmRemovePointAtIndex()`
1980 @*/
1981 PetscErrorCode DMSwarmRemovePoint(DM dm)
1982 {
1983   DM_Swarm *swarm = (DM_Swarm *)dm->data;
1984 
1985   PetscFunctionBegin;
1986   PetscCall(PetscLogEventBegin(DMSWARM_RemovePoints, 0, 0, 0, 0));
1987   PetscCall(DMSwarmDataBucketRemovePoint(swarm->db));
1988   PetscCall(PetscLogEventEnd(DMSWARM_RemovePoints, 0, 0, 0, 0));
1989   PetscFunctionReturn(PETSC_SUCCESS);
1990 }
1991 
1992 /*@
1993   DMSwarmRemovePointAtIndex - Removes a specific point from the `DMSWARM`
1994 
1995   Not Collective
1996 
1997   Input Parameters:
1998 + dm  - a `DMSWARM`
1999 - idx - index of point to remove
2000 
2001   Level: beginner
2002 
2003 .seealso: `DM`, `DMSWARM`, `DMSwarmRemovePoint()`
2004 @*/
2005 PetscErrorCode DMSwarmRemovePointAtIndex(DM dm, PetscInt idx)
2006 {
2007   DM_Swarm *swarm = (DM_Swarm *)dm->data;
2008 
2009   PetscFunctionBegin;
2010   PetscCall(PetscLogEventBegin(DMSWARM_RemovePoints, 0, 0, 0, 0));
2011   PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db, idx));
2012   PetscCall(PetscLogEventEnd(DMSWARM_RemovePoints, 0, 0, 0, 0));
2013   PetscFunctionReturn(PETSC_SUCCESS);
2014 }
2015 
2016 /*@
2017   DMSwarmCopyPoint - Copy point pj to point pi in the `DMSWARM`
2018 
2019   Not Collective
2020 
2021   Input Parameters:
2022 + dm - a `DMSWARM`
2023 . pi - the index of the point to copy
2024 - pj - the point index where the copy should be located
2025 
2026   Level: beginner
2027 
2028 .seealso: `DM`, `DMSWARM`, `DMSwarmRemovePoint()`
2029 @*/
2030 PetscErrorCode DMSwarmCopyPoint(DM dm, PetscInt pi, PetscInt pj)
2031 {
2032   DM_Swarm *swarm = (DM_Swarm *)dm->data;
2033 
2034   PetscFunctionBegin;
2035   if (!swarm->issetup) PetscCall(DMSetUp(dm));
2036   PetscCall(DMSwarmDataBucketCopyPoint(swarm->db, pi, swarm->db, pj));
2037   PetscFunctionReturn(PETSC_SUCCESS);
2038 }
2039 
2040 static PetscErrorCode DMSwarmMigrate_Basic(DM dm, PetscBool remove_sent_points)
2041 {
2042   PetscFunctionBegin;
2043   PetscCall(DMSwarmMigrate_Push_Basic(dm, remove_sent_points));
2044   PetscFunctionReturn(PETSC_SUCCESS);
2045 }
2046 
2047 /*@
2048   DMSwarmMigrate - Relocates points defined in the `DMSWARM` to other MPI-ranks
2049 
2050   Collective
2051 
2052   Input Parameters:
2053 + dm                 - the `DMSWARM`
2054 - remove_sent_points - flag indicating if sent points should be removed from the current MPI-rank
2055 
2056   Level: advanced
2057 
2058   Notes:
2059   The `DM` will be modified to accommodate received points.
2060   If `remove_sent_points` is `PETSC_TRUE`, any points that were sent will be removed from the `DM`.
2061   Different styles of migration are supported. See `DMSwarmSetMigrateType()`.
2062 
2063 .seealso: `DM`, `DMSWARM`, `DMSwarmSetMigrateType()`
2064 @*/
2065 PetscErrorCode DMSwarmMigrate(DM dm, PetscBool remove_sent_points)
2066 {
2067   DM_Swarm *swarm = (DM_Swarm *)dm->data;
2068 
2069   PetscFunctionBegin;
2070   PetscCall(PetscLogEventBegin(DMSWARM_Migrate, 0, 0, 0, 0));
2071   switch (swarm->migrate_type) {
2072   case DMSWARM_MIGRATE_BASIC:
2073     PetscCall(DMSwarmMigrate_Basic(dm, remove_sent_points));
2074     break;
2075   case DMSWARM_MIGRATE_DMCELLNSCATTER:
2076     PetscCall(DMSwarmMigrate_CellDMScatter(dm, remove_sent_points));
2077     break;
2078   case DMSWARM_MIGRATE_DMCELLEXACT:
2079     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_MIGRATE_DMCELLEXACT not implemented");
2080   case DMSWARM_MIGRATE_USER:
2081     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_MIGRATE_USER not implemented");
2082   default:
2083     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_MIGRATE type unknown");
2084   }
2085   PetscCall(PetscLogEventEnd(DMSWARM_Migrate, 0, 0, 0, 0));
2086   PetscCall(DMClearGlobalVectors(dm));
2087   PetscFunctionReturn(PETSC_SUCCESS);
2088 }
2089 
2090 PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm, PetscInt *globalsize);
2091 
2092 /*
2093  DMSwarmCollectViewCreate
2094 
2095  * Applies a collection method and gathers point neighbour points into dm
2096 
2097  Notes:
2098  Users should call DMSwarmCollectViewDestroy() after
2099  they have finished computations associated with the collected points
2100 */
2101 
2102 /*@
2103   DMSwarmCollectViewCreate - Applies a collection method and gathers points
2104   in neighbour ranks into the `DMSWARM`
2105 
2106   Collective
2107 
2108   Input Parameter:
2109 . dm - the `DMSWARM`
2110 
2111   Level: advanced
2112 
2113   Notes:
2114   Users should call `DMSwarmCollectViewDestroy()` after
2115   they have finished computations associated with the collected points
2116 
2117   Different collect methods are supported. See `DMSwarmSetCollectType()`.
2118 
2119   Developer Note:
2120   Create and Destroy routines create new objects that can get destroyed, they do not change the state
2121   of the current object.
2122 
2123 .seealso: `DM`, `DMSWARM`, `DMSwarmCollectViewDestroy()`, `DMSwarmSetCollectType()`
2124 @*/
2125 PetscErrorCode DMSwarmCollectViewCreate(DM dm)
2126 {
2127   DM_Swarm *swarm = (DM_Swarm *)dm->data;
2128   PetscInt  ng;
2129 
2130   PetscFunctionBegin;
2131   PetscCheck(!swarm->collect_view_active, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "CollectView currently active");
2132   PetscCall(DMSwarmGetLocalSize(dm, &ng));
2133   switch (swarm->collect_type) {
2134   case DMSWARM_COLLECT_BASIC:
2135     PetscCall(DMSwarmMigrate_GlobalToLocal_Basic(dm, &ng));
2136     break;
2137   case DMSWARM_COLLECT_DMDABOUNDINGBOX:
2138     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_COLLECT_DMDABOUNDINGBOX not implemented");
2139   case DMSWARM_COLLECT_GENERAL:
2140     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_COLLECT_GENERAL not implemented");
2141   default:
2142     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_COLLECT type unknown");
2143   }
2144   swarm->collect_view_active       = PETSC_TRUE;
2145   swarm->collect_view_reset_nlocal = ng;
2146   PetscFunctionReturn(PETSC_SUCCESS);
2147 }
2148 
2149 /*@
2150   DMSwarmCollectViewDestroy - Resets the `DMSWARM` to the size prior to calling `DMSwarmCollectViewCreate()`
2151 
2152   Collective
2153 
2154   Input Parameters:
2155 . dm - the `DMSWARM`
2156 
2157   Notes:
2158   Users should call `DMSwarmCollectViewCreate()` before this function is called.
2159 
2160   Level: advanced
2161 
2162   Developer Note:
2163   Create and Destroy routines create new objects that can get destroyed, they do not change the state
2164   of the current object.
2165 
2166 .seealso: `DM`, `DMSWARM`, `DMSwarmCollectViewCreate()`, `DMSwarmSetCollectType()`
2167 @*/
2168 PetscErrorCode DMSwarmCollectViewDestroy(DM dm)
2169 {
2170   DM_Swarm *swarm = (DM_Swarm *)dm->data;
2171 
2172   PetscFunctionBegin;
2173   PetscCheck(swarm->collect_view_active, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "CollectView is currently not active");
2174   PetscCall(DMSwarmSetLocalSizes(dm, swarm->collect_view_reset_nlocal, -1));
2175   swarm->collect_view_active = PETSC_FALSE;
2176   PetscFunctionReturn(PETSC_SUCCESS);
2177 }
2178 
2179 static PetscErrorCode DMSwarmSetUpPIC(DM dm)
2180 {
2181   PetscInt dim;
2182 
2183   PetscFunctionBegin;
2184   PetscCall(DMSwarmSetNumSpecies(dm, 1));
2185   PetscCall(DMGetDimension(dm, &dim));
2186   PetscCheck(dim >= 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Dimension must be 1,2,3 - found %" PetscInt_FMT, dim);
2187   PetscCheck(dim <= 3, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Dimension must be 1,2,3 - found %" PetscInt_FMT, dim);
2188   PetscFunctionReturn(PETSC_SUCCESS);
2189 }
2190 
2191 /*@
2192   DMSwarmSetPointCoordinatesRandom - Sets initial coordinates for particles in each cell
2193 
2194   Collective
2195 
2196   Input Parameters:
2197 + dm  - the `DMSWARM`
2198 - Npc - The number of particles per cell in the cell `DM`
2199 
2200   Level: intermediate
2201 
2202   Notes:
2203   The user must use `DMSwarmSetCellDM()` to set the cell `DM` first. The particles are placed randomly inside each cell. If only
2204   one particle is in each cell, it is placed at the centroid.
2205 
2206 .seealso: `DM`, `DMSWARM`, `DMSwarmSetCellDM()`
2207 @*/
2208 PetscErrorCode DMSwarmSetPointCoordinatesRandom(DM dm, PetscInt Npc)
2209 {
2210   DM             cdm;
2211   DMSwarmCellDM  celldm;
2212   PetscRandom    rnd;
2213   DMPolytopeType ct;
2214   PetscBool      simplex;
2215   PetscReal     *centroid, *coords, *xi0, *v0, *J, *invJ, detJ;
2216   PetscInt       dim, d, cStart, cEnd, c, p, Nfc;
2217   const char   **coordFields;
2218 
2219   PetscFunctionBeginUser;
2220   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rnd));
2221   PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0));
2222   PetscCall(PetscRandomSetType(rnd, PETSCRAND48));
2223 
2224   PetscCall(DMSwarmGetCellDMActive(dm, &celldm));
2225   PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Nfc, &coordFields));
2226   PetscCheck(Nfc == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "We only support a single coordinate field right now, not %" PetscInt_FMT, Nfc);
2227   PetscCall(DMSwarmGetCellDM(dm, &cdm));
2228   PetscCall(DMGetDimension(cdm, &dim));
2229   PetscCall(DMPlexGetHeightStratum(cdm, 0, &cStart, &cEnd));
2230   PetscCall(DMPlexGetCellType(cdm, cStart, &ct));
2231   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
2232 
2233   PetscCall(PetscMalloc5(dim, &centroid, dim, &xi0, dim, &v0, dim * dim, &J, dim * dim, &invJ));
2234   for (d = 0; d < dim; ++d) xi0[d] = -1.0;
2235   PetscCall(DMSwarmGetField(dm, coordFields[0], NULL, NULL, (void **)&coords));
2236   for (c = cStart; c < cEnd; ++c) {
2237     if (Npc == 1) {
2238       PetscCall(DMPlexComputeCellGeometryFVM(cdm, c, NULL, centroid, NULL));
2239       for (d = 0; d < dim; ++d) coords[c * dim + d] = centroid[d];
2240     } else {
2241       PetscCall(DMPlexComputeCellGeometryFEM(cdm, c, NULL, v0, J, invJ, &detJ)); /* affine */
2242       for (p = 0; p < Npc; ++p) {
2243         const PetscInt n   = c * Npc + p;
2244         PetscReal      sum = 0.0, refcoords[3];
2245 
2246         for (d = 0; d < dim; ++d) {
2247           PetscCall(PetscRandomGetValueReal(rnd, &refcoords[d]));
2248           sum += refcoords[d];
2249         }
2250         if (simplex && sum > 0.0)
2251           for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim) * sum;
2252         CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n * dim]);
2253       }
2254     }
2255   }
2256   PetscCall(DMSwarmRestoreField(dm, coordFields[0], NULL, NULL, (void **)&coords));
2257   PetscCall(PetscFree5(centroid, xi0, v0, J, invJ));
2258   PetscCall(PetscRandomDestroy(&rnd));
2259   PetscFunctionReturn(PETSC_SUCCESS);
2260 }
2261 
2262 /*@
2263   DMSwarmGetType - Get particular flavor of `DMSWARM`
2264 
2265   Collective
2266 
2267   Input Parameter:
2268 . sw - the `DMSWARM`
2269 
2270   Output Parameter:
2271 . stype - the `DMSWARM` type (e.g. `DMSWARM_PIC`)
2272 
2273   Level: advanced
2274 
2275 .seealso: `DM`, `DMSWARM`, `DMSwarmSetMigrateType()`, `DMSwarmSetCollectType()`, `DMSwarmType`, `DMSWARM_PIC`, `DMSWARM_BASIC`
2276 @*/
2277 PetscErrorCode DMSwarmGetType(DM sw, DMSwarmType *stype)
2278 {
2279   DM_Swarm *swarm = (DM_Swarm *)sw->data;
2280 
2281   PetscFunctionBegin;
2282   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
2283   PetscAssertPointer(stype, 2);
2284   *stype = swarm->swarm_type;
2285   PetscFunctionReturn(PETSC_SUCCESS);
2286 }
2287 
2288 /*@
2289   DMSwarmSetType - Set particular flavor of `DMSWARM`
2290 
2291   Collective
2292 
2293   Input Parameters:
2294 + sw    - the `DMSWARM`
2295 - stype - the `DMSWARM` type (e.g. `DMSWARM_PIC`)
2296 
2297   Level: advanced
2298 
2299 .seealso: `DM`, `DMSWARM`, `DMSwarmSetMigrateType()`, `DMSwarmSetCollectType()`, `DMSwarmType`, `DMSWARM_PIC`, `DMSWARM_BASIC`
2300 @*/
2301 PetscErrorCode DMSwarmSetType(DM sw, DMSwarmType stype)
2302 {
2303   DM_Swarm *swarm = (DM_Swarm *)sw->data;
2304 
2305   PetscFunctionBegin;
2306   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
2307   swarm->swarm_type = stype;
2308   if (swarm->swarm_type == DMSWARM_PIC) PetscCall(DMSwarmSetUpPIC(sw));
2309   PetscFunctionReturn(PETSC_SUCCESS);
2310 }
2311 
2312 static PetscErrorCode DMSwarmCreateRemapDM_Private(DM sw, DM *rdm)
2313 {
2314   PetscFE        fe;
2315   DMPolytopeType ct;
2316   PetscInt       dim, cStart;
2317   const char    *prefix = "remap_";
2318 
2319   PetscFunctionBegin;
2320   PetscCall(DMCreate(PetscObjectComm((PetscObject)sw), rdm));
2321   PetscCall(DMSetType(*rdm, DMPLEX));
2322   PetscCall(DMPlexSetOptionsPrefix(*rdm, prefix));
2323   PetscCall(DMSetFromOptions(*rdm));
2324   PetscCall(PetscObjectSetName((PetscObject)*rdm, "remap"));
2325   PetscCall(DMViewFromOptions(*rdm, NULL, "-dm_view"));
2326 
2327   PetscCall(DMGetDimension(*rdm, &dim));
2328   PetscCall(DMPlexGetHeightStratum(*rdm, 0, &cStart, NULL));
2329   PetscCall(DMPlexGetCellType(*rdm, cStart, &ct));
2330   PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, 1, ct, prefix, PETSC_DETERMINE, &fe));
2331   PetscCall(PetscObjectSetName((PetscObject)fe, "distribution"));
2332   PetscCall(DMSetField(*rdm, 0, NULL, (PetscObject)fe));
2333   PetscCall(DMCreateDS(*rdm));
2334   PetscCall(PetscFEDestroy(&fe));
2335   PetscFunctionReturn(PETSC_SUCCESS);
2336 }
2337 
2338 static PetscErrorCode DMSetup_Swarm(DM sw)
2339 {
2340   DM_Swarm *swarm = (DM_Swarm *)sw->data;
2341 
2342   PetscFunctionBegin;
2343   if (swarm->issetup) PetscFunctionReturn(PETSC_SUCCESS);
2344   swarm->issetup = PETSC_TRUE;
2345 
2346   if (swarm->remap_type != DMSWARM_REMAP_NONE) {
2347     DMSwarmCellDM celldm;
2348     DM            rdm;
2349     const char   *fieldnames[2]  = {DMSwarmPICField_coor, "velocity"};
2350     const char   *vfieldnames[1] = {"w_q"};
2351 
2352     PetscCall(DMSwarmCreateRemapDM_Private(sw, &rdm));
2353     PetscCall(DMSwarmCellDMCreate(rdm, 1, vfieldnames, 2, fieldnames, &celldm));
2354     PetscCall(DMSwarmAddCellDM(sw, celldm));
2355     PetscCall(DMSwarmCellDMDestroy(&celldm));
2356     PetscCall(DMDestroy(&rdm));
2357   }
2358 
2359   if (swarm->swarm_type == DMSWARM_PIC) {
2360     DMSwarmCellDM celldm;
2361 
2362     PetscCall(DMSwarmGetCellDMActive(sw, &celldm));
2363     PetscCheck(celldm, PetscObjectComm((PetscObject)sw), PETSC_ERR_USER, "No active cell DM. DMSWARM_PIC requires you call DMSwarmSetCellDM() or DMSwarmAddCellDM()");
2364     if (celldm->dm->ops->locatepointssubdomain) {
2365       /* check methods exists for exact ownership identificiation */
2366       PetscCall(PetscInfo(sw, "DMSWARM_PIC: Using method CellDM->ops->LocatePointsSubdomain\n"));
2367       swarm->migrate_type = DMSWARM_MIGRATE_DMCELLEXACT;
2368     } else {
2369       /* check methods exist for point location AND rank neighbor identification */
2370       if (celldm->dm->ops->locatepoints) {
2371         PetscCall(PetscInfo(sw, "DMSWARM_PIC: Using method CellDM->LocatePoints\n"));
2372       } else SETERRQ(PetscObjectComm((PetscObject)sw), PETSC_ERR_USER, "DMSWARM_PIC requires the method CellDM->ops->locatepoints be defined");
2373 
2374       if (celldm->dm->ops->getneighbors) {
2375         PetscCall(PetscInfo(sw, "DMSWARM_PIC: Using method CellDM->GetNeigbors\n"));
2376       } else SETERRQ(PetscObjectComm((PetscObject)sw), PETSC_ERR_USER, "DMSWARM_PIC requires the method CellDM->ops->getneighbors be defined");
2377 
2378       swarm->migrate_type = DMSWARM_MIGRATE_DMCELLNSCATTER;
2379     }
2380   }
2381 
2382   PetscCall(DMSwarmFinalizeFieldRegister(sw));
2383 
2384   /* check some fields were registered */
2385   PetscCheck(swarm->db->nfields > 2, PetscObjectComm((PetscObject)sw), PETSC_ERR_USER, "At least one field user must be registered via DMSwarmRegisterXXX()");
2386   PetscFunctionReturn(PETSC_SUCCESS);
2387 }
2388 
2389 static PetscErrorCode DMDestroy_Swarm(DM dm)
2390 {
2391   DM_Swarm *swarm = (DM_Swarm *)dm->data;
2392 
2393   PetscFunctionBegin;
2394   if (--swarm->refct > 0) PetscFunctionReturn(PETSC_SUCCESS);
2395   PetscCall(PetscObjectListDestroy(&swarm->cellDMs));
2396   PetscCall(PetscFree(swarm->activeCellDM));
2397   PetscCall(DMSwarmDataBucketDestroy(&swarm->db));
2398   PetscCall(PetscFree(swarm));
2399   PetscFunctionReturn(PETSC_SUCCESS);
2400 }
2401 
2402 static PetscErrorCode DMSwarmView_Draw(DM dm, PetscViewer viewer)
2403 {
2404   DM            cdm;
2405   DMSwarmCellDM celldm;
2406   PetscDraw     draw;
2407   PetscReal    *coords, oldPause, radius = 0.01;
2408   PetscInt      Np, p, bs, Nfc;
2409   const char  **coordFields;
2410 
2411   PetscFunctionBegin;
2412   PetscCall(PetscOptionsGetReal(NULL, ((PetscObject)dm)->prefix, "-dm_view_swarm_radius", &radius, NULL));
2413   PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
2414   PetscCall(DMSwarmGetCellDM(dm, &cdm));
2415   PetscCall(PetscDrawGetPause(draw, &oldPause));
2416   PetscCall(PetscDrawSetPause(draw, 0.0));
2417   PetscCall(DMView(cdm, viewer));
2418   PetscCall(PetscDrawSetPause(draw, oldPause));
2419 
2420   PetscCall(DMSwarmGetCellDMActive(dm, &celldm));
2421   PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Nfc, &coordFields));
2422   PetscCheck(Nfc == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "We only support a single coordinate field right now, not %" PetscInt_FMT, Nfc);
2423   PetscCall(DMSwarmGetLocalSize(dm, &Np));
2424   PetscCall(DMSwarmGetField(dm, coordFields[0], &bs, NULL, (void **)&coords));
2425   for (p = 0; p < Np; ++p) {
2426     const PetscInt i = p * bs;
2427 
2428     PetscCall(PetscDrawEllipse(draw, coords[i], coords[i + 1], radius, radius, PETSC_DRAW_BLUE));
2429   }
2430   PetscCall(DMSwarmRestoreField(dm, coordFields[0], &bs, NULL, (void **)&coords));
2431   PetscCall(PetscDrawFlush(draw));
2432   PetscCall(PetscDrawPause(draw));
2433   PetscCall(PetscDrawSave(draw));
2434   PetscFunctionReturn(PETSC_SUCCESS);
2435 }
2436 
2437 static PetscErrorCode DMView_Swarm_Ascii(DM dm, PetscViewer viewer)
2438 {
2439   PetscViewerFormat format;
2440   PetscInt         *sizes;
2441   PetscInt          dim, Np, maxSize = 17;
2442   MPI_Comm          comm;
2443   PetscMPIInt       rank, size, p;
2444   const char       *name, *cellid;
2445 
2446   PetscFunctionBegin;
2447   PetscCall(PetscViewerGetFormat(viewer, &format));
2448   PetscCall(DMGetDimension(dm, &dim));
2449   PetscCall(DMSwarmGetLocalSize(dm, &Np));
2450   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
2451   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2452   PetscCallMPI(MPI_Comm_size(comm, &size));
2453   PetscCall(PetscObjectGetName((PetscObject)dm, &name));
2454   if (name) PetscCall(PetscViewerASCIIPrintf(viewer, "%s in %" PetscInt_FMT " dimension%s:\n", name, dim, dim == 1 ? "" : "s"));
2455   else PetscCall(PetscViewerASCIIPrintf(viewer, "Swarm in %" PetscInt_FMT " dimension%s:\n", dim, dim == 1 ? "" : "s"));
2456   if (size < maxSize) PetscCall(PetscCalloc1(size, &sizes));
2457   else PetscCall(PetscCalloc1(3, &sizes));
2458   if (size < maxSize) {
2459     PetscCallMPI(MPI_Gather(&Np, 1, MPIU_INT, sizes, 1, MPIU_INT, 0, comm));
2460     PetscCall(PetscViewerASCIIPrintf(viewer, "  Number of particles per rank:"));
2461     for (p = 0; p < size; ++p) {
2462       if (rank == 0) PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, sizes[p]));
2463     }
2464   } else {
2465     PetscInt locMinMax[2] = {Np, Np};
2466 
2467     PetscCall(PetscGlobalMinMaxInt(comm, locMinMax, sizes));
2468     PetscCall(PetscViewerASCIIPrintf(viewer, "  Min/Max of particles per rank: %" PetscInt_FMT "/%" PetscInt_FMT, sizes[0], sizes[1]));
2469   }
2470   PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
2471   PetscCall(PetscFree(sizes));
2472   if (format == PETSC_VIEWER_ASCII_INFO) {
2473     DMSwarmCellDM celldm;
2474     PetscInt     *cell;
2475 
2476     PetscCall(PetscViewerASCIIPrintf(viewer, "  Cells containing each particle:\n"));
2477     PetscCall(PetscViewerASCIIPushSynchronized(viewer));
2478     PetscCall(DMSwarmGetCellDMActive(dm, &celldm));
2479     PetscCall(DMSwarmCellDMGetCellID(celldm, &cellid));
2480     PetscCall(DMSwarmGetField(dm, cellid, NULL, NULL, (void **)&cell));
2481     for (p = 0; p < Np; ++p) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "  p%d: %" PetscInt_FMT "\n", p, cell[p]));
2482     PetscCall(PetscViewerFlush(viewer));
2483     PetscCall(PetscViewerASCIIPopSynchronized(viewer));
2484     PetscCall(DMSwarmRestoreField(dm, cellid, NULL, NULL, (void **)&cell));
2485   }
2486   PetscFunctionReturn(PETSC_SUCCESS);
2487 }
2488 
2489 static PetscErrorCode DMView_Swarm(DM dm, PetscViewer viewer)
2490 {
2491   DM_Swarm *swarm = (DM_Swarm *)dm->data;
2492   PetscBool iascii, ibinary, isvtk, isdraw, ispython;
2493 #if defined(PETSC_HAVE_HDF5)
2494   PetscBool ishdf5;
2495 #endif
2496 
2497   PetscFunctionBegin;
2498   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2499   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
2500   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
2501   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &ibinary));
2502   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk));
2503 #if defined(PETSC_HAVE_HDF5)
2504   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
2505 #endif
2506   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
2507   PetscCall(PetscObjectHasFunction((PetscObject)viewer, "PetscViewerPythonViewObject_C", &ispython));
2508   if (iascii) {
2509     PetscViewerFormat format;
2510 
2511     PetscCall(PetscViewerGetFormat(viewer, &format));
2512     switch (format) {
2513     case PETSC_VIEWER_ASCII_INFO_DETAIL:
2514       PetscCall(DMSwarmDataBucketView(PetscObjectComm((PetscObject)dm), swarm->db, NULL, DATABUCKET_VIEW_STDOUT));
2515       break;
2516     default:
2517       PetscCall(DMView_Swarm_Ascii(dm, viewer));
2518     }
2519   } else {
2520 #if defined(PETSC_HAVE_HDF5)
2521     if (ishdf5) PetscCall(DMSwarmView_HDF5(dm, viewer));
2522 #endif
2523     if (isdraw) PetscCall(DMSwarmView_Draw(dm, viewer));
2524     if (ispython) PetscCall(PetscViewerPythonViewObject(viewer, (PetscObject)dm));
2525   }
2526   PetscFunctionReturn(PETSC_SUCCESS);
2527 }
2528 
2529 /*@
2530   DMSwarmGetCellSwarm - Extracts a single cell from the `DMSWARM` object, returns it as a single cell `DMSWARM`.
2531   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.
2532 
2533   Noncollective
2534 
2535   Input Parameters:
2536 + sw        - the `DMSWARM`
2537 . cellID    - the integer id of the cell to be extracted and filtered
2538 - cellswarm - The `DMSWARM` to receive the cell
2539 
2540   Level: beginner
2541 
2542   Notes:
2543   This presently only supports `DMSWARM_PIC` type
2544 
2545   Should be restored with `DMSwarmRestoreCellSwarm()`
2546 
2547   Changes to this cell of the swarm will be lost if they are made prior to restoring this cell.
2548 
2549 .seealso: `DM`, `DMSWARM`, `DMSwarmRestoreCellSwarm()`
2550 @*/
2551 PetscErrorCode DMSwarmGetCellSwarm(DM sw, PetscInt cellID, DM cellswarm)
2552 {
2553   DM_Swarm   *original = (DM_Swarm *)sw->data;
2554   DMLabel     label;
2555   DM          dmc, subdmc;
2556   PetscInt   *pids, particles, dim;
2557   const char *name;
2558 
2559   PetscFunctionBegin;
2560   /* Configure new swarm */
2561   PetscCall(DMSetType(cellswarm, DMSWARM));
2562   PetscCall(DMGetDimension(sw, &dim));
2563   PetscCall(DMSetDimension(cellswarm, dim));
2564   PetscCall(DMSwarmSetType(cellswarm, DMSWARM_PIC));
2565   /* Destroy the unused, unconfigured data bucket to prevent stragglers in memory */
2566   PetscCall(DMSwarmDataBucketDestroy(&((DM_Swarm *)cellswarm->data)->db));
2567   PetscCall(DMSwarmSortGetAccess(sw));
2568   PetscCall(DMSwarmSortGetNumberOfPointsPerCell(sw, cellID, &particles));
2569   PetscCall(DMSwarmSortGetPointsPerCell(sw, cellID, &particles, &pids));
2570   PetscCall(DMSwarmDataBucketCreateFromSubset(original->db, particles, pids, &((DM_Swarm *)cellswarm->data)->db));
2571   PetscCall(DMSwarmSortRestoreAccess(sw));
2572   PetscCall(DMSwarmSortRestorePointsPerCell(sw, cellID, &particles, &pids));
2573   PetscCall(DMSwarmGetCellDM(sw, &dmc));
2574   PetscCall(DMLabelCreate(PetscObjectComm((PetscObject)sw), "singlecell", &label));
2575   PetscCall(DMAddLabel(dmc, label));
2576   PetscCall(DMLabelSetValue(label, cellID, 1));
2577   PetscCall(DMPlexFilter(dmc, label, 1, PETSC_FALSE, PETSC_FALSE, NULL, &subdmc));
2578   PetscCall(PetscObjectGetName((PetscObject)dmc, &name));
2579   PetscCall(PetscObjectSetName((PetscObject)subdmc, name));
2580   PetscCall(DMSwarmSetCellDM(cellswarm, subdmc));
2581   PetscCall(DMLabelDestroy(&label));
2582   PetscFunctionReturn(PETSC_SUCCESS);
2583 }
2584 
2585 /*@
2586   DMSwarmRestoreCellSwarm - Restores a `DMSWARM` object obtained with `DMSwarmGetCellSwarm()`. All fields are copied back into the parent swarm.
2587 
2588   Noncollective
2589 
2590   Input Parameters:
2591 + sw        - the parent `DMSWARM`
2592 . cellID    - the integer id of the cell to be copied back into the parent swarm
2593 - cellswarm - the cell swarm object
2594 
2595   Level: beginner
2596 
2597   Note:
2598   This only supports `DMSWARM_PIC` types of `DMSWARM`s
2599 
2600 .seealso: `DM`, `DMSWARM`, `DMSwarmGetCellSwarm()`
2601 @*/
2602 PetscErrorCode DMSwarmRestoreCellSwarm(DM sw, PetscInt cellID, DM cellswarm)
2603 {
2604   DM        dmc;
2605   PetscInt *pids, particles, p;
2606 
2607   PetscFunctionBegin;
2608   PetscCall(DMSwarmSortGetAccess(sw));
2609   PetscCall(DMSwarmSortGetPointsPerCell(sw, cellID, &particles, &pids));
2610   PetscCall(DMSwarmSortRestoreAccess(sw));
2611   /* Pointwise copy of each particle based on pid. The parent swarm may not be altered during this process. */
2612   for (p = 0; p < particles; ++p) PetscCall(DMSwarmDataBucketCopyPoint(((DM_Swarm *)cellswarm->data)->db, pids[p], ((DM_Swarm *)sw->data)->db, pids[p]));
2613   /* Free memory, destroy cell dm */
2614   PetscCall(DMSwarmGetCellDM(cellswarm, &dmc));
2615   PetscCall(DMDestroy(&dmc));
2616   PetscCall(DMSwarmSortRestorePointsPerCell(sw, cellID, &particles, &pids));
2617   PetscFunctionReturn(PETSC_SUCCESS);
2618 }
2619 
2620 /*@
2621   DMSwarmComputeMoments - Compute the first three particle moments for a given field
2622 
2623   Noncollective
2624 
2625   Input Parameters:
2626 + sw         - the `DMSWARM`
2627 . coordinate - the coordinate field name
2628 - weight     - the weight field name
2629 
2630   Output Parameter:
2631 . moments - the field moments
2632 
2633   Level: intermediate
2634 
2635   Notes:
2636   The `moments` array should be of length bs + 2, where bs is the block size of the coordinate field.
2637 
2638   The weight field must be a scalar, having blocksize 1.
2639 
2640 .seealso: `DM`, `DMSWARM`, `DMPlexComputeMoments()`
2641 @*/
2642 PetscErrorCode DMSwarmComputeMoments(DM sw, const char coordinate[], const char weight[], PetscReal moments[])
2643 {
2644   const PetscReal *coords;
2645   const PetscReal *w;
2646   PetscReal       *mom;
2647   PetscDataType    dtc, dtw;
2648   PetscInt         bsc, bsw, Np;
2649   MPI_Comm         comm;
2650 
2651   PetscFunctionBegin;
2652   PetscValidHeaderSpecific(sw, DM_CLASSID, 1);
2653   PetscAssertPointer(coordinate, 2);
2654   PetscAssertPointer(weight, 3);
2655   PetscAssertPointer(moments, 4);
2656   PetscCall(PetscObjectGetComm((PetscObject)sw, &comm));
2657   PetscCall(DMSwarmGetField(sw, coordinate, &bsc, &dtc, (void **)&coords));
2658   PetscCall(DMSwarmGetField(sw, weight, &bsw, &dtw, (void **)&w));
2659   PetscCheck(dtc == PETSC_REAL, comm, PETSC_ERR_ARG_WRONG, "Coordinate field %s must be real, not %s", coordinate, PetscDataTypes[dtc]);
2660   PetscCheck(dtw == PETSC_REAL, comm, PETSC_ERR_ARG_WRONG, "Weight field %s must be real, not %s", weight, PetscDataTypes[dtw]);
2661   PetscCheck(bsw == 1, comm, PETSC_ERR_ARG_WRONG, "Weight field %s must be a scalar, not blocksize %" PetscInt_FMT, weight, bsw);
2662   PetscCall(DMSwarmGetLocalSize(sw, &Np));
2663   PetscCall(DMGetWorkArray(sw, bsc + 2, MPIU_REAL, &mom));
2664   PetscCall(PetscArrayzero(mom, bsc + 2));
2665   for (PetscInt p = 0; p < Np; ++p) {
2666     const PetscReal *c  = &coords[p * bsc];
2667     const PetscReal  wp = w[p];
2668 
2669     mom[0] += wp;
2670     for (PetscInt d = 0; d < bsc; ++d) {
2671       mom[d + 1] += wp * c[d];
2672       mom[d + bsc + 1] += wp * PetscSqr(c[d]);
2673     }
2674   }
2675   PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&coords));
2676   PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&w));
2677   PetscCallMPI(MPIU_Allreduce(mom, moments, bsc + 2, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)sw)));
2678   PetscCall(DMRestoreWorkArray(sw, bsc + 2, MPIU_REAL, &mom));
2679   PetscFunctionReturn(PETSC_SUCCESS);
2680 }
2681 
2682 static PetscErrorCode DMSetFromOptions_Swarm(DM dm, PetscOptionItems PetscOptionsObject)
2683 {
2684   DM_Swarm *swarm = (DM_Swarm *)dm->data;
2685 
2686   PetscFunctionBegin;
2687   PetscOptionsHeadBegin(PetscOptionsObject, "DMSwarm Options");
2688   PetscCall(PetscOptionsEnum("-dm_swarm_remap_type", "Remap algorithm", "DMSwarmSetRemapType", DMSwarmRemapTypeNames, (PetscEnum)swarm->remap_type, (PetscEnum *)&swarm->remap_type, NULL));
2689   PetscOptionsHeadEnd();
2690   PetscFunctionReturn(PETSC_SUCCESS);
2691 }
2692 
2693 PETSC_INTERN PetscErrorCode DMClone_Swarm(DM, DM *);
2694 
2695 static PetscErrorCode DMInitialize_Swarm(DM sw)
2696 {
2697   PetscFunctionBegin;
2698   sw->ops->view                     = DMView_Swarm;
2699   sw->ops->load                     = NULL;
2700   sw->ops->setfromoptions           = DMSetFromOptions_Swarm;
2701   sw->ops->clone                    = DMClone_Swarm;
2702   sw->ops->setup                    = DMSetup_Swarm;
2703   sw->ops->createlocalsection       = NULL;
2704   sw->ops->createsectionpermutation = NULL;
2705   sw->ops->createdefaultconstraints = NULL;
2706   sw->ops->createglobalvector       = DMCreateGlobalVector_Swarm;
2707   sw->ops->createlocalvector        = DMCreateLocalVector_Swarm;
2708   sw->ops->getlocaltoglobalmapping  = NULL;
2709   sw->ops->createfieldis            = NULL;
2710   sw->ops->createcoordinatedm       = NULL;
2711   sw->ops->createcellcoordinatedm   = NULL;
2712   sw->ops->getcoloring              = NULL;
2713   sw->ops->creatematrix             = DMCreateMatrix_Swarm;
2714   sw->ops->createinterpolation      = NULL;
2715   sw->ops->createinjection          = NULL;
2716   sw->ops->createmassmatrix         = DMCreateMassMatrix_Swarm;
2717   sw->ops->creategradientmatrix     = DMCreateGradientMatrix_Swarm;
2718   sw->ops->refine                   = NULL;
2719   sw->ops->coarsen                  = NULL;
2720   sw->ops->refinehierarchy          = NULL;
2721   sw->ops->coarsenhierarchy         = NULL;
2722   sw->ops->globaltolocalbegin       = DMGlobalToLocalBegin_Swarm;
2723   sw->ops->globaltolocalend         = DMGlobalToLocalEnd_Swarm;
2724   sw->ops->localtoglobalbegin       = DMLocalToGlobalBegin_Swarm;
2725   sw->ops->localtoglobalend         = DMLocalToGlobalEnd_Swarm;
2726   sw->ops->destroy                  = DMDestroy_Swarm;
2727   sw->ops->createsubdm              = NULL;
2728   sw->ops->getdimpoints             = NULL;
2729   sw->ops->locatepoints             = NULL;
2730   sw->ops->projectfieldlocal        = DMProjectFieldLocal_Swarm;
2731   PetscFunctionReturn(PETSC_SUCCESS);
2732 }
2733 
2734 PETSC_INTERN PetscErrorCode DMClone_Swarm(DM dm, DM *newdm)
2735 {
2736   DM_Swarm *swarm = (DM_Swarm *)dm->data;
2737 
2738   PetscFunctionBegin;
2739   swarm->refct++;
2740   (*newdm)->data = swarm;
2741   PetscCall(PetscObjectChangeTypeName((PetscObject)*newdm, DMSWARM));
2742   PetscCall(DMInitialize_Swarm(*newdm));
2743   (*newdm)->dim = dm->dim;
2744   PetscFunctionReturn(PETSC_SUCCESS);
2745 }
2746 
2747 /*MC
2748  DMSWARM = "swarm" - A `DM` object for particle methods, such as particle-in-cell (PIC), in which the underlying
2749            data is both (i) dynamic in length, (ii) and of arbitrary data type.
2750 
2751  Level: intermediate
2752 
2753  Notes:
2754  User data can be represented by `DMSWARM` through a registering "fields" which are to be stored on particles.
2755  To register a field, the user must provide;
2756  (a) a unique name;
2757  (b) the data type (or size in bytes);
2758  (c) the block size of the data.
2759 
2760  For example, suppose the application requires a unique id, energy, momentum and density to be stored
2761  on a set of particles. Then the following code could be used
2762 .vb
2763     DMSwarmInitializeFieldRegister(dm)
2764     DMSwarmRegisterPetscDatatypeField(dm,"uid",1,PETSC_LONG);
2765     DMSwarmRegisterPetscDatatypeField(dm,"energy",1,PETSC_REAL);
2766     DMSwarmRegisterPetscDatatypeField(dm,"momentum",3,PETSC_REAL);
2767     DMSwarmRegisterPetscDatatypeField(dm,"density",1,PETSC_FLOAT);
2768     DMSwarmFinalizeFieldRegister(dm)
2769 .ve
2770 
2771  The fields represented by `DMSWARM` are dynamic and can be re-sized at any time.
2772  The only restriction imposed by `DMSWARM` is that all fields contain the same number of particles.
2773 
2774  To support particle methods, "migration" techniques are provided. These methods migrate data
2775  between ranks.
2776 
2777  `DMSWARM` supports the methods `DMCreateGlobalVector()` and `DMCreateLocalVector()`.
2778  As a `DMSWARM` may internally define and store values of different data types,
2779  before calling `DMCreateGlobalVector()` or `DMCreateLocalVector()`, the user must inform `DMSWARM` which
2780  fields should be used to define a `Vec` object via `DMSwarmVectorDefineField()`
2781  The specified field can be changed at any time - thereby permitting vectors
2782  compatible with different fields to be created.
2783 
2784  A dual representation of fields in the `DMSWARM` and a Vec object is permitted via `DMSwarmCreateGlobalVectorFromField()`
2785  Here the data defining the field in the `DMSWARM` is shared with a `Vec`.
2786  This is inherently unsafe if you alter the size of the field at any time between
2787  calls to `DMSwarmCreateGlobalVectorFromField()` and `DMSwarmDestroyGlobalVectorFromField()`.
2788  If the local size of the `DMSWARM` does not match the local size of the global vector
2789  when `DMSwarmDestroyGlobalVectorFromField()` is called, an error is thrown.
2790 
2791  Additional high-level support is provided for Particle-In-Cell methods. Refer to `DMSwarmSetType()`.
2792 
2793 .seealso: `DM`, `DMSWARM`, `DMType`, `DMCreate()`, `DMSetType()`, `DMSwarmSetType()`, `DMSwarmType`, `DMSwarmCreateGlobalVectorFromField()`,
2794          `DMCreateGlobalVector()`, `DMCreateLocalVector()`
2795 M*/
2796 
2797 PETSC_EXTERN PetscErrorCode DMCreate_Swarm(DM dm)
2798 {
2799   DM_Swarm *swarm;
2800 
2801   PetscFunctionBegin;
2802   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2803   PetscCall(PetscNew(&swarm));
2804   dm->data = swarm;
2805   PetscCall(DMSwarmDataBucketCreate(&swarm->db));
2806   PetscCall(DMSwarmInitializeFieldRegister(dm));
2807   dm->dim                               = 0;
2808   swarm->refct                          = 1;
2809   swarm->issetup                        = PETSC_FALSE;
2810   swarm->swarm_type                     = DMSWARM_BASIC;
2811   swarm->migrate_type                   = DMSWARM_MIGRATE_BASIC;
2812   swarm->collect_type                   = DMSWARM_COLLECT_BASIC;
2813   swarm->migrate_error_on_missing_point = PETSC_FALSE;
2814   swarm->collect_view_active            = PETSC_FALSE;
2815   swarm->collect_view_reset_nlocal      = -1;
2816   PetscCall(DMInitialize_Swarm(dm));
2817   if (SwarmDataFieldId == -1) PetscCall(PetscObjectComposedDataRegister(&SwarmDataFieldId));
2818   PetscFunctionReturn(PETSC_SUCCESS);
2819 }
2820 
2821 /* Replace dm with the contents of ndm, and then destroy ndm
2822    - Share the DM_Swarm structure
2823 */
2824 PetscErrorCode DMSwarmReplace(DM dm, DM *ndm)
2825 {
2826   DM               dmNew = *ndm;
2827   const PetscReal *maxCell, *Lstart, *L;
2828   PetscInt         dim;
2829 
2830   PetscFunctionBegin;
2831   if (dm == dmNew) {
2832     PetscCall(DMDestroy(ndm));
2833     PetscFunctionReturn(PETSC_SUCCESS);
2834   }
2835   dm->setupcalled = dmNew->setupcalled;
2836   if (!dm->hdr.name) {
2837     const char *name;
2838 
2839     PetscCall(PetscObjectGetName((PetscObject)*ndm, &name));
2840     PetscCall(PetscObjectSetName((PetscObject)dm, name));
2841   }
2842   PetscCall(DMGetDimension(dmNew, &dim));
2843   PetscCall(DMSetDimension(dm, dim));
2844   PetscCall(DMGetPeriodicity(dmNew, &maxCell, &Lstart, &L));
2845   PetscCall(DMSetPeriodicity(dm, maxCell, Lstart, L));
2846   PetscCall(DMDestroy_Swarm(dm));
2847   PetscCall(DMInitialize_Swarm(dm));
2848   dm->data = dmNew->data;
2849   ((DM_Swarm *)dmNew->data)->refct++;
2850   PetscCall(DMDestroy(ndm));
2851   PetscFunctionReturn(PETSC_SUCCESS);
2852 }
2853 
2854 /*@
2855   DMSwarmDuplicate - Creates a new `DMSWARM` with the same fields and cell `DM`s but no particles
2856 
2857   Collective
2858 
2859   Input Parameter:
2860 . sw - the `DMSWARM`
2861 
2862   Output Parameter:
2863 . nsw - the new `DMSWARM`
2864 
2865   Level: beginner
2866 
2867 .seealso: `DM`, `DMSWARM`, `DMSwarmCreate()`, `DMClone()`
2868 @*/
2869 PetscErrorCode DMSwarmDuplicate(DM sw, DM *nsw)
2870 {
2871   DM_Swarm         *swarm = (DM_Swarm *)sw->data;
2872   DMSwarmDataField *fields;
2873   DMSwarmCellDM     celldm, ncelldm;
2874   DMSwarmType       stype;
2875   const char       *name, **celldmnames;
2876   void             *ctx;
2877   PetscInt          dim, Nf, Ndm;
2878   PetscBool         flg;
2879 
2880   PetscFunctionBegin;
2881   PetscCall(DMCreate(PetscObjectComm((PetscObject)sw), nsw));
2882   PetscCall(DMSetType(*nsw, DMSWARM));
2883   PetscCall(PetscObjectGetName((PetscObject)sw, &name));
2884   PetscCall(PetscObjectSetName((PetscObject)*nsw, name));
2885   PetscCall(DMGetDimension(sw, &dim));
2886   PetscCall(DMSetDimension(*nsw, dim));
2887   PetscCall(DMSwarmGetType(sw, &stype));
2888   PetscCall(DMSwarmSetType(*nsw, stype));
2889   PetscCall(DMGetApplicationContext(sw, &ctx));
2890   PetscCall(DMSetApplicationContext(*nsw, ctx));
2891 
2892   PetscCall(DMSwarmDataBucketGetDMSwarmDataFields(swarm->db, &Nf, &fields));
2893   for (PetscInt f = 0; f < Nf; ++f) {
2894     PetscCall(DMSwarmDataFieldStringInList(fields[f]->name, ((DM_Swarm *)(*nsw)->data)->db->nfields, (const DMSwarmDataField *)((DM_Swarm *)(*nsw)->data)->db->field, &flg));
2895     if (!flg) PetscCall(DMSwarmRegisterPetscDatatypeField(*nsw, fields[f]->name, fields[f]->bs, fields[f]->petsc_type));
2896   }
2897 
2898   PetscCall(DMSwarmGetCellDMNames(sw, &Ndm, &celldmnames));
2899   for (PetscInt c = 0; c < Ndm; ++c) {
2900     DM           dm;
2901     PetscInt     Ncf;
2902     const char **coordfields, **fields;
2903 
2904     PetscCall(DMSwarmGetCellDMByName(sw, celldmnames[c], &celldm));
2905     PetscCall(DMSwarmCellDMGetDM(celldm, &dm));
2906     PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Ncf, &coordfields));
2907     PetscCall(DMSwarmCellDMGetFields(celldm, &Nf, &fields));
2908     PetscCall(DMSwarmCellDMCreate(dm, Nf, fields, Ncf, coordfields, &ncelldm));
2909     PetscCall(DMSwarmAddCellDM(*nsw, ncelldm));
2910     PetscCall(DMSwarmCellDMDestroy(&ncelldm));
2911   }
2912   PetscCall(PetscFree(celldmnames));
2913 
2914   PetscCall(DMSetFromOptions(*nsw));
2915   PetscCall(DMSetUp(*nsw));
2916   PetscCall(DMSwarmGetCellDMActive(sw, &celldm));
2917   PetscCall(PetscObjectGetName((PetscObject)celldm, &name));
2918   PetscCall(DMSwarmSetCellDMActive(*nsw, name));
2919   PetscFunctionReturn(PETSC_SUCCESS);
2920 }
2921 
2922 PetscErrorCode DMLocalToGlobalBegin_Swarm(DM dm, Vec l, InsertMode mode, Vec g)
2923 {
2924   PetscFunctionBegin;
2925   PetscFunctionReturn(PETSC_SUCCESS);
2926 }
2927 
2928 PetscErrorCode DMLocalToGlobalEnd_Swarm(DM dm, Vec l, InsertMode mode, Vec g)
2929 {
2930   PetscFunctionBegin;
2931   switch (mode) {
2932   case INSERT_VALUES:
2933     PetscCall(VecCopy(l, g));
2934     break;
2935   case ADD_VALUES:
2936     PetscCall(VecAXPY(g, 1., l));
2937     break;
2938   default:
2939     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mode not supported: %d", mode);
2940   }
2941   PetscFunctionReturn(PETSC_SUCCESS);
2942 }
2943 
2944 PetscErrorCode DMGlobalToLocalBegin_Swarm(DM dm, Vec g, InsertMode mode, Vec l)
2945 {
2946   PetscFunctionBegin;
2947   PetscFunctionReturn(PETSC_SUCCESS);
2948 }
2949 
2950 PetscErrorCode DMGlobalToLocalEnd_Swarm(DM dm, Vec g, InsertMode mode, Vec l)
2951 {
2952   PetscFunctionBegin;
2953   PetscCall(DMLocalToGlobalEnd_Swarm(dm, g, mode, l));
2954   PetscFunctionReturn(PETSC_SUCCESS);
2955 }
2956