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