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