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