xref: /petsc/src/dm/impls/swarm/swarm.c (revision f2ed2dc71a2ab9ffda85eae8afa0cbea9ed570de)
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 "../src/dm/impls/swarm/data_bucket.h"
9 
10 PetscLogEvent DMSWARM_Migrate, DMSWARM_SetSizes, DMSWARM_AddPoints, DMSWARM_RemovePoints, DMSWARM_Sort;
11 PetscLogEvent DMSWARM_DataExchangerTopologySetup, DMSWARM_DataExchangerBegin, DMSWARM_DataExchangerEnd;
12 PetscLogEvent DMSWARM_DataExchangerSendCount, DMSWARM_DataExchangerPack;
13 
14 const char* DMSwarmTypeNames[] = { "basic", "pic", 0 };
15 const char* DMSwarmMigrateTypeNames[] = { "basic", "dmcellnscatter", "dmcellexact", "user", 0 };
16 const char* DMSwarmCollectTypeNames[] = { "basic", "boundingbox", "general", "user", 0 };
17 const char* DMSwarmPICLayoutTypeNames[] = { "regular", "gauss", "subdivision", 0 };
18 
19 const char DMSwarmField_pid[] = "DMSwarm_pid";
20 const char DMSwarmField_rank[] = "DMSwarm_rank";
21 const char DMSwarmPICField_coor[] = "DMSwarmPIC_coor";
22 const char DMSwarmPICField_cellid[] = "DMSwarm_cellid";
23 
24 PETSC_EXTERN PetscErrorCode VecView_Seq(Vec, PetscViewer);
25 PETSC_EXTERN PetscErrorCode VecView_MPI(Vec, PetscViewer);
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   PetscErrorCode ierr;
37 
38   PetscFunctionBegin;
39   ierr = VecGetDM(v, &dm);CHKERRQ(ierr);
40   ierr = VecGetBlockSize(v, &bs);CHKERRQ(ierr);
41   ierr = PetscViewerHDF5PushGroup(viewer, "/particle_fields");CHKERRQ(ierr);
42   ierr = PetscObjectTypeCompare((PetscObject) v, VECSEQ, &isseq);CHKERRQ(ierr);
43   ierr = DMGetOutputSequenceNumber(dm, &seqnum, &seqval);CHKERRQ(ierr);
44   ierr = PetscViewerHDF5SetTimestep(viewer, seqnum);CHKERRQ(ierr);
45   /* ierr = DMSequenceView_HDF5(dm, "time", seqnum, (PetscScalar) seqval, viewer);CHKERRQ(ierr); */
46   if (isseq) {ierr = VecView_Seq(v, viewer);CHKERRQ(ierr);}
47   else       {ierr = VecView_MPI(v, viewer);CHKERRQ(ierr);}
48   ierr = PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject) v, "Nc", PETSC_INT, (void *) &bs);CHKERRQ(ierr);
49   ierr = PetscViewerHDF5PopGroup(viewer);CHKERRQ(ierr);
50   PetscFunctionReturn(0);
51 }
52 
53 PetscErrorCode DMSwarmView_HDF5(DM dm, PetscViewer viewer)
54 {
55   Vec            coordinates;
56   PetscInt       Np;
57   PetscBool      isseq;
58   PetscErrorCode ierr;
59 
60   PetscFunctionBegin;
61   ierr = DMSwarmGetSize(dm, &Np);CHKERRQ(ierr);
62   ierr = DMSwarmCreateGlobalVectorFromField(dm, DMSwarmPICField_coor, &coordinates);CHKERRQ(ierr);
63   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
64   ierr = PetscViewerHDF5PushGroup(viewer, "/particles");CHKERRQ(ierr);
65   ierr = PetscObjectTypeCompare((PetscObject) coordinates, VECSEQ, &isseq);CHKERRQ(ierr);
66   if (isseq) {ierr = VecView_Seq(coordinates, viewer);CHKERRQ(ierr);}
67   else       {ierr = VecView_MPI(coordinates, viewer);CHKERRQ(ierr);}
68   ierr = PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject) coordinates, "Np", PETSC_INT, (void *) &Np);CHKERRQ(ierr);
69   ierr = PetscViewerHDF5PopGroup(viewer);CHKERRQ(ierr);
70   ierr = DMSwarmDestroyGlobalVectorFromField(dm, DMSwarmPICField_coor, &coordinates);CHKERRQ(ierr);
71   PetscFunctionReturn(0);
72 }
73 #endif
74 
75 PetscErrorCode VecView_Swarm(Vec v, PetscViewer viewer)
76 {
77   DM             dm;
78   PetscBool      ishdf5;
79   PetscErrorCode ierr;
80 
81   PetscFunctionBegin;
82   ierr = VecGetDM(v, &dm);CHKERRQ(ierr);
83   if (!dm) SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
84   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5,  &ishdf5);CHKERRQ(ierr);
85   if (ishdf5) {
86 #if defined(PETSC_HAVE_HDF5)
87     ierr = VecView_Swarm_HDF5_Internal(v, viewer);CHKERRQ(ierr);
88 #else
89     SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
90 #endif
91   } else {
92     PetscBool isseq;
93 
94     ierr = PetscObjectTypeCompare((PetscObject) v, VECSEQ, &isseq);CHKERRQ(ierr);
95     if (isseq) {ierr = VecView_Seq(v, viewer);CHKERRQ(ierr);}
96     else       {ierr = VecView_MPI(v, viewer);CHKERRQ(ierr);}
97   }
98   PetscFunctionReturn(0);
99 }
100 
101 /*@C
102    DMSwarmVectorDefineField - Sets the field from which to define a Vec object
103                              when DMCreateLocalVector(), or DMCreateGlobalVector() is called
104 
105    Collective on dm
106 
107    Input parameters:
108 +  dm - a DMSwarm
109 -  fieldname - the textual name given to a registered field
110 
111    Level: beginner
112 
113    Notes:
114 
115    The field with name fieldname must be defined as having a data type of PetscScalar.
116 
117    This function must be called prior to calling DMCreateLocalVector(), DMCreateGlobalVector().
118    Mutiple calls to DMSwarmVectorDefineField() are permitted.
119 
120 .seealso: DMSwarmRegisterPetscDatatypeField(), DMCreateGlobalVector(), DMCreateLocalVector()
121 @*/
122 PetscErrorCode DMSwarmVectorDefineField(DM dm,const char fieldname[])
123 {
124   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
125   PetscErrorCode ierr;
126   PetscInt       bs,n;
127   PetscScalar    *array;
128   PetscDataType  type;
129 
130   PetscFunctionBegin;
131   if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); }
132   ierr = DMSwarmDataBucketGetSizes(swarm->db,&n,NULL,NULL);CHKERRQ(ierr);
133   ierr = DMSwarmGetField(dm,fieldname,&bs,&type,(void**)&array);CHKERRQ(ierr);
134 
135   /* Check all fields are of type PETSC_REAL or PETSC_SCALAR */
136   if (type != PETSC_REAL) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid for PETSC_REAL");
137   ierr = PetscSNPrintf(swarm->vec_field_name,PETSC_MAX_PATH_LEN-1,"%s",fieldname);CHKERRQ(ierr);
138   swarm->vec_field_set = PETSC_TRUE;
139   swarm->vec_field_bs = bs;
140   swarm->vec_field_nlocal = n;
141   ierr = DMSwarmRestoreField(dm,fieldname,&bs,&type,(void**)&array);CHKERRQ(ierr);
142   PetscFunctionReturn(0);
143 }
144 
145 /* requires DMSwarmDefineFieldVector has been called */
146 PetscErrorCode DMCreateGlobalVector_Swarm(DM dm,Vec *vec)
147 {
148   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
149   PetscErrorCode ierr;
150   Vec            x;
151   char           name[PETSC_MAX_PATH_LEN];
152 
153   PetscFunctionBegin;
154   if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); }
155   if (!swarm->vec_field_set) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmVectorDefineField first");
156   if (swarm->vec_field_nlocal != swarm->db->L) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSwarm sizes have changed since last call to VectorDefineField first"); /* Stale data */
157 
158   ierr = PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarmField_%s",swarm->vec_field_name);CHKERRQ(ierr);
159   ierr = VecCreate(PetscObjectComm((PetscObject)dm),&x);CHKERRQ(ierr);
160   ierr = PetscObjectSetName((PetscObject)x,name);CHKERRQ(ierr);
161   ierr = VecSetSizes(x,swarm->db->L*swarm->vec_field_bs,PETSC_DETERMINE);CHKERRQ(ierr);
162   ierr = VecSetBlockSize(x,swarm->vec_field_bs);CHKERRQ(ierr);
163   ierr = VecSetDM(x,dm);CHKERRQ(ierr);
164   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
165   ierr = VecSetDM(x, dm);CHKERRQ(ierr);
166   ierr = VecSetOperation(x, VECOP_VIEW, (void (*)(void)) VecView_Swarm);CHKERRQ(ierr);
167   *vec = x;
168   PetscFunctionReturn(0);
169 }
170 
171 /* requires DMSwarmDefineFieldVector has been called */
172 PetscErrorCode DMCreateLocalVector_Swarm(DM dm,Vec *vec)
173 {
174   DM_Swarm *swarm = (DM_Swarm*)dm->data;
175   PetscErrorCode ierr;
176   Vec x;
177   char name[PETSC_MAX_PATH_LEN];
178 
179   PetscFunctionBegin;
180   if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); }
181   if (!swarm->vec_field_set) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmVectorDefineField first");
182   if (swarm->vec_field_nlocal != swarm->db->L) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSwarm sizes have changed since last call to VectorDefineField first"); /* Stale data */
183 
184   ierr = PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarmField_%s",swarm->vec_field_name);CHKERRQ(ierr);
185   ierr = VecCreate(PETSC_COMM_SELF,&x);CHKERRQ(ierr);
186   ierr = PetscObjectSetName((PetscObject)x,name);CHKERRQ(ierr);
187   ierr = VecSetSizes(x,swarm->db->L*swarm->vec_field_bs,PETSC_DETERMINE);CHKERRQ(ierr);
188   ierr = VecSetBlockSize(x,swarm->vec_field_bs);CHKERRQ(ierr);
189   ierr = VecSetDM(x,dm);CHKERRQ(ierr);
190   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
191   *vec = x;
192   PetscFunctionReturn(0);
193 }
194 
195 static PetscErrorCode DMSwarmDestroyVectorFromField_Private(DM dm, const char fieldname[], Vec *vec)
196 {
197   DM_Swarm      *swarm = (DM_Swarm *) dm->data;
198   DMSwarmDataField      gfield;
199   void         (*fptr)(void);
200   PetscInt       bs, nlocal;
201   char           name[PETSC_MAX_PATH_LEN];
202   PetscErrorCode ierr;
203 
204   PetscFunctionBegin;
205   ierr = VecGetLocalSize(*vec, &nlocal);CHKERRQ(ierr);
206   ierr = VecGetBlockSize(*vec, &bs);CHKERRQ(ierr);
207   if (nlocal/bs != swarm->db->L) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSwarm sizes have changed since vector was created - cannot ensure pointers are valid"); /* Stale data */
208   ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield);CHKERRQ(ierr);
209   /* check vector is an inplace array */
210   ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "DMSwarm_VecFieldInPlace_%s", fieldname);CHKERRQ(ierr);
211   ierr = PetscObjectQueryFunction((PetscObject) *vec, name, &fptr);CHKERRQ(ierr);
212   if (!fptr) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_USER, "Vector being destroyed was not created from DMSwarm field(%s)", fieldname);
213   ierr = DMSwarmDataFieldRestoreAccess(gfield);CHKERRQ(ierr);
214   ierr = VecDestroy(vec);CHKERRQ(ierr);
215   PetscFunctionReturn(0);
216 }
217 
218 static PetscErrorCode DMSwarmCreateVectorFromField_Private(DM dm, const char fieldname[], MPI_Comm comm, Vec *vec)
219 {
220   DM_Swarm      *swarm = (DM_Swarm *) dm->data;
221   PetscDataType  type;
222   PetscScalar   *array;
223   PetscInt       bs, n;
224   char           name[PETSC_MAX_PATH_LEN];
225   PetscMPIInt    size;
226   PetscErrorCode ierr;
227 
228   PetscFunctionBegin;
229   if (!swarm->issetup) {ierr = DMSetUp(dm);CHKERRQ(ierr);}
230   ierr = DMSwarmDataBucketGetSizes(swarm->db, &n, NULL, NULL);CHKERRQ(ierr);
231   ierr = DMSwarmGetField(dm, fieldname, &bs, &type, (void **) &array);CHKERRQ(ierr);
232   if (type != PETSC_REAL) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Only valid for PETSC_REAL");
233 
234   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
235   if (size == 1) {
236     ierr = VecCreateSeqWithArray(comm, bs, n*bs, array, vec);CHKERRQ(ierr);
237   } else {
238     ierr = VecCreateMPIWithArray(comm, bs, n*bs, PETSC_DETERMINE, array, vec);CHKERRQ(ierr);
239   }
240   ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "DMSwarmSharedField_%s", fieldname);CHKERRQ(ierr);
241   ierr = PetscObjectSetName((PetscObject) *vec, name);CHKERRQ(ierr);
242 
243   /* Set guard */
244   ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "DMSwarm_VecFieldInPlace_%s", fieldname);CHKERRQ(ierr);
245   ierr = PetscObjectComposeFunction((PetscObject) *vec, name, DMSwarmDestroyVectorFromField_Private);CHKERRQ(ierr);
246 
247   ierr = VecSetDM(*vec, dm);CHKERRQ(ierr);
248   ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Swarm);CHKERRQ(ierr);
249   PetscFunctionReturn(0);
250 }
251 
252 /* This creates a "mass matrix" between a finite element and particle space. If a finite element interpolant is given by
253 
254      \hat f = \sum_i f_i \phi_i
255 
256    and a particle function is given by
257 
258      f = \sum_i w_i \delta(x - x_i)
259 
260    then we want to require that
261 
262      M \hat f = M_p f
263 
264    where the particle mass matrix is given by
265 
266      (M_p)_{ij} = \int \phi_i \delta(x - x_j)
267 
268    The way Dave May does particles, they amount to quadratue weights rather than delta functions, so he has |J| is in
269    his integral. We allow this with the boolean flag.
270 */
271 static PetscErrorCode DMSwarmComputeMassMatrix_Private(DM dmc, DM dmf, Mat mass, PetscBool useDeltaFunction, void *ctx)
272 {
273   const char    *name = "Mass Matrix";
274   MPI_Comm       comm;
275   PetscDS        prob;
276   PetscSection   fsection, globalFSection;
277   PetscHSetIJ    ht;
278   PetscLayout    rLayout, colLayout;
279   PetscInt      *dnz, *onz;
280   PetscInt       locRows, locCols, rStart, colStart, colEnd, *rowIDXs;
281   PetscReal     *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0};
282   PetscScalar   *elemMat;
283   PetscInt       dim, Nf, field, cStart, cEnd, cell, totDim, maxC = 0;
284   PetscErrorCode ierr;
285 
286   PetscFunctionBegin;
287   ierr = PetscObjectGetComm((PetscObject) mass, &comm);CHKERRQ(ierr);
288   ierr = DMGetCoordinateDim(dmf, &dim);CHKERRQ(ierr);
289   ierr = DMGetDS(dmf, &prob);CHKERRQ(ierr);
290   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
291   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
292   ierr = PetscMalloc3(dim, &v0, dim*dim, &J, dim*dim,&invJ);CHKERRQ(ierr);
293   ierr = DMGetLocalSection(dmf, &fsection);CHKERRQ(ierr);
294   ierr = DMGetGlobalSection(dmf, &globalFSection);CHKERRQ(ierr);
295   ierr = DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd);CHKERRQ(ierr);
296   ierr = MatGetLocalSize(mass, &locRows, &locCols);CHKERRQ(ierr);
297 
298   ierr = PetscLayoutCreate(comm, &colLayout);CHKERRQ(ierr);
299   ierr = PetscLayoutSetLocalSize(colLayout, locCols);CHKERRQ(ierr);
300   ierr = PetscLayoutSetBlockSize(colLayout, 1);CHKERRQ(ierr);
301   ierr = PetscLayoutSetUp(colLayout);CHKERRQ(ierr);
302   ierr = PetscLayoutGetRange(colLayout, &colStart, &colEnd);CHKERRQ(ierr);
303   ierr = PetscLayoutDestroy(&colLayout);CHKERRQ(ierr);
304 
305   ierr = PetscLayoutCreate(comm, &rLayout);CHKERRQ(ierr);
306   ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr);
307   ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr);
308   ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr);
309   ierr = PetscLayoutGetRange(rLayout, &rStart, NULL);CHKERRQ(ierr);
310   ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr);
311 
312   ierr = PetscCalloc2(locRows, &dnz, locRows, &onz);CHKERRQ(ierr);
313   ierr = PetscHSetIJCreate(&ht);CHKERRQ(ierr);
314 
315   ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr);
316   /* count non-zeros */
317   ierr = DMSwarmSortGetAccess(dmc);CHKERRQ(ierr);
318   for (field = 0; field < Nf; ++field) {
319     for (cell = cStart; cell < cEnd; ++cell) {
320       PetscInt  c, i;
321       PetscInt *findices,   *cindices; /* fine is vertices, coarse is particles */
322       PetscInt  numFIndices, numCIndices;
323 
324       ierr = DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL);CHKERRQ(ierr);
325       ierr = DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices);CHKERRQ(ierr);
326       maxC = PetscMax(maxC, numCIndices);
327       {
328         PetscHashIJKey key;
329         PetscBool      missing;
330         for (i = 0; i < numFIndices; ++i) {
331           key.j = findices[i]; /* global column (from Plex) */
332           if (key.j >= 0) {
333             /* Get indices for coarse elements */
334             for (c = 0; c < numCIndices; ++c) {
335               key.i = cindices[c] + rStart; /* global cols (from Swarm) */
336               if (key.i < 0) continue;
337               ierr = PetscHSetIJQueryAdd(ht, key, &missing);CHKERRQ(ierr);
338               if (missing) {
339                 if ((key.j >= colStart) && (key.j < colEnd)) ++dnz[key.i - rStart];
340                 else                                         ++onz[key.i - rStart];
341               } else SETERRQ2(PetscObjectComm((PetscObject) dmf), PETSC_ERR_SUP, "Set new value at %D,%D", key.i, key.j);
342             }
343           }
344         }
345         ierr = PetscFree(cindices);CHKERRQ(ierr);
346       }
347       ierr = DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL);CHKERRQ(ierr);
348     }
349   }
350   ierr = PetscHSetIJDestroy(&ht);CHKERRQ(ierr);
351   ierr = MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL);CHKERRQ(ierr);
352   ierr = MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE);CHKERRQ(ierr);
353   ierr = PetscFree2(dnz, onz);CHKERRQ(ierr);
354   ierr = PetscMalloc3(maxC*totDim, &elemMat, maxC, &rowIDXs, maxC*dim, &xi);CHKERRQ(ierr);
355   for (field = 0; field < Nf; ++field) {
356     PetscTabulation Tcoarse;
357     PetscObject       obj;
358     PetscReal        *coords;
359     PetscInt          Nc, i;
360 
361     ierr = PetscDSGetDiscretization(prob, field, &obj);CHKERRQ(ierr);
362     ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);
363     if (Nc != 1) SETERRQ1(PetscObjectComm((PetscObject) dmf), PETSC_ERR_SUP, "Can only interpolate a scalar field from particles, Nc = %D", Nc);
364     ierr = DMSwarmGetField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
365     for (cell = cStart; cell < cEnd; ++cell) {
366       PetscInt *findices  , *cindices;
367       PetscInt  numFIndices, numCIndices;
368       PetscInt  p, c;
369 
370       /* TODO: Use DMField instead of assuming affine */
371       ierr = DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
372       ierr = DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL);CHKERRQ(ierr);
373       ierr = DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices);CHKERRQ(ierr);
374       for (p = 0; p < numCIndices; ++p) {
375         CoordinatesRealToRef(dim, dim, v0ref, v0, invJ, &coords[cindices[p]*dim], &xi[p*dim]);
376       }
377       ierr = PetscFECreateTabulation((PetscFE) obj, 1, numCIndices, xi, 0, &Tcoarse);CHKERRQ(ierr);
378       /* Get elemMat entries by multiplying by weight */
379       ierr = PetscArrayzero(elemMat, numCIndices*totDim);CHKERRQ(ierr);
380       for (i = 0; i < numFIndices; ++i) {
381         for (p = 0; p < numCIndices; ++p) {
382           for (c = 0; c < Nc; ++c) {
383             /* B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c */
384             elemMat[p*numFIndices+i] += Tcoarse->T[0][(p*numFIndices + i)*Nc + c]*(useDeltaFunction ? 1.0 : detJ);
385           }
386         }
387       }
388       for (p = 0; p < numCIndices; ++p) rowIDXs[p] = cindices[p] + rStart;
389       if (0) {ierr = DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat);CHKERRQ(ierr);}
390       ierr = MatSetValues(mass, numCIndices, rowIDXs, numFIndices, findices, elemMat, ADD_VALUES);CHKERRQ(ierr);
391       ierr = PetscFree(cindices);CHKERRQ(ierr);
392       ierr = DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL);CHKERRQ(ierr);
393       ierr = PetscTabulationDestroy(&Tcoarse);CHKERRQ(ierr);
394     }
395     ierr = DMSwarmRestoreField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
396   }
397   ierr = PetscFree3(elemMat, rowIDXs, xi);CHKERRQ(ierr);
398   ierr = DMSwarmSortRestoreAccess(dmc);CHKERRQ(ierr);
399   ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr);
400   ierr = MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
401   ierr = MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
402   PetscFunctionReturn(0);
403 }
404 
405 /* FEM cols, Particle rows */
406 static PetscErrorCode DMCreateMassMatrix_Swarm(DM dmCoarse, DM dmFine, Mat *mass)
407 {
408   PetscSection   gsf;
409   PetscInt       m, n;
410   void          *ctx;
411   PetscErrorCode ierr;
412 
413   PetscFunctionBegin;
414   ierr = DMGetGlobalSection(dmFine, &gsf);CHKERRQ(ierr);
415   ierr = PetscSectionGetConstrainedStorageSize(gsf, &m);CHKERRQ(ierr);
416   ierr = DMSwarmGetLocalSize(dmCoarse, &n);CHKERRQ(ierr);
417   ierr = MatCreate(PetscObjectComm((PetscObject) dmCoarse), mass);CHKERRQ(ierr);
418   ierr = MatSetSizes(*mass, n, m, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
419   ierr = MatSetType(*mass, dmCoarse->mattype);CHKERRQ(ierr);
420   ierr = DMGetApplicationContext(dmFine, &ctx);CHKERRQ(ierr);
421 
422   ierr = DMSwarmComputeMassMatrix_Private(dmCoarse, dmFine, *mass, PETSC_TRUE, ctx);CHKERRQ(ierr);
423   ierr = MatViewFromOptions(*mass, NULL, "-mass_mat_view");CHKERRQ(ierr);
424   PetscFunctionReturn(0);
425 }
426 
427 /*@C
428    DMSwarmCreateGlobalVectorFromField - Creates a Vec object sharing the array associated with a given field
429 
430    Collective on dm
431 
432    Input parameters:
433 +  dm - a DMSwarm
434 -  fieldname - the textual name given to a registered field
435 
436    Output parameter:
437 .  vec - the vector
438 
439    Level: beginner
440 
441    Notes:
442    The vector must be returned using a matching call to DMSwarmDestroyGlobalVectorFromField().
443 
444 .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmDestroyGlobalVectorFromField()
445 @*/
446 PetscErrorCode DMSwarmCreateGlobalVectorFromField(DM dm,const char fieldname[],Vec *vec)
447 {
448   MPI_Comm       comm = PetscObjectComm((PetscObject) dm);
449   PetscErrorCode ierr;
450 
451   PetscFunctionBegin;
452   ierr = DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec);CHKERRQ(ierr);
453   PetscFunctionReturn(0);
454 }
455 
456 /*@C
457    DMSwarmDestroyGlobalVectorFromField - Destroys the Vec object which share the array associated with a given field
458 
459    Collective on dm
460 
461    Input parameters:
462 +  dm - a DMSwarm
463 -  fieldname - the textual name given to a registered field
464 
465    Output parameter:
466 .  vec - the vector
467 
468    Level: beginner
469 
470 .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmCreateGlobalVectorFromField()
471 @*/
472 PetscErrorCode DMSwarmDestroyGlobalVectorFromField(DM dm,const char fieldname[],Vec *vec)
473 {
474   PetscErrorCode ierr;
475 
476   PetscFunctionBegin;
477   ierr = DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec);CHKERRQ(ierr);
478   PetscFunctionReturn(0);
479 }
480 
481 /*@C
482    DMSwarmCreateLocalVectorFromField - Creates a Vec object sharing the array associated with a given field
483 
484    Collective on dm
485 
486    Input parameters:
487 +  dm - a DMSwarm
488 -  fieldname - the textual name given to a registered field
489 
490    Output parameter:
491 .  vec - the vector
492 
493    Level: beginner
494 
495    Notes:
496    The vector must be returned using a matching call to DMSwarmDestroyLocalVectorFromField().
497 
498 .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmDestroyLocalVectorFromField()
499 @*/
500 PetscErrorCode DMSwarmCreateLocalVectorFromField(DM dm,const char fieldname[],Vec *vec)
501 {
502   MPI_Comm       comm = PETSC_COMM_SELF;
503   PetscErrorCode ierr;
504 
505   PetscFunctionBegin;
506   ierr = DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec);CHKERRQ(ierr);
507   PetscFunctionReturn(0);
508 }
509 
510 /*@C
511    DMSwarmDestroyLocalVectorFromField - Destroys the Vec object which share the array associated with a given field
512 
513    Collective on dm
514 
515    Input parameters:
516 +  dm - a DMSwarm
517 -  fieldname - the textual name given to a registered field
518 
519    Output parameter:
520 .  vec - the vector
521 
522    Level: beginner
523 
524 .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmCreateLocalVectorFromField()
525 @*/
526 PetscErrorCode DMSwarmDestroyLocalVectorFromField(DM dm,const char fieldname[],Vec *vec)
527 {
528   PetscErrorCode ierr;
529 
530   PetscFunctionBegin;
531   ierr = DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec);CHKERRQ(ierr);
532   PetscFunctionReturn(0);
533 }
534 
535 /*
536 PetscErrorCode DMSwarmCreateGlobalVectorFromFields(DM dm,const PetscInt nf,const char *fieldnames[],Vec *vec)
537 {
538   PetscFunctionReturn(0);
539 }
540 
541 PetscErrorCode DMSwarmRestoreGlobalVectorFromFields(DM dm,Vec *vec)
542 {
543   PetscFunctionReturn(0);
544 }
545 */
546 
547 /*@C
548    DMSwarmInitializeFieldRegister - Initiates the registration of fields to a DMSwarm
549 
550    Collective on dm
551 
552    Input parameter:
553 .  dm - a DMSwarm
554 
555    Level: beginner
556 
557    Notes:
558    After all fields have been registered, you must call DMSwarmFinalizeFieldRegister().
559 
560 .seealso: DMSwarmFinalizeFieldRegister(), DMSwarmRegisterPetscDatatypeField(),
561  DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField()
562 @*/
563 PetscErrorCode DMSwarmInitializeFieldRegister(DM dm)
564 {
565   DM_Swarm      *swarm = (DM_Swarm *) dm->data;
566   PetscErrorCode ierr;
567 
568   PetscFunctionBegin;
569   if (!swarm->field_registration_initialized) {
570     swarm->field_registration_initialized = PETSC_TRUE;
571     ierr = DMSwarmRegisterPetscDatatypeField(dm,DMSwarmField_pid,1,PETSC_INT64);CHKERRQ(ierr); /* unique identifer */
572     ierr = DMSwarmRegisterPetscDatatypeField(dm,DMSwarmField_rank,1,PETSC_INT);CHKERRQ(ierr); /* used for communication */
573   }
574   PetscFunctionReturn(0);
575 }
576 
577 /*@
578    DMSwarmFinalizeFieldRegister - Finalizes the registration of fields to a DMSwarm
579 
580    Collective on dm
581 
582    Input parameter:
583 .  dm - a DMSwarm
584 
585    Level: beginner
586 
587    Notes:
588    After DMSwarmFinalizeFieldRegister() has been called, no new fields can be defined on the DMSwarm.
589 
590 .seealso: DMSwarmInitializeFieldRegister(), DMSwarmRegisterPetscDatatypeField(),
591  DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField()
592 @*/
593 PetscErrorCode DMSwarmFinalizeFieldRegister(DM dm)
594 {
595   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
596   PetscErrorCode ierr;
597 
598   PetscFunctionBegin;
599   if (!swarm->field_registration_finalized) {
600     ierr = DMSwarmDataBucketFinalize(swarm->db);CHKERRQ(ierr);
601   }
602   swarm->field_registration_finalized = PETSC_TRUE;
603   PetscFunctionReturn(0);
604 }
605 
606 /*@
607    DMSwarmSetLocalSizes - Sets the length of all registered fields on the DMSwarm
608 
609    Not collective
610 
611    Input parameters:
612 +  dm - a DMSwarm
613 .  nlocal - the length of each registered field
614 -  buffer - the length of the buffer used to efficient dynamic re-sizing
615 
616    Level: beginner
617 
618 .seealso: DMSwarmGetLocalSize()
619 @*/
620 PetscErrorCode DMSwarmSetLocalSizes(DM dm,PetscInt nlocal,PetscInt buffer)
621 {
622   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
623   PetscErrorCode ierr;
624 
625   PetscFunctionBegin;
626   ierr = PetscLogEventBegin(DMSWARM_SetSizes,0,0,0,0);CHKERRQ(ierr);
627   ierr = DMSwarmDataBucketSetSizes(swarm->db,nlocal,buffer);CHKERRQ(ierr);
628   ierr = PetscLogEventEnd(DMSWARM_SetSizes,0,0,0,0);CHKERRQ(ierr);
629   PetscFunctionReturn(0);
630 }
631 
632 /*@
633    DMSwarmSetCellDM - Attachs a DM to a DMSwarm
634 
635    Collective on dm
636 
637    Input parameters:
638 +  dm - a DMSwarm
639 -  dmcell - the DM to attach to the DMSwarm
640 
641    Level: beginner
642 
643    Notes:
644    The attached DM (dmcell) will be queried for point location and
645    neighbor MPI-rank information if DMSwarmMigrate() is called.
646 
647 .seealso: DMSwarmSetType(), DMSwarmGetCellDM(), DMSwarmMigrate()
648 @*/
649 PetscErrorCode DMSwarmSetCellDM(DM dm,DM dmcell)
650 {
651   DM_Swarm *swarm = (DM_Swarm*)dm->data;
652 
653   PetscFunctionBegin;
654   swarm->dmcell = dmcell;
655   PetscFunctionReturn(0);
656 }
657 
658 /*@
659    DMSwarmGetCellDM - Fetches the attached cell DM
660 
661    Collective on dm
662 
663    Input parameter:
664 .  dm - a DMSwarm
665 
666    Output parameter:
667 .  dmcell - the DM which was attached to the DMSwarm
668 
669    Level: beginner
670 
671 .seealso: DMSwarmSetCellDM()
672 @*/
673 PetscErrorCode DMSwarmGetCellDM(DM dm,DM *dmcell)
674 {
675   DM_Swarm *swarm = (DM_Swarm*)dm->data;
676 
677   PetscFunctionBegin;
678   *dmcell = swarm->dmcell;
679   PetscFunctionReturn(0);
680 }
681 
682 /*@
683    DMSwarmGetLocalSize - Retrives the local length of fields registered
684 
685    Not collective
686 
687    Input parameter:
688 .  dm - a DMSwarm
689 
690    Output parameter:
691 .  nlocal - the length of each registered field
692 
693    Level: beginner
694 
695 .seealso: DMSwarmGetSize(), DMSwarmSetLocalSizes()
696 @*/
697 PetscErrorCode DMSwarmGetLocalSize(DM dm,PetscInt *nlocal)
698 {
699   DM_Swarm *swarm = (DM_Swarm*)dm->data;
700   PetscErrorCode ierr;
701 
702   PetscFunctionBegin;
703   if (nlocal) {ierr = DMSwarmDataBucketGetSizes(swarm->db,nlocal,NULL,NULL);CHKERRQ(ierr);}
704   PetscFunctionReturn(0);
705 }
706 
707 /*@
708    DMSwarmGetSize - Retrives the total length of fields registered
709 
710    Collective on dm
711 
712    Input parameter:
713 .  dm - a DMSwarm
714 
715    Output parameter:
716 .  n - the total length of each registered field
717 
718    Level: beginner
719 
720    Note:
721    This calls MPI_Allreduce upon each call (inefficient but safe)
722 
723 .seealso: DMSwarmGetLocalSize(), DMSwarmSetLocalSizes()
724 @*/
725 PetscErrorCode DMSwarmGetSize(DM dm,PetscInt *n)
726 {
727   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
728   PetscErrorCode ierr;
729   PetscInt       nlocal,ng;
730 
731   PetscFunctionBegin;
732   ierr = DMSwarmDataBucketGetSizes(swarm->db,&nlocal,NULL,NULL);CHKERRQ(ierr);
733   ierr = MPI_Allreduce(&nlocal,&ng,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
734   if (n) { *n = ng; }
735   PetscFunctionReturn(0);
736 }
737 
738 /*@C
739    DMSwarmRegisterPetscDatatypeField - Register a field to a DMSwarm with a native PETSc data type
740 
741    Collective on dm
742 
743    Input parameters:
744 +  dm - a DMSwarm
745 .  fieldname - the textual name to identify this field
746 .  blocksize - the number of each data type
747 -  type - a valid PETSc data type (PETSC_CHAR, PETSC_SHORT, PETSC_INT, PETSC_FLOAT, PETSC_REAL, PETSC_LONG)
748 
749    Level: beginner
750 
751    Notes:
752    The textual name for each registered field must be unique.
753 
754 .seealso: DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField()
755 @*/
756 PetscErrorCode DMSwarmRegisterPetscDatatypeField(DM dm,const char fieldname[],PetscInt blocksize,PetscDataType type)
757 {
758   PetscErrorCode ierr;
759   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
760   size_t         size;
761 
762   PetscFunctionBegin;
763   if (!swarm->field_registration_initialized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmInitializeFieldRegister() first");
764   if (swarm->field_registration_finalized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Cannot register additional fields after calling DMSwarmFinalizeFieldRegister() first");
765 
766   if (type == PETSC_OBJECT) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
767   if (type == PETSC_FUNCTION) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
768   if (type == PETSC_STRING) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
769   if (type == PETSC_STRUCT) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
770   if (type == PETSC_DATATYPE_UNKNOWN) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
771 
772   ierr = PetscDataTypeGetSize(type, &size);CHKERRQ(ierr);
773   /* Load a specific data type into data bucket, specifying textual name and its size in bytes */
774   ierr = DMSwarmDataBucketRegisterField(swarm->db,"DMSwarmRegisterPetscDatatypeField",fieldname,blocksize*size,NULL);CHKERRQ(ierr);
775   {
776     DMSwarmDataField gfield;
777 
778     ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr);
779     ierr = DMSwarmDataFieldSetBlockSize(gfield,blocksize);CHKERRQ(ierr);
780   }
781   swarm->db->field[swarm->db->nfields-1]->petsc_type = type;
782   PetscFunctionReturn(0);
783 }
784 
785 /*@C
786    DMSwarmRegisterUserStructField - Register a user defined struct to a DMSwarm
787 
788    Collective on dm
789 
790    Input parameters:
791 +  dm - a DMSwarm
792 .  fieldname - the textual name to identify this field
793 -  size - the size in bytes of the user struct of each data type
794 
795    Level: beginner
796 
797    Notes:
798    The textual name for each registered field must be unique.
799 
800 .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmRegisterUserDatatypeField()
801 @*/
802 PetscErrorCode DMSwarmRegisterUserStructField(DM dm,const char fieldname[],size_t size)
803 {
804   PetscErrorCode ierr;
805   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
806 
807   PetscFunctionBegin;
808   ierr = DMSwarmDataBucketRegisterField(swarm->db,"DMSwarmRegisterUserStructField",fieldname,size,NULL);CHKERRQ(ierr);
809   swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_STRUCT ;
810   PetscFunctionReturn(0);
811 }
812 
813 /*@C
814    DMSwarmRegisterUserDatatypeField - Register a user defined data type to a DMSwarm
815 
816    Collective on dm
817 
818    Input parameters:
819 +  dm - a DMSwarm
820 .  fieldname - the textual name to identify this field
821 .  size - the size in bytes of the user data type
822 -  blocksize - the number of each data type
823 
824    Level: beginner
825 
826    Notes:
827    The textual name for each registered field must be unique.
828 
829 .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField()
830 @*/
831 PetscErrorCode DMSwarmRegisterUserDatatypeField(DM dm,const char fieldname[],size_t size,PetscInt blocksize)
832 {
833   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
834   PetscErrorCode ierr;
835 
836   PetscFunctionBegin;
837   ierr = DMSwarmDataBucketRegisterField(swarm->db,"DMSwarmRegisterUserDatatypeField",fieldname,blocksize*size,NULL);CHKERRQ(ierr);
838   {
839     DMSwarmDataField gfield;
840 
841     ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr);
842     ierr = DMSwarmDataFieldSetBlockSize(gfield,blocksize);CHKERRQ(ierr);
843   }
844   swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_DATATYPE_UNKNOWN;
845   PetscFunctionReturn(0);
846 }
847 
848 /*@C
849    DMSwarmGetField - Get access to the underlying array storing all entries associated with a registered field
850 
851    Not collective
852 
853    Input parameters:
854 +  dm - a DMSwarm
855 -  fieldname - the textual name to identify this field
856 
857    Output parameters:
858 +  blocksize - the number of each data type
859 .  type - the data type
860 -  data - pointer to raw array
861 
862    Level: beginner
863 
864    Notes:
865    The array must be returned using a matching call to DMSwarmRestoreField().
866 
867 .seealso: DMSwarmRestoreField()
868 @*/
869 PetscErrorCode DMSwarmGetField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data)
870 {
871   DM_Swarm         *swarm = (DM_Swarm*)dm->data;
872   DMSwarmDataField gfield;
873   PetscErrorCode ierr;
874 
875   PetscFunctionBegin;
876   if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); }
877   ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr);
878   ierr = DMSwarmDataFieldGetAccess(gfield);CHKERRQ(ierr);
879   ierr = DMSwarmDataFieldGetEntries(gfield,data);CHKERRQ(ierr);
880   if (blocksize) {*blocksize = gfield->bs; }
881   if (type) { *type = gfield->petsc_type; }
882   PetscFunctionReturn(0);
883 }
884 
885 /*@C
886    DMSwarmRestoreField - Restore access to the underlying array storing all entries associated with a registered field
887 
888    Not collective
889 
890    Input parameters:
891 +  dm - a DMSwarm
892 -  fieldname - the textual name to identify this field
893 
894    Output parameters:
895 +  blocksize - the number of each data type
896 .  type - the data type
897 -  data - pointer to raw array
898 
899    Level: beginner
900 
901    Notes:
902    The user must call DMSwarmGetField() prior to calling DMSwarmRestoreField().
903 
904 .seealso: DMSwarmGetField()
905 @*/
906 PetscErrorCode DMSwarmRestoreField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data)
907 {
908   DM_Swarm         *swarm = (DM_Swarm*)dm->data;
909   DMSwarmDataField gfield;
910   PetscErrorCode ierr;
911 
912   PetscFunctionBegin;
913   ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr);
914   ierr = DMSwarmDataFieldRestoreAccess(gfield);CHKERRQ(ierr);
915   if (data) *data = NULL;
916   PetscFunctionReturn(0);
917 }
918 
919 /*@
920    DMSwarmAddPoint - Add space for one new point in the DMSwarm
921 
922    Not collective
923 
924    Input parameter:
925 .  dm - a DMSwarm
926 
927    Level: beginner
928 
929    Notes:
930    The new point will have all fields initialized to zero.
931 
932 .seealso: DMSwarmAddNPoints()
933 @*/
934 PetscErrorCode DMSwarmAddPoint(DM dm)
935 {
936   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
937   PetscErrorCode ierr;
938 
939   PetscFunctionBegin;
940   if (!swarm->issetup) {ierr = DMSetUp(dm);CHKERRQ(ierr);}
941   ierr = PetscLogEventBegin(DMSWARM_AddPoints,0,0,0,0);CHKERRQ(ierr);
942   ierr = DMSwarmDataBucketAddPoint(swarm->db);CHKERRQ(ierr);
943   ierr = PetscLogEventEnd(DMSWARM_AddPoints,0,0,0,0);CHKERRQ(ierr);
944   PetscFunctionReturn(0);
945 }
946 
947 /*@
948    DMSwarmAddNPoints - Add space for a number of new points in the DMSwarm
949 
950    Not collective
951 
952    Input parameters:
953 +  dm - a DMSwarm
954 -  npoints - the number of new points to add
955 
956    Level: beginner
957 
958    Notes:
959    The new point will have all fields initialized to zero.
960 
961 .seealso: DMSwarmAddPoint()
962 @*/
963 PetscErrorCode DMSwarmAddNPoints(DM dm,PetscInt npoints)
964 {
965   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
966   PetscErrorCode ierr;
967   PetscInt       nlocal;
968 
969   PetscFunctionBegin;
970   ierr = PetscLogEventBegin(DMSWARM_AddPoints,0,0,0,0);CHKERRQ(ierr);
971   ierr = DMSwarmDataBucketGetSizes(swarm->db,&nlocal,NULL,NULL);CHKERRQ(ierr);
972   nlocal = nlocal + npoints;
973   ierr = DMSwarmDataBucketSetSizes(swarm->db,nlocal,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr);
974   ierr = PetscLogEventEnd(DMSWARM_AddPoints,0,0,0,0);CHKERRQ(ierr);
975   PetscFunctionReturn(0);
976 }
977 
978 /*@
979    DMSwarmRemovePoint - Remove the last point from the DMSwarm
980 
981    Not collective
982 
983    Input parameter:
984 .  dm - a DMSwarm
985 
986    Level: beginner
987 
988 .seealso: DMSwarmRemovePointAtIndex()
989 @*/
990 PetscErrorCode DMSwarmRemovePoint(DM dm)
991 {
992   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
993   PetscErrorCode ierr;
994 
995   PetscFunctionBegin;
996   ierr = PetscLogEventBegin(DMSWARM_RemovePoints,0,0,0,0);CHKERRQ(ierr);
997   ierr = DMSwarmDataBucketRemovePoint(swarm->db);CHKERRQ(ierr);
998   ierr = PetscLogEventEnd(DMSWARM_RemovePoints,0,0,0,0);CHKERRQ(ierr);
999   PetscFunctionReturn(0);
1000 }
1001 
1002 /*@
1003    DMSwarmRemovePointAtIndex - Removes a specific point from the DMSwarm
1004 
1005    Not collective
1006 
1007    Input parameters:
1008 +  dm - a DMSwarm
1009 -  idx - index of point to remove
1010 
1011    Level: beginner
1012 
1013 .seealso: DMSwarmRemovePoint()
1014 @*/
1015 PetscErrorCode DMSwarmRemovePointAtIndex(DM dm,PetscInt idx)
1016 {
1017   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1018   PetscErrorCode ierr;
1019 
1020   PetscFunctionBegin;
1021   ierr = PetscLogEventBegin(DMSWARM_RemovePoints,0,0,0,0);CHKERRQ(ierr);
1022   ierr = DMSwarmDataBucketRemovePointAtIndex(swarm->db,idx);CHKERRQ(ierr);
1023   ierr = PetscLogEventEnd(DMSWARM_RemovePoints,0,0,0,0);CHKERRQ(ierr);
1024   PetscFunctionReturn(0);
1025 }
1026 
1027 /*@
1028    DMSwarmCopyPoint - Copy point pj to point pi in the DMSwarm
1029 
1030    Not collective
1031 
1032    Input parameters:
1033 +  dm - a DMSwarm
1034 .  pi - the index of the point to copy
1035 -  pj - the point index where the copy should be located
1036 
1037  Level: beginner
1038 
1039 .seealso: DMSwarmRemovePoint()
1040 @*/
1041 PetscErrorCode DMSwarmCopyPoint(DM dm,PetscInt pi,PetscInt pj)
1042 {
1043   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1044   PetscErrorCode ierr;
1045 
1046   PetscFunctionBegin;
1047   if (!swarm->issetup) {ierr = DMSetUp(dm);CHKERRQ(ierr);}
1048   ierr = DMSwarmDataBucketCopyPoint(swarm->db,pi,swarm->db,pj);CHKERRQ(ierr);
1049   PetscFunctionReturn(0);
1050 }
1051 
1052 PetscErrorCode DMSwarmMigrate_Basic(DM dm,PetscBool remove_sent_points)
1053 {
1054   PetscErrorCode ierr;
1055 
1056   PetscFunctionBegin;
1057   ierr = DMSwarmMigrate_Push_Basic(dm,remove_sent_points);CHKERRQ(ierr);
1058   PetscFunctionReturn(0);
1059 }
1060 
1061 /*@
1062    DMSwarmMigrate - Relocates points defined in the DMSwarm to other MPI-ranks
1063 
1064    Collective on dm
1065 
1066    Input parameters:
1067 +  dm - the DMSwarm
1068 -  remove_sent_points - flag indicating if sent points should be removed from the current MPI-rank
1069 
1070    Notes:
1071    The DM will be modified to accomodate received points.
1072    If remove_sent_points = PETSC_TRUE, any points that were sent will be removed from the DM.
1073    Different styles of migration are supported. See DMSwarmSetMigrateType().
1074 
1075    Level: advanced
1076 
1077 .seealso: DMSwarmSetMigrateType()
1078 @*/
1079 PetscErrorCode DMSwarmMigrate(DM dm,PetscBool remove_sent_points)
1080 {
1081   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1082   PetscErrorCode ierr;
1083 
1084   PetscFunctionBegin;
1085   ierr = PetscLogEventBegin(DMSWARM_Migrate,0,0,0,0);CHKERRQ(ierr);
1086   switch (swarm->migrate_type) {
1087     case DMSWARM_MIGRATE_BASIC:
1088       ierr = DMSwarmMigrate_Basic(dm,remove_sent_points);CHKERRQ(ierr);
1089       break;
1090     case DMSWARM_MIGRATE_DMCELLNSCATTER:
1091       ierr = DMSwarmMigrate_CellDMScatter(dm,remove_sent_points);CHKERRQ(ierr);
1092       break;
1093     case DMSWARM_MIGRATE_DMCELLEXACT:
1094       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_MIGRATE_DMCELLEXACT not implemented");
1095       break;
1096     case DMSWARM_MIGRATE_USER:
1097       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_MIGRATE_USER not implemented");
1098       break;
1099     default:
1100       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_MIGRATE type unknown");
1101       break;
1102   }
1103   ierr = PetscLogEventEnd(DMSWARM_Migrate,0,0,0,0);CHKERRQ(ierr);
1104   PetscFunctionReturn(0);
1105 }
1106 
1107 PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm,PetscInt *globalsize);
1108 
1109 /*
1110  DMSwarmCollectViewCreate
1111 
1112  * Applies a collection method and gathers point neighbour points into dm
1113 
1114  Notes:
1115  Users should call DMSwarmCollectViewDestroy() after
1116  they have finished computations associated with the collected points
1117 */
1118 
1119 /*@
1120    DMSwarmCollectViewCreate - Applies a collection method and gathers points
1121    in neighbour MPI-ranks into the DMSwarm
1122 
1123    Collective on dm
1124 
1125    Input parameter:
1126 .  dm - the DMSwarm
1127 
1128    Notes:
1129    Users should call DMSwarmCollectViewDestroy() after
1130    they have finished computations associated with the collected points
1131    Different collect methods are supported. See DMSwarmSetCollectType().
1132 
1133    Level: advanced
1134 
1135 .seealso: DMSwarmCollectViewDestroy(), DMSwarmSetCollectType()
1136 @*/
1137 PetscErrorCode DMSwarmCollectViewCreate(DM dm)
1138 {
1139   PetscErrorCode ierr;
1140   DM_Swarm *swarm = (DM_Swarm*)dm->data;
1141   PetscInt ng;
1142 
1143   PetscFunctionBegin;
1144   if (swarm->collect_view_active) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"CollectView currently active");
1145   ierr = DMSwarmGetLocalSize(dm,&ng);CHKERRQ(ierr);
1146   switch (swarm->collect_type) {
1147 
1148     case DMSWARM_COLLECT_BASIC:
1149       ierr = DMSwarmMigrate_GlobalToLocal_Basic(dm,&ng);CHKERRQ(ierr);
1150       break;
1151     case DMSWARM_COLLECT_DMDABOUNDINGBOX:
1152       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_COLLECT_DMDABOUNDINGBOX not implemented");
1153       break;
1154     case DMSWARM_COLLECT_GENERAL:
1155       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_COLLECT_GENERAL not implemented");
1156       break;
1157     default:
1158       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_COLLECT type unknown");
1159       break;
1160   }
1161   swarm->collect_view_active = PETSC_TRUE;
1162   swarm->collect_view_reset_nlocal = ng;
1163   PetscFunctionReturn(0);
1164 }
1165 
1166 /*@
1167    DMSwarmCollectViewDestroy - Resets the DMSwarm to the size prior to calling DMSwarmCollectViewCreate()
1168 
1169    Collective on dm
1170 
1171    Input parameters:
1172 .  dm - the DMSwarm
1173 
1174    Notes:
1175    Users should call DMSwarmCollectViewCreate() before this function is called.
1176 
1177    Level: advanced
1178 
1179 .seealso: DMSwarmCollectViewCreate(), DMSwarmSetCollectType()
1180 @*/
1181 PetscErrorCode DMSwarmCollectViewDestroy(DM dm)
1182 {
1183   PetscErrorCode ierr;
1184   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1185 
1186   PetscFunctionBegin;
1187   if (!swarm->collect_view_active) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"CollectView is currently not active");
1188   ierr = DMSwarmSetLocalSizes(dm,swarm->collect_view_reset_nlocal,-1);CHKERRQ(ierr);
1189   swarm->collect_view_active = PETSC_FALSE;
1190   PetscFunctionReturn(0);
1191 }
1192 
1193 PetscErrorCode DMSwarmSetUpPIC(DM dm)
1194 {
1195   PetscInt       dim;
1196   PetscErrorCode ierr;
1197 
1198   PetscFunctionBegin;
1199   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
1200   if (dim < 1) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Dimension must be 1,2,3 - found %D",dim);
1201   if (dim > 3) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Dimension must be 1,2,3 - found %D",dim);
1202   ierr = DMSwarmRegisterPetscDatatypeField(dm,DMSwarmPICField_coor,dim,PETSC_DOUBLE);CHKERRQ(ierr);
1203   ierr = DMSwarmRegisterPetscDatatypeField(dm,DMSwarmPICField_cellid,1,PETSC_INT);CHKERRQ(ierr);
1204   PetscFunctionReturn(0);
1205 }
1206 
1207 /*@
1208    DMSwarmSetType - Set particular flavor of DMSwarm
1209 
1210    Collective on dm
1211 
1212    Input parameters:
1213 +  dm - the DMSwarm
1214 -  stype - the DMSwarm type (e.g. DMSWARM_PIC)
1215 
1216    Level: advanced
1217 
1218 .seealso: DMSwarmSetMigrateType(), DMSwarmSetCollectType()
1219 @*/
1220 PetscErrorCode DMSwarmSetType(DM dm,DMSwarmType stype)
1221 {
1222   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1223   PetscErrorCode ierr;
1224 
1225   PetscFunctionBegin;
1226   swarm->swarm_type = stype;
1227   if (swarm->swarm_type == DMSWARM_PIC) {
1228     ierr = DMSwarmSetUpPIC(dm);CHKERRQ(ierr);
1229   }
1230   PetscFunctionReturn(0);
1231 }
1232 
1233 PetscErrorCode DMSetup_Swarm(DM dm)
1234 {
1235   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1236   PetscErrorCode ierr;
1237   PetscMPIInt    rank;
1238   PetscInt       p,npoints,*rankval;
1239 
1240   PetscFunctionBegin;
1241   if (swarm->issetup) PetscFunctionReturn(0);
1242 
1243   swarm->issetup = PETSC_TRUE;
1244 
1245   if (swarm->swarm_type == DMSWARM_PIC) {
1246     /* check dmcell exists */
1247     if (!swarm->dmcell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSWARM_PIC requires you call DMSwarmSetCellDM");
1248 
1249     if (swarm->dmcell->ops->locatepointssubdomain) {
1250       /* check methods exists for exact ownership identificiation */
1251       ierr = PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->ops->LocatePointsSubdomain\n");CHKERRQ(ierr);
1252       swarm->migrate_type = DMSWARM_MIGRATE_DMCELLEXACT;
1253     } else {
1254       /* check methods exist for point location AND rank neighbor identification */
1255       if (swarm->dmcell->ops->locatepoints) {
1256         ierr = PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->LocatePoints\n");CHKERRQ(ierr);
1257       } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSWARM_PIC requires the method CellDM->ops->locatepoints be defined");
1258 
1259       if (swarm->dmcell->ops->getneighbors) {
1260         ierr = PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->GetNeigbors\n");CHKERRQ(ierr);
1261       } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSWARM_PIC requires the method CellDM->ops->getneighbors be defined");
1262 
1263       swarm->migrate_type = DMSWARM_MIGRATE_DMCELLNSCATTER;
1264     }
1265   }
1266 
1267   ierr = DMSwarmFinalizeFieldRegister(dm);CHKERRQ(ierr);
1268 
1269   /* check some fields were registered */
1270   if (swarm->db->nfields <= 2) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"At least one field user must be registered via DMSwarmRegisterXXX()");
1271 
1272   /* check local sizes were set */
1273   if (swarm->db->L == -1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Local sizes must be set via DMSwarmSetLocalSizes()");
1274 
1275   /* initialize values in pid and rank placeholders */
1276   /* TODO: [pid - use MPI_Scan] */
1277   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
1278   ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
1279   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
1280   for (p=0; p<npoints; p++) {
1281     rankval[p] = (PetscInt)rank;
1282   }
1283   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
1284   PetscFunctionReturn(0);
1285 }
1286 
1287 extern PetscErrorCode DMSwarmSortDestroy(DMSwarmSort *_ctx);
1288 
1289 PetscErrorCode DMDestroy_Swarm(DM dm)
1290 {
1291   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1292   PetscErrorCode ierr;
1293 
1294   PetscFunctionBegin;
1295   ierr = DMSwarmDataBucketDestroy(&swarm->db);CHKERRQ(ierr);
1296   if (swarm->sort_context) {
1297     ierr = DMSwarmSortDestroy(&swarm->sort_context);CHKERRQ(ierr);
1298   }
1299   ierr = PetscFree(swarm);CHKERRQ(ierr);
1300   PetscFunctionReturn(0);
1301 }
1302 
1303 PetscErrorCode DMSwarmView_Draw(DM dm, PetscViewer viewer)
1304 {
1305   DM             cdm;
1306   PetscDraw      draw;
1307   PetscReal     *coords, oldPause;
1308   PetscInt       Np, p, bs;
1309   PetscErrorCode ierr;
1310 
1311   PetscFunctionBegin;
1312   ierr = PetscViewerDrawGetDraw(viewer, 0, &draw);CHKERRQ(ierr);
1313   ierr = DMSwarmGetCellDM(dm, &cdm);CHKERRQ(ierr);
1314   ierr = PetscDrawGetPause(draw, &oldPause);CHKERRQ(ierr);
1315   ierr = PetscDrawSetPause(draw, 0.0);CHKERRQ(ierr);
1316   ierr = DMView(cdm, viewer);CHKERRQ(ierr);
1317   ierr = PetscDrawSetPause(draw, oldPause);CHKERRQ(ierr);
1318 
1319   ierr = DMSwarmGetLocalSize(dm, &Np);CHKERRQ(ierr);
1320   ierr = DMSwarmGetField(dm, DMSwarmPICField_coor, &bs, NULL, (void **) &coords);CHKERRQ(ierr);
1321   for (p = 0; p < Np; ++p) {
1322     const PetscInt i = p*bs;
1323 
1324     ierr = PetscDrawEllipse(draw, coords[i], coords[i+1], 0.01, 0.01, PETSC_DRAW_BLUE);CHKERRQ(ierr);
1325   }
1326   ierr = DMSwarmRestoreField(dm, DMSwarmPICField_coor, &bs, NULL, (void **) &coords);CHKERRQ(ierr);
1327   ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
1328   ierr = PetscDrawPause(draw);CHKERRQ(ierr);
1329   ierr = PetscDrawSave(draw);CHKERRQ(ierr);
1330   PetscFunctionReturn(0);
1331 }
1332 
1333 PetscErrorCode DMView_Swarm(DM dm, PetscViewer viewer)
1334 {
1335   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1336   PetscBool      iascii,ibinary,ishdf5,isvtk,isdraw;
1337   PetscErrorCode ierr;
1338 
1339   PetscFunctionBegin;
1340   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1341   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1342   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
1343   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY,&ibinary);CHKERRQ(ierr);
1344   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK,   &isvtk);CHKERRQ(ierr);
1345   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5,  &ishdf5);CHKERRQ(ierr);
1346   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW,  &isdraw);CHKERRQ(ierr);
1347   if (iascii) {
1348     ierr = DMSwarmDataBucketView(PetscObjectComm((PetscObject)dm),swarm->db,NULL,DATABUCKET_VIEW_STDOUT);CHKERRQ(ierr);
1349   } else if (ibinary) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO Binary support");
1350 #if defined(PETSC_HAVE_HDF5)
1351   else if (ishdf5) {
1352     ierr = DMSwarmView_HDF5(dm, viewer);CHKERRQ(ierr);
1353   }
1354 #else
1355   else if (ishdf5) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"HDF5 not supported. Please reconfigure using --download-hdf5");
1356 #endif
1357   else if (isvtk) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO VTK support");
1358   else if (isdraw) {
1359     ierr = DMSwarmView_Draw(dm, viewer);CHKERRQ(ierr);
1360   }
1361   PetscFunctionReturn(0);
1362 }
1363 
1364 /*MC
1365 
1366  DMSWARM = "swarm" - A DM object used to represent arrays of data (fields) of arbitrary data type.
1367  This implementation was designed for particle methods in which the underlying
1368  data required to be represented is both (i) dynamic in length, (ii) and of arbitrary data type.
1369 
1370  User data can be represented by DMSwarm through a registering "fields".
1371  To register a field, the user must provide;
1372  (a) a unique name;
1373  (b) the data type (or size in bytes);
1374  (c) the block size of the data.
1375 
1376  For example, suppose the application requires a unique id, energy, momentum and density to be stored
1377  on a set of particles. Then the following code could be used
1378 
1379 $    DMSwarmInitializeFieldRegister(dm)
1380 $    DMSwarmRegisterPetscDatatypeField(dm,"uid",1,PETSC_LONG);
1381 $    DMSwarmRegisterPetscDatatypeField(dm,"energy",1,PETSC_REAL);
1382 $    DMSwarmRegisterPetscDatatypeField(dm,"momentum",3,PETSC_REAL);
1383 $    DMSwarmRegisterPetscDatatypeField(dm,"density",1,PETSC_FLOAT);
1384 $    DMSwarmFinalizeFieldRegister(dm)
1385 
1386  The fields represented by DMSwarm are dynamic and can be re-sized at any time.
1387  The only restriction imposed by DMSwarm is that all fields contain the same number of points.
1388 
1389  To support particle methods, "migration" techniques are provided. These methods migrate data
1390  between MPI-ranks.
1391 
1392  DMSwarm supports the methods DMCreateGlobalVector() and DMCreateLocalVector().
1393  As a DMSwarm may internally define and store values of different data types,
1394  before calling DMCreateGlobalVector() or DMCreateLocalVector(), the user must inform DMSwarm which
1395  fields should be used to define a Vec object via
1396    DMSwarmVectorDefineField()
1397  The specified field can be changed at any time - thereby permitting vectors
1398  compatible with different fields to be created.
1399 
1400  A dual representation of fields in the DMSwarm and a Vec object is permitted via
1401    DMSwarmCreateGlobalVectorFromField()
1402  Here the data defining the field in the DMSwarm is shared with a Vec.
1403  This is inherently unsafe if you alter the size of the field at any time between
1404  calls to DMSwarmCreateGlobalVectorFromField() and DMSwarmDestroyGlobalVectorFromField().
1405  If the local size of the DMSwarm does not match the local size of the global vector
1406  when DMSwarmDestroyGlobalVectorFromField() is called, an error is thrown.
1407 
1408  Additional high-level support is provided for Particle-In-Cell methods.
1409  Please refer to the man page for DMSwarmSetType().
1410 
1411  Level: beginner
1412 
1413 .seealso: DMType, DMCreate(), DMSetType()
1414 M*/
1415 PETSC_EXTERN PetscErrorCode DMCreate_Swarm(DM dm)
1416 {
1417   DM_Swarm      *swarm;
1418   PetscErrorCode ierr;
1419 
1420   PetscFunctionBegin;
1421   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1422   ierr     = PetscNewLog(dm,&swarm);CHKERRQ(ierr);
1423   dm->data = swarm;
1424 
1425   ierr = DMSwarmDataBucketCreate(&swarm->db);CHKERRQ(ierr);
1426   ierr = DMSwarmInitializeFieldRegister(dm);CHKERRQ(ierr);
1427 
1428   swarm->vec_field_set = PETSC_FALSE;
1429   swarm->issetup = PETSC_FALSE;
1430   swarm->swarm_type = DMSWARM_BASIC;
1431   swarm->migrate_type = DMSWARM_MIGRATE_BASIC;
1432   swarm->collect_type = DMSWARM_COLLECT_BASIC;
1433   swarm->migrate_error_on_missing_point = PETSC_FALSE;
1434 
1435   swarm->dmcell = NULL;
1436   swarm->collect_view_active = PETSC_FALSE;
1437   swarm->collect_view_reset_nlocal = -1;
1438 
1439   dm->dim  = 0;
1440   dm->ops->view                            = DMView_Swarm;
1441   dm->ops->load                            = NULL;
1442   dm->ops->setfromoptions                  = NULL;
1443   dm->ops->clone                           = NULL;
1444   dm->ops->setup                           = DMSetup_Swarm;
1445   dm->ops->createlocalsection              = NULL;
1446   dm->ops->createdefaultconstraints        = NULL;
1447   dm->ops->createglobalvector              = DMCreateGlobalVector_Swarm;
1448   dm->ops->createlocalvector               = DMCreateLocalVector_Swarm;
1449   dm->ops->getlocaltoglobalmapping         = NULL;
1450   dm->ops->createfieldis                   = NULL;
1451   dm->ops->createcoordinatedm              = NULL;
1452   dm->ops->getcoloring                     = NULL;
1453   dm->ops->creatematrix                    = NULL;
1454   dm->ops->createinterpolation             = NULL;
1455   dm->ops->createinjection                 = NULL;
1456   dm->ops->createmassmatrix                = DMCreateMassMatrix_Swarm;
1457   dm->ops->refine                          = NULL;
1458   dm->ops->coarsen                         = NULL;
1459   dm->ops->refinehierarchy                 = NULL;
1460   dm->ops->coarsenhierarchy                = NULL;
1461   dm->ops->globaltolocalbegin              = NULL;
1462   dm->ops->globaltolocalend                = NULL;
1463   dm->ops->localtoglobalbegin              = NULL;
1464   dm->ops->localtoglobalend                = NULL;
1465   dm->ops->destroy                         = DMDestroy_Swarm;
1466   dm->ops->createsubdm                     = NULL;
1467   dm->ops->getdimpoints                    = NULL;
1468   dm->ops->locatepoints                    = NULL;
1469   PetscFunctionReturn(0);
1470 }
1471