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