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