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