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