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