xref: /petsc/src/dm/impls/swarm/swarm.c (revision 7f978cf1625ca657898bed9392afcbd8113e6639)
1 
2 #define PETSCDM_DLL
3 #include <petsc/private/dmswarmimpl.h>    /*I   "petscdmswarm.h"   I*/
4 #include "data_bucket.h"
5 
6 //typedef PetscErrorCode (*swarm_project)(DM,DM,Vec) DMSwarmProjectMethod; /* swarm, geometry, result */
7 
8 //typedef enum { PROJECT_DMDA_AQ1=0, PROJECT_DMDA_P0 } DMSwarmDMDAProjectionType;
9 
10 #if 0
11 
12 /* Defines what the local space will be */
13 PetscErrorCode DMSwarmSetOverlap(void)
14 {
15 
16   PetscFunctionReturn(0);
17 }
18 
19 
20 /* coordinates */
21 /*
22 DMGetCoordinateDM returns self
23 DMGetCoordinates and DMGetCoordinatesLocal return same thing
24 Local view could be used to define overlapping information
25 */
26 
27 #endif
28 
29 #undef __FUNCT__
30 #define __FUNCT__ "DMSwarmVectorDefineField"
31 PETSC_EXTERN PetscErrorCode DMSwarmVectorDefineField(DM dm,const char fieldname[])
32 {
33   DM_Swarm *swarm = (DM_Swarm*)dm->data;
34   PetscErrorCode ierr;
35   PetscInt bs,n;
36   PetscScalar *array;
37   PetscDataType type;
38 
39   if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); }
40   ierr = DataBucketGetSizes(swarm->db,&n,NULL,NULL);CHKERRQ(ierr);
41   ierr = DMSwarmGetField(dm,fieldname,&bs,&type,(void**)&array);CHKERRQ(ierr);
42 
43   /* Check all fields are of type PETSC_REAL or PETSC_SCALAR */
44   if (type != PETSC_REAL) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid for PETSC_REAL");
45 
46   PetscSNPrintf(swarm->vec_field_name,PETSC_MAX_PATH_LEN-1,"%s",fieldname);
47   swarm->vec_field_set = PETSC_TRUE;
48   swarm->vec_field_bs = 1;//bs;
49   swarm->vec_field_nlocal = n;
50 
51   PetscFunctionReturn(0);
52 }
53 
54 /* requires DMSwarmDefineFieldVector has been called */
55 #undef __FUNCT__
56 #define __FUNCT__ "DMCreateGlobalVector_Swarm"
57 PetscErrorCode DMCreateGlobalVector_Swarm(DM dm,Vec *vec)
58 {
59   DM_Swarm *swarm = (DM_Swarm*)dm->data;
60   PetscErrorCode ierr;
61   Vec x;
62   char name[PETSC_MAX_PATH_LEN];
63 
64   if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); }
65   if (!swarm->vec_field_set) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmVectorDefineField first");
66   PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarmField_%s",swarm->vec_field_name);
67   ierr = VecCreate(PetscObjectComm((PetscObject)dm),&x);CHKERRQ(ierr);
68   ierr = PetscObjectSetName((PetscObject)x,name);CHKERRQ(ierr);
69   ierr = VecSetSizes(x,swarm->vec_field_nlocal,PETSC_DETERMINE);CHKERRQ(ierr);
70   ierr = VecSetBlockSize(x,swarm->vec_field_bs);CHKERRQ(ierr);
71   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
72   *vec = x;
73 
74   PetscFunctionReturn(0);
75 }
76 
77 /* requires DMSwarmDefineFieldVector has been called */
78 #undef __FUNCT__
79 #define __FUNCT__ "DMCreateLocalVector_Swarm"
80 PetscErrorCode DMCreateLocalVector_Swarm(DM dm,Vec *vec)
81 {
82   DM_Swarm *swarm = (DM_Swarm*)dm->data;
83   PetscErrorCode ierr;
84   Vec x;
85   char name[PETSC_MAX_PATH_LEN];
86 
87   if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); }
88   if (!swarm->vec_field_set) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmVectorDefineField first");
89   PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarmField_%s",swarm->vec_field_name);
90   ierr = VecCreate(PETSC_COMM_SELF,&x);CHKERRQ(ierr);
91   ierr = PetscObjectSetName((PetscObject)x,name);CHKERRQ(ierr);
92   ierr = VecSetSizes(x,swarm->vec_field_nlocal,swarm->vec_field_nlocal);CHKERRQ(ierr);
93   ierr = VecSetBlockSize(x,swarm->vec_field_bs);CHKERRQ(ierr);
94   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
95   *vec = x;
96 
97   PetscFunctionReturn(0);
98 }
99 
100 #undef __FUNCT__
101 #define __FUNCT__ "DMSwarmCreateGlobalVectorFromField"
102 PETSC_EXTERN PetscErrorCode DMSwarmCreateGlobalVectorFromField(DM dm,const char fieldname[],Vec *vec)
103 {
104   DM_Swarm *swarm = (DM_Swarm*)dm->data;
105   PetscErrorCode ierr;
106   PetscInt bs,n;
107   PetscScalar *array;
108   Vec x;
109   PetscDataType type;
110   char name[PETSC_MAX_PATH_LEN];
111   PetscMPIInt commsize;
112 
113   if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); }
114   ierr = DataBucketGetSizes(swarm->db,&n,NULL,NULL);CHKERRQ(ierr);
115   ierr = DMSwarmGetField(dm,fieldname,&bs,&type,(void**)&array);CHKERRQ(ierr);
116 
117   /* Check all fields are of type PETSC_REAL or PETSC_SCALAR */
118   if (type != PETSC_REAL) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid for PETSC_REAL");
119 
120   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&commsize);CHKERRQ(ierr);
121   if (commsize == 1) {
122     ierr = VecCreateSeqWithArray(PetscObjectComm((PetscObject)dm),1,n,array,&x);CHKERRQ(ierr);
123   } else {
124     ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)dm),1,n,PETSC_DETERMINE,array,&x);CHKERRQ(ierr);
125   }
126   PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarm_VecFieldInPlace_%s",fieldname);
127   ierr = PetscObjectComposeFunction((PetscObject)x,name,DMSwarmDestroyGlobalVectorFromField);CHKERRQ(ierr);
128 
129   /* Set guard */
130 
131   *vec = x;
132   PetscFunctionReturn(0);
133 }
134 
135 #undef __FUNCT__
136 #define __FUNCT__ "DMSwarmDestroyGlobalVectorFromField"
137 PETSC_EXTERN PetscErrorCode DMSwarmDestroyGlobalVectorFromField(DM dm,const char fieldname[],Vec *vec)
138 {
139   DM_Swarm *swarm = (DM_Swarm*)dm->data;
140   PetscErrorCode ierr;
141   DataField gfield;
142   char name[PETSC_MAX_PATH_LEN];
143   void (*fptr)(void);
144 
145   /* get data field */
146   ierr = DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr);
147 
148   /* check vector is an inplace array */
149   PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarm_VecFieldInPlace_%s",fieldname);
150   ierr = PetscObjectQueryFunction((PetscObject)(*vec),name,&fptr);CHKERRQ(ierr);
151   if (!fptr) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Vector being destroyed was not created from DMSwarm field(%s)",fieldname);
152 
153   /* restore data field */
154   ierr = DataFieldRestoreAccess(gfield);CHKERRQ(ierr);
155 
156   ierr = VecDestroy(vec);CHKERRQ(ierr);
157 
158   PetscFunctionReturn(0);
159 }
160 
161 /*
162 PETSC_EXTERN PetscErrorCode DMSwarmCreateGlobalVectorFromFields(DM dm,const PetscInt nf,const char *fieldnames[],Vec *vec)
163 {
164   PetscFunctionReturn(0);
165 }
166 
167 PETSC_EXTERN PetscErrorCode DMSwarmRestoreGlobalVectorFromFields(DM dm,Vec *vec)
168 {
169   PetscFunctionReturn(0);
170 }
171 */
172 
173 #undef __FUNCT__
174 #define __FUNCT__ "DMSwarmInitializeFieldRegister"
175 PETSC_EXTERN PetscErrorCode DMSwarmInitializeFieldRegister(DM dm)
176 {
177   DM_Swarm *swarm = (DM_Swarm*)dm->data;
178   PetscErrorCode ierr;
179 
180   swarm->field_registration_initialized = PETSC_TRUE;
181 
182   ierr = DMSwarmRegisterPetscDatatypeField(dm,"DMSwarm_pid",1,PETSC_LONG);CHKERRQ(ierr); /* unique identifer */
183   ierr = DMSwarmRegisterPetscDatatypeField(dm,"DMSwarm_rank",1,PETSC_INT);CHKERRQ(ierr); /* used for communication */
184 
185   PetscFunctionReturn(0);
186 }
187 
188 #undef __FUNCT__
189 #define __FUNCT__ "DMSwarmFinalizeFieldRegister"
190 PETSC_EXTERN PetscErrorCode DMSwarmFinalizeFieldRegister(DM dm)
191 {
192   DM_Swarm *swarm = (DM_Swarm*)dm->data;
193   PetscErrorCode ierr;
194 
195   swarm->field_registration_finalized = PETSC_TRUE;
196   ierr = DataBucketFinalize(swarm->db);CHKERRQ(ierr);
197   PetscFunctionReturn(0);
198 }
199 
200 #undef __FUNCT__
201 #define __FUNCT__ "DMSwarmSetLocalSizes"
202 PETSC_EXTERN PetscErrorCode DMSwarmSetLocalSizes(DM dm,PetscInt nlocal,PetscInt buffer)
203 {
204   DM_Swarm *swarm = (DM_Swarm*)dm->data;
205   PetscErrorCode ierr;
206 
207   ierr = DataBucketSetSizes(swarm->db,nlocal,buffer);CHKERRQ(ierr);
208 
209   PetscFunctionReturn(0);
210 }
211 
212 #undef __FUNCT__
213 #define __FUNCT__ "DMSwarmRegisterPetscDatatypeField"
214 PETSC_EXTERN PetscErrorCode DMSwarmRegisterPetscDatatypeField(DM dm,const char fieldname[],PetscInt blocksize,PetscDataType type)
215 {
216   PetscErrorCode ierr;
217   DM_Swarm *swarm = (DM_Swarm*)dm->data;
218   size_t size;
219 
220   if (!swarm->field_registration_initialized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmInitializeFieldRegister() first");
221   if (swarm->field_registration_finalized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Cannot register additional fields after calling DMSwarmFinalizeFieldRegister() first");
222 
223   if (type == PETSC_OBJECT) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
224   if (type == PETSC_FUNCTION) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
225   if (type == PETSC_STRING) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
226   if (type == PETSC_STRUCT) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
227   if (type == PETSC_DATATYPE_UNKNOWN) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
228 
229   switch (type) {
230     case PETSC_CHAR:
231       size = sizeof(PetscChar);
232       break;
233     case PETSC_SHORT:
234       size = sizeof(PetscShort);
235       break;
236     case PETSC_INT:
237       size = sizeof(PetscInt);
238       break;
239     case PETSC_LONG:
240       size = sizeof(Petsc64bitInt);
241       break;
242     case PETSC_FLOAT:
243       size = sizeof(PetscFloat);
244       break;
245     case PETSC_DOUBLE:
246       size = sizeof(PetscReal);
247       break;
248 
249     default:
250       SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
251       break;
252   }
253 
254   /* Load a specific data type into data bucket, specifying textual name and its size in bytes */
255 	ierr = DataBucketRegisterField(swarm->db,"DMSwarmRegisterPetscDatatypeField",fieldname,size,NULL);CHKERRQ(ierr);
256   swarm->db->field[swarm->db->nfields-1]->petsc_type = type;
257 
258   PetscFunctionReturn(0);
259 }
260 
261 #undef __FUNCT__
262 #define __FUNCT__ "DMSwarmRegisterUserStructField"
263 PETSC_EXTERN PetscErrorCode DMSwarmRegisterUserStructField(DM dm,const char fieldname[],size_t size)
264 {
265   PetscErrorCode ierr;
266   DM_Swarm *swarm = (DM_Swarm*)dm->data;
267 
268 	ierr = DataBucketRegisterField(swarm->db,"DMSwarmRegisterUserStructField",fieldname,size,NULL);CHKERRQ(ierr);
269   swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_STRUCT ;
270 
271   PetscFunctionReturn(0);
272 }
273 
274 #undef __FUNCT__
275 #define __FUNCT__ "DMSwarmRegisterUserDatatypeField"
276 PETSC_EXTERN PetscErrorCode DMSwarmRegisterUserDatatypeField(DM dm,const char fieldname[],size_t size)
277 {
278   DM_Swarm *swarm = (DM_Swarm*)dm->data;
279   PetscErrorCode ierr;
280 
281 	ierr = DataBucketRegisterField(swarm->db,"DMSwarmRegisterUserDatatypeField",fieldname,size,NULL);CHKERRQ(ierr);
282   swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_DATATYPE_UNKNOWN;
283 
284   PetscFunctionReturn(0);
285 }
286 
287 #undef __FUNCT__
288 #define __FUNCT__ "DMSwarmGetField"
289 PETSC_EXTERN PetscErrorCode DMSwarmGetField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data)
290 {
291   DM_Swarm *swarm = (DM_Swarm*)dm->data;
292   DataField gfield;
293   PetscErrorCode ierr;
294 
295   if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); }
296 
297   ierr = DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr);
298   ierr = DataFieldGetAccess(gfield);CHKERRQ(ierr);
299   ierr = DataFieldGetEntries(gfield,data);CHKERRQ(ierr);
300   if (blocksize) {}
301   if (type) { *type = gfield->petsc_type; }
302 
303   PetscFunctionReturn(0);
304 }
305 
306 #undef __FUNCT__
307 #define __FUNCT__ "DMSwarmRestoreField"
308 PETSC_EXTERN PetscErrorCode DMSwarmRestoreField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data)
309 {
310   DM_Swarm *swarm = (DM_Swarm*)dm->data;
311   DataField gfield;
312   PetscErrorCode ierr;
313 
314   ierr = DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr);
315   ierr = DataFieldRestoreAccess(gfield);CHKERRQ(ierr);
316   if (data) *data = NULL;
317 
318   PetscFunctionReturn(0);
319 }
320 
321 #undef __FUNCT__
322 #define __FUNCT__ "DMSwarmAddPoint"
323 PETSC_EXTERN PetscErrorCode DMSwarmAddPoint(DM dm)
324 {
325   DM_Swarm *swarm = (DM_Swarm*)dm->data;
326   PetscErrorCode ierr;
327 
328   if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); }
329   ierr = DataBucketAddPoint(swarm->db);CHKERRQ(ierr);
330   PetscFunctionReturn(0);
331 }
332 
333 #undef __FUNCT__
334 #define __FUNCT__ "DMSwarmAddNPoints"
335 PETSC_EXTERN PetscErrorCode DMSwarmAddNPoints(DM dm,PetscInt npoints)
336 {
337   DM_Swarm *swarm = (DM_Swarm*)dm->data;
338   PetscErrorCode ierr;
339   PetscInt nlocal;
340 
341   ierr = DataBucketGetSizes(swarm->db,&nlocal,NULL,NULL);CHKERRQ(ierr);
342   nlocal = nlocal + npoints;
343   ierr = DataBucketSetSizes(swarm->db,nlocal,-1);CHKERRQ(ierr);
344   PetscFunctionReturn(0);
345 }
346 
347 #undef __FUNCT__
348 #define __FUNCT__ "DMSwarmRemovePoint"
349 PETSC_EXTERN PetscErrorCode DMSwarmRemovePoint(DM dm)
350 {
351   DM_Swarm *swarm = (DM_Swarm*)dm->data;
352   PetscErrorCode ierr;
353 
354   ierr = DataBucketRemovePoint(swarm->db);CHKERRQ(ierr);
355   PetscFunctionReturn(0);
356 }
357 
358 #undef __FUNCT__
359 #define __FUNCT__ "DMSwarmRemovePointAtIndex"
360 PETSC_EXTERN PetscErrorCode DMSwarmRemovePointAtIndex(DM dm,PetscInt idx)
361 {
362   DM_Swarm *swarm = (DM_Swarm*)dm->data;
363   PetscErrorCode ierr;
364 
365   ierr = DataBucketRemovePointAtIndex(swarm->db,idx);CHKERRQ(ierr);
366   PetscFunctionReturn(0);
367 }
368 
369 #undef __FUNCT__
370 #define __FUNCT__ "DMSwarmPush_Basic"
371 PetscErrorCode DMSwarmPush_Basic(DM dm,PetscBool remove_sent_points)
372 {
373   PetscFunctionReturn(0);
374 }
375 
376 #undef __FUNCT__
377 #define __FUNCT__ "DMSwarmPush"
378 PETSC_EXTERN PetscErrorCode DMSwarmPush(DM dm,PetscBool remove_sent_points)
379 {
380   PetscErrorCode ierr;
381   ierr = DMSwarmPush_Basic(dm,remove_sent_points);CHKERRQ(ierr);
382   PetscFunctionReturn(0);
383 }
384 
385 
386 #undef __FUNCT__
387 #define __FUNCT__ "DMSetup_Swarm"
388 PetscErrorCode DMSetup_Swarm(DM dm)
389 {
390   DM_Swarm *swarm = (DM_Swarm*)dm->data;
391   PetscErrorCode ierr;
392   PetscMPIInt rank;
393   PetscInt p,npoints,*rankval;
394 
395   if (swarm->issetup) PetscFunctionReturn(0);
396 
397   PetscPrintf(PETSC_COMM_SELF,"Swarm setup \n");
398   swarm->issetup = PETSC_TRUE;
399 
400   /* check some fields were registered */
401   if (swarm->db->nfields <= 2) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"At least one field user must be registered via DMSwarmRegisterXXX()");
402 
403   /* check local sizes were set */
404   if (swarm->db->L == -1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Local sizes must be set via DMSwarmSetLocalSizes()");
405 
406   /* initialize values in pid and rank placeholders */
407   /* TODO: [pid - use MPI_Scan] */
408 
409   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
410   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
411   ierr = DMSwarmGetField(dm,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
412   for (p=0; p<npoints; p++) {
413     rankval[p] = (PetscInt)rank;
414   }
415   ierr = DMSwarmRestoreField(dm,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
416 
417   PetscFunctionReturn(0);
418 }
419 
420 #undef __FUNCT__
421 #define __FUNCT__ "DMDestroy_Swarm"
422 PetscErrorCode DMDestroy_Swarm(DM dm)
423 {
424   DM_Swarm *swarm = (DM_Swarm*)dm->data;
425   PetscErrorCode ierr;
426 
427   PetscFunctionBegin;
428   ierr = DataBucketDestroy(&swarm->db);CHKERRQ(ierr);
429   ierr = PetscFree(swarm);CHKERRQ(ierr);
430   PetscFunctionReturn(0);
431 }
432 
433 #undef __FUNCT__
434 #define __FUNCT__ "DMView_Swarm"
435 PetscErrorCode DMView_Swarm(DM dm, PetscViewer viewer)
436 {
437   DM_Swarm *swarm = (DM_Swarm*)dm->data;
438   PetscBool      iascii,ibinary,ishdf5,isvtk;
439   PetscErrorCode ierr;
440 
441   PetscFunctionBegin;
442   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
443   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
444   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
445   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY,&ibinary);CHKERRQ(ierr);
446   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK,   &isvtk);CHKERRQ(ierr);
447   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5,  &ishdf5);CHKERRQ(ierr);
448   if (iascii) {
449     ierr = DataBucketView(PetscObjectComm((PetscObject)dm),swarm->db,NULL,DATABUCKET_VIEW_STDOUT);CHKERRQ(ierr);
450   } else if (ibinary) {
451     SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO VTK support");
452   } else if (ishdf5) {
453 #if defined(PETSC_HAVE_HDF5)
454     SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO HDF5 support");
455 #else
456     SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"HDF5 not supported. Please reconfigure using --download-hdf5");
457 #endif
458   } else if (isvtk) {
459     SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO VTK support");
460   }
461   PetscFunctionReturn(0);
462 }
463 
464 #undef __FUNCT__
465 #define __FUNCT__ "DMCreate_Swarm"
466 PETSC_EXTERN PetscErrorCode DMCreate_Swarm(DM dm)
467 {
468   DM_Swarm      *swarm;
469   PetscErrorCode ierr;
470 
471   PetscFunctionBegin;
472   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
473   ierr     = PetscNewLog(dm,&swarm);CHKERRQ(ierr);
474 
475   ierr = DataBucketCreate(&swarm->db);CHKERRQ(ierr);
476   swarm->vec_field_set = PETSC_FALSE;
477   swarm->issetup = PETSC_FALSE;
478 
479   dm->dim  = 0;
480   dm->data = swarm;
481 
482   dm->ops->view                            = DMView_Swarm;
483   dm->ops->load                            = NULL;
484   dm->ops->setfromoptions                  = NULL;
485   dm->ops->clone                           = NULL;
486   dm->ops->setup                           = DMSetup_Swarm;
487   dm->ops->createdefaultsection            = NULL;
488   dm->ops->createdefaultconstraints        = NULL;
489   dm->ops->createglobalvector              = DMCreateGlobalVector_Swarm;
490   dm->ops->createlocalvector               = DMCreateLocalVector_Swarm;
491   dm->ops->getlocaltoglobalmapping         = NULL;
492   dm->ops->createfieldis                   = NULL;
493   dm->ops->createcoordinatedm              = NULL;
494   dm->ops->getcoloring                     = NULL;
495   dm->ops->creatematrix                    = NULL;
496   dm->ops->createinterpolation             = NULL;
497   dm->ops->getaggregates                   = NULL;
498   dm->ops->getinjection                    = NULL;
499   dm->ops->refine                          = NULL;
500   dm->ops->coarsen                         = NULL;
501   dm->ops->refinehierarchy                 = NULL;
502   dm->ops->coarsenhierarchy                = NULL;
503   dm->ops->globaltolocalbegin              = NULL;
504   dm->ops->globaltolocalend                = NULL;
505   dm->ops->localtoglobalbegin              = NULL;
506   dm->ops->localtoglobalend                = NULL;
507   dm->ops->destroy                         = DMDestroy_Swarm;
508   dm->ops->createsubdm                     = NULL;
509   dm->ops->getdimpoints                    = NULL;
510   dm->ops->locatepoints                    = NULL;
511 
512   PetscFunctionReturn(0);
513 }