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