Actual source code: swarm.c

petsc-master 2020-02-22
Report Typos and Errors
  1: #define PETSCDM_DLL
  2:  #include <petsc/private/dmswarmimpl.h>
  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"

 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;

 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 };

 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";

 24: /*@C
 25:    DMSwarmVectorDefineField - Sets the field from which to define a Vec object
 26:                              when DMCreateLocalVector(), or DMCreateGlobalVector() is called

 28:    Collective on dm

 30:    Input parameters:
 31: +  dm - a DMSwarm
 32: -  fieldname - the textual name given to a registered field

 34:    Level: beginner

 36:    Notes:

 38:    The field with name fieldname must be defined as having a data type of PetscScalar.

 40:    This function must be called prior to calling DMCreateLocalVector(), DMCreateGlobalVector().
 41:    Mutiple calls to DMSwarmVectorDefineField() are permitted.

 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;
 49:   PetscInt       bs,n;
 50:   PetscScalar    *array;
 51:   PetscDataType  type;

 54:   if (!swarm->issetup) { DMSetUp(dm); }
 55:   DMSwarmDataBucketGetSizes(swarm->db,&n,NULL,NULL);
 56:   DMSwarmGetField(dm,fieldname,&bs,&type,(void**)&array);

 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:   PetscSNPrintf(swarm->vec_field_name,PETSC_MAX_PATH_LEN-1,"%s",fieldname);
 61:   swarm->vec_field_set = PETSC_TRUE;
 62:   swarm->vec_field_bs = bs;
 63:   swarm->vec_field_nlocal = n;
 64:   DMSwarmRestoreField(dm,fieldname,&bs,&type,(void**)&array);
 65:   return(0);
 66: }

 68: /* requires DMSwarmDefineFieldVector has been called */
 69: PetscErrorCode DMCreateGlobalVector_Swarm(DM dm,Vec *vec)
 70: {
 71:   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
 73:   Vec            x;
 74:   char           name[PETSC_MAX_PATH_LEN];

 77:   if (!swarm->issetup) { DMSetUp(dm); }
 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 */

 81:   PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarmField_%s",swarm->vec_field_name);
 82:   VecCreate(PetscObjectComm((PetscObject)dm),&x);
 83:   PetscObjectSetName((PetscObject)x,name);
 84:   VecSetSizes(x,swarm->db->L*swarm->vec_field_bs,PETSC_DETERMINE);
 85:   VecSetBlockSize(x,swarm->vec_field_bs);
 86:   VecSetFromOptions(x);
 87:   *vec = x;
 88:   return(0);
 89: }

 91: /* requires DMSwarmDefineFieldVector has been called */
 92: PetscErrorCode DMCreateLocalVector_Swarm(DM dm,Vec *vec)
 93: {
 94:   DM_Swarm *swarm = (DM_Swarm*)dm->data;
 96:   Vec x;
 97:   char name[PETSC_MAX_PATH_LEN];

100:   if (!swarm->issetup) { DMSetUp(dm); }
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 */

104:   PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarmField_%s",swarm->vec_field_name);
105:   VecCreate(PETSC_COMM_SELF,&x);
106:   PetscObjectSetName((PetscObject)x,name);
107:   VecSetSizes(x,swarm->db->L*swarm->vec_field_bs,PETSC_DETERMINE);
108:   VecSetBlockSize(x,swarm->vec_field_bs);
109:   VecSetFromOptions(x);
110:   *vec = x;
111:   return(0);
112: }

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];

124:   VecGetLocalSize(*vec, &nlocal);
125:   VecGetBlockSize(*vec, &bs);
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:   DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield);
128:   /* check vector is an inplace array */
129:   PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "DMSwarm_VecFieldInPlace_%s", fieldname);
130:   PetscObjectQueryFunction((PetscObject) *vec, name, &fptr);
131:   if (!fptr) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_USER, "Vector being destroyed was not created from DMSwarm field(%s)", fieldname);
132:   DMSwarmDataFieldRestoreAccess(gfield);
133:   VecDestroy(vec);
134:   return(0);
135: }

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;

148:   if (!swarm->issetup) {DMSetUp(dm);}
149:   DMSwarmDataBucketGetSizes(swarm->db, &n, NULL, NULL);
150:   DMSwarmGetField(dm, fieldname, &bs, &type, (void **) &array);
151:   if (type != PETSC_REAL) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Only valid for PETSC_REAL");

153:   MPI_Comm_size(comm, &size);
154:   if (size == 1) {
155:     VecCreateSeqWithArray(comm, bs, n*bs, array, vec);
156:   } else {
157:     VecCreateMPIWithArray(comm, bs, n*bs, PETSC_DETERMINE, array, vec);
158:   }
159:   PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "DMSwarmSharedField_%s", fieldname);
160:   PetscObjectSetName((PetscObject) *vec, name);

162:   /* Set guard */
163:   PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "DMSwarm_VecFieldInPlace_%s", fieldname);
164:   PetscObjectComposeFunction((PetscObject) *vec, name, DMSwarmDestroyVectorFromField_Private);
165:   return(0);
166: }

168: /* This creates a "mass matrix" between a finite element and particle space. If a finite element interpolant is given by

170:      \hat f = \sum_i f_i \phi_i

172:    and a particle function is given by

174:      f = \sum_i w_i \delta(x - x_i)

176:    then we want to require that

178:      M \hat f = M_p f

180:    where the particle mass matrix is given by

182:      (M_p)_{ij} = \int \phi_i \delta(x - x_j)

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;

203:   PetscObjectGetComm((PetscObject) mass, &comm);
204:   DMGetCoordinateDim(dmf, &dim);
205:   DMGetDS(dmf, &prob);
206:   PetscDSGetNumFields(prob, &Nf);
207:   PetscDSGetTotalDimension(prob, &totDim);
208:   PetscMalloc3(dim, &v0, dim*dim, &J, dim*dim,&invJ);
209:   DMGetLocalSection(dmf, &fsection);
210:   DMGetGlobalSection(dmf, &globalFSection);
211:   DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd);
212:   MatGetLocalSize(mass, &locRows, &locCols);

214:   PetscLayoutCreate(comm, &colLayout);
215:   PetscLayoutSetLocalSize(colLayout, locCols);
216:   PetscLayoutSetBlockSize(colLayout, 1);
217:   PetscLayoutSetUp(colLayout);
218:   PetscLayoutGetRange(colLayout, &colStart, &colEnd);
219:   PetscLayoutDestroy(&colLayout);

221:   PetscLayoutCreate(comm, &rLayout);
222:   PetscLayoutSetLocalSize(rLayout, locRows);
223:   PetscLayoutSetBlockSize(rLayout, 1);
224:   PetscLayoutSetUp(rLayout);
225:   PetscLayoutGetRange(rLayout, &rStart, NULL);
226:   PetscLayoutDestroy(&rLayout);

228:   PetscCalloc2(locRows, &dnz, locRows, &onz);
229:   PetscHSetIJCreate(&ht);

231:   PetscSynchronizedFlush(comm, NULL);
232:   /* count non-zeros */
233:   DMSwarmSortGetAccess(dmc);
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;

240:       DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);
241:       DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices);
242:       maxC = PetscMax(maxC, numCIndices);
243:       {
244:         PetscHashIJKey key;
245:         PetscBool      missing;
246:         for (i = 0; i < numFIndices; ++i) {
247:           key.j = findices[i]; /* global column (from Plex) */
248:           if (key.j >= 0) {
249:             /* Get indices for coarse elements */
250:             for (c = 0; c < numCIndices; ++c) {
251:               key.i = cindices[c] + rStart; /* global cols (from Swarm) */
252:               if (key.i < 0) continue;
253:               PetscHSetIJQueryAdd(ht, key, &missing);
254:               if (missing) {
255:                 if ((key.j >= colStart) && (key.j < colEnd)) ++dnz[key.i - rStart];
256:                 else                                         ++onz[key.i - rStart];
257:               } else SETERRQ2(PetscObjectComm((PetscObject) dmf), PETSC_ERR_SUP, "Set new value at %D,%D", key.i, key.j);
258:             }
259:           }
260:         }
261:         PetscFree(cindices);
262:       }
263:       DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);
264:     }
265:   }
266:   PetscHSetIJDestroy(&ht);
267:   MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL);
268:   MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE);
269:   PetscFree2(dnz, onz);
270:   PetscMalloc3(maxC*totDim, &elemMat, maxC, &rowIDXs, maxC*dim, &xi);
271:   for (field = 0; field < Nf; ++field) {
272:     PetscTabulation Tcoarse;
273:     PetscObject       obj;
274:     PetscReal        *coords;
275:     PetscInt          Nc, i;

277:     PetscDSGetDiscretization(prob, field, &obj);
278:     PetscFEGetNumComponents((PetscFE) obj, &Nc);
279:     if (Nc != 1) SETERRQ1(PetscObjectComm((PetscObject) dmf), PETSC_ERR_SUP, "Can only interpolate a scalar field from particles, Nc = %D", Nc);
280:     DMSwarmGetField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);
281:     for (cell = cStart; cell < cEnd; ++cell) {
282:       PetscInt *findices  , *cindices;
283:       PetscInt  numFIndices, numCIndices;
284:       PetscInt  p, c;

286:       /* TODO: Use DMField instead of assuming affine */
287:       DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);
288:       DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);
289:       DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices);
290:       for (p = 0; p < numCIndices; ++p) {
291:         CoordinatesRealToRef(dim, dim, v0ref, v0, invJ, &coords[cindices[p]*dim], &xi[p*dim]);
292:       }
293:       PetscFECreateTabulation((PetscFE) obj, 1, numCIndices, xi, 0, &Tcoarse);
294:       /* Get elemMat entries by multiplying by weight */
295:       PetscArrayzero(elemMat, numCIndices*totDim);
296:       for (i = 0; i < numFIndices; ++i) {
297:         for (p = 0; p < numCIndices; ++p) {
298:           for (c = 0; c < Nc; ++c) {
299:             /* B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c */
300:             elemMat[p*numFIndices+i] += Tcoarse->T[0][(p*numFIndices + i)*Nc + c]*(useDeltaFunction ? 1.0 : detJ);
301:           }
302:         }
303:       }
304:       for (p = 0; p < numCIndices; ++p) rowIDXs[p] = cindices[p] + rStart;
305:       if (0) {DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat);}
306:       MatSetValues(mass, numCIndices, rowIDXs, numFIndices, findices, elemMat, ADD_VALUES);
307:       PetscFree(cindices);
308:       DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);
309:       PetscTabulationDestroy(&Tcoarse);
310:     }
311:     DMSwarmRestoreField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);
312:   }
313:   PetscFree3(elemMat, rowIDXs, xi);
314:   DMSwarmSortRestoreAccess(dmc);
315:   PetscFree3(v0, J, invJ);
316:   MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY);
317:   MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY);
318:   return(0);
319: }

321: /* FEM cols, Particle rows */
322: static PetscErrorCode DMCreateMassMatrix_Swarm(DM dmCoarse, DM dmFine, Mat *mass)
323: {
324:   PetscSection   gsf;
325:   PetscInt       m, n;
326:   void          *ctx;

330:   DMGetGlobalSection(dmFine, &gsf);
331:   PetscSectionGetConstrainedStorageSize(gsf, &m);
332:   DMSwarmGetLocalSize(dmCoarse, &n);
333:   MatCreate(PetscObjectComm((PetscObject) dmCoarse), mass);
334:   MatSetSizes(*mass, n, m, PETSC_DETERMINE, PETSC_DETERMINE);
335:   MatSetType(*mass, dmCoarse->mattype);
336:   DMGetApplicationContext(dmFine, &ctx);

338:   DMSwarmComputeMassMatrix_Private(dmCoarse, dmFine, *mass, PETSC_TRUE, ctx);
339:   MatViewFromOptions(*mass, NULL, "-mass_mat_view");
340:   return(0);
341: }

343: /*@C
344:    DMSwarmCreateGlobalVectorFromField - Creates a Vec object sharing the array associated with a given field

346:    Collective on dm

348:    Input parameters:
349: +  dm - a DMSwarm
350: -  fieldname - the textual name given to a registered field

352:    Output parameter:
353: .  vec - the vector

355:    Level: beginner

357:    Notes:
358:    The vector must be returned using a matching call to DMSwarmDestroyGlobalVectorFromField().

360: .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmDestroyGlobalVectorFromField()
361: @*/
362: PETSC_EXTERN PetscErrorCode DMSwarmCreateGlobalVectorFromField(DM dm,const char fieldname[],Vec *vec)
363: {
364:   MPI_Comm       comm = PetscObjectComm((PetscObject) dm);

368:   DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec);
369:   return(0);
370: }

372: /*@C
373:    DMSwarmDestroyGlobalVectorFromField - Destroys the Vec object which share the array associated with a given field

375:    Collective on dm

377:    Input parameters:
378: +  dm - a DMSwarm
379: -  fieldname - the textual name given to a registered field

381:    Output parameter:
382: .  vec - the vector

384:    Level: beginner

386: .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmCreateGlobalVectorFromField()
387: @*/
388: PETSC_EXTERN PetscErrorCode DMSwarmDestroyGlobalVectorFromField(DM dm,const char fieldname[],Vec *vec)
389: {

393:   DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec);
394:   return(0);
395: }

397: /*@C
398:    DMSwarmCreateLocalVectorFromField - Creates a Vec object sharing the array associated with a given field

400:    Collective on dm

402:    Input parameters:
403: +  dm - a DMSwarm
404: -  fieldname - the textual name given to a registered field

406:    Output parameter:
407: .  vec - the vector

409:    Level: beginner

411:    Notes:
412:    The vector must be returned using a matching call to DMSwarmDestroyLocalVectorFromField().

414: .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmDestroyLocalVectorFromField()
415: @*/
416: PETSC_EXTERN PetscErrorCode DMSwarmCreateLocalVectorFromField(DM dm,const char fieldname[],Vec *vec)
417: {
418:   MPI_Comm       comm = PETSC_COMM_SELF;

422:   DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec);
423:   return(0);
424: }

426: /*@C
427:    DMSwarmDestroyLocalVectorFromField - Destroys the Vec object which share the array associated with a given field

429:    Collective on dm

431:    Input parameters:
432: +  dm - a DMSwarm
433: -  fieldname - the textual name given to a registered field

435:    Output parameter:
436: .  vec - the vector

438:    Level: beginner

440: .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmCreateLocalVectorFromField()
441: @*/
442: PETSC_EXTERN PetscErrorCode DMSwarmDestroyLocalVectorFromField(DM dm,const char fieldname[],Vec *vec)
443: {

447:   DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec);
448:   return(0);
449: }

451: /*
452: PETSC_EXTERN PetscErrorCode DMSwarmCreateGlobalVectorFromFields(DM dm,const PetscInt nf,const char *fieldnames[],Vec *vec)
453: {
454:   return(0);
455: }

457: PETSC_EXTERN PetscErrorCode DMSwarmRestoreGlobalVectorFromFields(DM dm,Vec *vec)
458: {
459:   return(0);
460: }
461: */

463: /*@C
464:    DMSwarmInitializeFieldRegister - Initiates the registration of fields to a DMSwarm

466:    Collective on dm

468:    Input parameter:
469: .  dm - a DMSwarm

471:    Level: beginner

473:    Notes:
474:    After all fields have been registered, you must call DMSwarmFinalizeFieldRegister().

476: .seealso: DMSwarmFinalizeFieldRegister(), DMSwarmRegisterPetscDatatypeField(),
477:  DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField()
478: @*/
479: PETSC_EXTERN PetscErrorCode DMSwarmInitializeFieldRegister(DM dm)
480: {
481:   DM_Swarm      *swarm = (DM_Swarm *) dm->data;

485:   if (!swarm->field_registration_initialized) {
486:     swarm->field_registration_initialized = PETSC_TRUE;
487:     DMSwarmRegisterPetscDatatypeField(dm,DMSwarmField_pid,1,PETSC_INT64); /* unique identifer */
488:     DMSwarmRegisterPetscDatatypeField(dm,DMSwarmField_rank,1,PETSC_INT); /* used for communication */
489:   }
490:   return(0);
491: }

493: /*@C
494:    DMSwarmFinalizeFieldRegister - Finalizes the registration of fields to a DMSwarm

496:    Collective on dm

498:    Input parameter:
499: .  dm - a DMSwarm

501:    Level: beginner

503:    Notes:
504:    After DMSwarmFinalizeFieldRegister() has been called, no new fields can be defined on the DMSwarm.

506: .seealso: DMSwarmInitializeFieldRegister(), DMSwarmRegisterPetscDatatypeField(),
507:  DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField()
508: @*/
509: PETSC_EXTERN PetscErrorCode DMSwarmFinalizeFieldRegister(DM dm)
510: {
511:   DM_Swarm       *swarm = (DM_Swarm*)dm->data;

515:   if (!swarm->field_registration_finalized) {
516:     DMSwarmDataBucketFinalize(swarm->db);
517:   }
518:   swarm->field_registration_finalized = PETSC_TRUE;
519:   return(0);
520: }

522: /*@C
523:    DMSwarmSetLocalSizes - Sets the length of all registered fields on the DMSwarm

525:    Not collective

527:    Input parameters:
528: +  dm - a DMSwarm
529: .  nlocal - the length of each registered field
530: -  buffer - the length of the buffer used to efficient dynamic re-sizing

532:    Level: beginner

534: .seealso: DMSwarmGetLocalSize()
535: @*/
536: PETSC_EXTERN PetscErrorCode DMSwarmSetLocalSizes(DM dm,PetscInt nlocal,PetscInt buffer)
537: {
538:   DM_Swarm       *swarm = (DM_Swarm*)dm->data;

542:   PetscLogEventBegin(DMSWARM_SetSizes,0,0,0,0);
543:   DMSwarmDataBucketSetSizes(swarm->db,nlocal,buffer);
544:   PetscLogEventEnd(DMSWARM_SetSizes,0,0,0,0);
545:   return(0);
546: }

548: /*@C
549:    DMSwarmSetCellDM - Attachs a DM to a DMSwarm

551:    Collective on dm

553:    Input parameters:
554: +  dm - a DMSwarm
555: -  dmcell - the DM to attach to the DMSwarm

557:    Level: beginner

559:    Notes:
560:    The attached DM (dmcell) will be queried for point location and
561:    neighbor MPI-rank information if DMSwarmMigrate() is called.

563: .seealso: DMSwarmSetType(), DMSwarmGetCellDM(), DMSwarmMigrate()
564: @*/
565: PETSC_EXTERN PetscErrorCode DMSwarmSetCellDM(DM dm,DM dmcell)
566: {
567:   DM_Swarm *swarm = (DM_Swarm*)dm->data;

570:   swarm->dmcell = dmcell;
571:   return(0);
572: }

574: /*@C
575:    DMSwarmGetCellDM - Fetches the attached cell DM

577:    Collective on dm

579:    Input parameter:
580: .  dm - a DMSwarm

582:    Output parameter:
583: .  dmcell - the DM which was attached to the DMSwarm

585:    Level: beginner

587: .seealso: DMSwarmSetCellDM()
588: @*/
589: PETSC_EXTERN PetscErrorCode DMSwarmGetCellDM(DM dm,DM *dmcell)
590: {
591:   DM_Swarm *swarm = (DM_Swarm*)dm->data;

594:   *dmcell = swarm->dmcell;
595:   return(0);
596: }

598: /*@C
599:    DMSwarmGetLocalSize - Retrives the local length of fields registered

601:    Not collective

603:    Input parameter:
604: .  dm - a DMSwarm

606:    Output parameter:
607: .  nlocal - the length of each registered field

609:    Level: beginner

611: .seealso: DMSwarmGetSize(), DMSwarmSetLocalSizes()
612: @*/
613: PETSC_EXTERN PetscErrorCode DMSwarmGetLocalSize(DM dm,PetscInt *nlocal)
614: {
615:   DM_Swarm *swarm = (DM_Swarm*)dm->data;

619:   if (nlocal) {DMSwarmDataBucketGetSizes(swarm->db,nlocal,NULL,NULL);}
620:   return(0);
621: }

623: /*@C
624:    DMSwarmGetSize - Retrives the total length of fields registered

626:    Collective on dm

628:    Input parameter:
629: .  dm - a DMSwarm

631:    Output parameter:
632: .  n - the total length of each registered field

634:    Level: beginner

636:    Note:
637:    This calls MPI_Allreduce upon each call (inefficient but safe)

639: .seealso: DMSwarmGetLocalSize(), DMSwarmSetLocalSizes()
640: @*/
641: PETSC_EXTERN PetscErrorCode DMSwarmGetSize(DM dm,PetscInt *n)
642: {
643:   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
645:   PetscInt       nlocal,ng;

648:   DMSwarmDataBucketGetSizes(swarm->db,&nlocal,NULL,NULL);
649:   MPI_Allreduce(&nlocal,&ng,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));
650:   if (n) { *n = ng; }
651:   return(0);
652: }

654: /*@C
655:    DMSwarmRegisterPetscDatatypeField - Register a field to a DMSwarm with a native PETSc data type

657:    Collective on dm

659:    Input parameters:
660: +  dm - a DMSwarm
661: .  fieldname - the textual name to identify this field
662: .  blocksize - the number of each data type
663: -  type - a valid PETSc data type (PETSC_CHAR, PETSC_SHORT, PETSC_INT, PETSC_FLOAT, PETSC_REAL, PETSC_LONG)

665:    Level: beginner

667:    Notes:
668:    The textual name for each registered field must be unique.

670: .seealso: DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField()
671: @*/
672: PETSC_EXTERN PetscErrorCode DMSwarmRegisterPetscDatatypeField(DM dm,const char fieldname[],PetscInt blocksize,PetscDataType type)
673: {
675:   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
676:   size_t         size;

679:   if (!swarm->field_registration_initialized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmInitializeFieldRegister() first");
680:   if (swarm->field_registration_finalized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Cannot register additional fields after calling DMSwarmFinalizeFieldRegister() first");

682:   if (type == PETSC_OBJECT) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
683:   if (type == PETSC_FUNCTION) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
684:   if (type == PETSC_STRING) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
685:   if (type == PETSC_STRUCT) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
686:   if (type == PETSC_DATATYPE_UNKNOWN) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");

688:   PetscDataTypeGetSize(type, &size);
689:   /* Load a specific data type into data bucket, specifying textual name and its size in bytes */
690:   DMSwarmDataBucketRegisterField(swarm->db,"DMSwarmRegisterPetscDatatypeField",fieldname,blocksize*size,NULL);
691:   {
692:     DMSwarmDataField gfield;

694:     DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield);
695:     DMSwarmDataFieldSetBlockSize(gfield,blocksize);
696:   }
697:   swarm->db->field[swarm->db->nfields-1]->petsc_type = type;
698:   return(0);
699: }

701: /*@C
702:    DMSwarmRegisterUserStructField - Register a user defined struct to a DMSwarm

704:    Collective on dm

706:    Input parameters:
707: +  dm - a DMSwarm
708: .  fieldname - the textual name to identify this field
709: -  size - the size in bytes of the user struct of each data type

711:    Level: beginner

713:    Notes:
714:    The textual name for each registered field must be unique.

716: .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmRegisterUserDatatypeField()
717: @*/
718: PETSC_EXTERN PetscErrorCode DMSwarmRegisterUserStructField(DM dm,const char fieldname[],size_t size)
719: {
721:   DM_Swarm       *swarm = (DM_Swarm*)dm->data;

724:   DMSwarmDataBucketRegisterField(swarm->db,"DMSwarmRegisterUserStructField",fieldname,size,NULL);
725:   swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_STRUCT ;
726:   return(0);
727: }

729: /*@C
730:    DMSwarmRegisterUserDatatypeField - Register a user defined data type to a DMSwarm

732:    Collective on dm

734:    Input parameters:
735: +  dm - a DMSwarm
736: .  fieldname - the textual name to identify this field
737: .  size - the size in bytes of the user data type
738: -  blocksize - the number of each data type

740:    Level: beginner

742:    Notes:
743:    The textual name for each registered field must be unique.

745: .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField()
746: @*/
747: PETSC_EXTERN PetscErrorCode DMSwarmRegisterUserDatatypeField(DM dm,const char fieldname[],size_t size,PetscInt blocksize)
748: {
749:   DM_Swarm       *swarm = (DM_Swarm*)dm->data;

753:   DMSwarmDataBucketRegisterField(swarm->db,"DMSwarmRegisterUserDatatypeField",fieldname,blocksize*size,NULL);
754:   {
755:     DMSwarmDataField gfield;

757:     DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield);
758:     DMSwarmDataFieldSetBlockSize(gfield,blocksize);
759:   }
760:   swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_DATATYPE_UNKNOWN;
761:   return(0);
762: }

764: /*@C
765:    DMSwarmGetField - Get access to the underlying array storing all entries associated with a registered field

767:    Not collective

769:    Input parameters:
770: +  dm - a DMSwarm
771: -  fieldname - the textual name to identify this field

773:    Output parameters:
774: +  blocksize - the number of each data type
775: .  type - the data type
776: -  data - pointer to raw array

778:    Level: beginner

780:    Notes:
781:    The array must be returned using a matching call to DMSwarmRestoreField().

783: .seealso: DMSwarmRestoreField()
784: @*/
785: PETSC_EXTERN PetscErrorCode DMSwarmGetField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data)
786: {
787:   DM_Swarm         *swarm = (DM_Swarm*)dm->data;
788:   DMSwarmDataField gfield;

792:   if (!swarm->issetup) { DMSetUp(dm); }
793:   DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield);
794:   DMSwarmDataFieldGetAccess(gfield);
795:   DMSwarmDataFieldGetEntries(gfield,data);
796:   if (blocksize) {*blocksize = gfield->bs; }
797:   if (type) { *type = gfield->petsc_type; }
798:   return(0);
799: }

801: /*@C
802:    DMSwarmRestoreField - Restore access to the underlying array storing all entries associated with a registered field

804:    Not collective

806:    Input parameters:
807: +  dm - a DMSwarm
808: -  fieldname - the textual name to identify this field

810:    Output parameters:
811: +  blocksize - the number of each data type
812: .  type - the data type
813: -  data - pointer to raw array

815:    Level: beginner

817:    Notes:
818:    The user must call DMSwarmGetField() prior to calling DMSwarmRestoreField().

820: .seealso: DMSwarmGetField()
821: @*/
822: PETSC_EXTERN PetscErrorCode DMSwarmRestoreField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data)
823: {
824:   DM_Swarm         *swarm = (DM_Swarm*)dm->data;
825:   DMSwarmDataField gfield;

829:   DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield);
830:   DMSwarmDataFieldRestoreAccess(gfield);
831:   if (data) *data = NULL;
832:   return(0);
833: }

835: /*@C
836:    DMSwarmAddPoint - Add space for one new point in the DMSwarm

838:    Not collective

840:    Input parameter:
841: .  dm - a DMSwarm

843:    Level: beginner

845:    Notes:
846:    The new point will have all fields initialized to zero.

848: .seealso: DMSwarmAddNPoints()
849: @*/
850: PETSC_EXTERN PetscErrorCode DMSwarmAddPoint(DM dm)
851: {
852:   DM_Swarm       *swarm = (DM_Swarm*)dm->data;

856:   if (!swarm->issetup) {DMSetUp(dm);}
857:   PetscLogEventBegin(DMSWARM_AddPoints,0,0,0,0);
858:   DMSwarmDataBucketAddPoint(swarm->db);
859:   PetscLogEventEnd(DMSWARM_AddPoints,0,0,0,0);
860:   return(0);
861: }

863: /*@C
864:    DMSwarmAddNPoints - Add space for a number of new points in the DMSwarm

866:    Not collective

868:    Input parameters:
869: +  dm - a DMSwarm
870: -  npoints - the number of new points to add

872:    Level: beginner

874:    Notes:
875:    The new point will have all fields initialized to zero.

877: .seealso: DMSwarmAddPoint()
878: @*/
879: PETSC_EXTERN PetscErrorCode DMSwarmAddNPoints(DM dm,PetscInt npoints)
880: {
881:   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
883:   PetscInt       nlocal;

886:   PetscLogEventBegin(DMSWARM_AddPoints,0,0,0,0);
887:   DMSwarmDataBucketGetSizes(swarm->db,&nlocal,NULL,NULL);
888:   nlocal = nlocal + npoints;
889:   DMSwarmDataBucketSetSizes(swarm->db,nlocal,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT);
890:   PetscLogEventEnd(DMSWARM_AddPoints,0,0,0,0);
891:   return(0);
892: }

894: /*@C
895:    DMSwarmRemovePoint - Remove the last point from the DMSwarm

897:    Not collective

899:    Input parameter:
900: .  dm - a DMSwarm

902:    Level: beginner

904: .seealso: DMSwarmRemovePointAtIndex()
905: @*/
906: PETSC_EXTERN PetscErrorCode DMSwarmRemovePoint(DM dm)
907: {
908:   DM_Swarm       *swarm = (DM_Swarm*)dm->data;

912:   PetscLogEventBegin(DMSWARM_RemovePoints,0,0,0,0);
913:   DMSwarmDataBucketRemovePoint(swarm->db);
914:   PetscLogEventEnd(DMSWARM_RemovePoints,0,0,0,0);
915:   return(0);
916: }

918: /*@C
919:    DMSwarmRemovePointAtIndex - Removes a specific point from the DMSwarm

921:    Not collective

923:    Input parameters:
924: +  dm - a DMSwarm
925: -  idx - index of point to remove

927:    Level: beginner

929: .seealso: DMSwarmRemovePoint()
930: @*/
931: PETSC_EXTERN PetscErrorCode DMSwarmRemovePointAtIndex(DM dm,PetscInt idx)
932: {
933:   DM_Swarm       *swarm = (DM_Swarm*)dm->data;

937:   PetscLogEventBegin(DMSWARM_RemovePoints,0,0,0,0);
938:   DMSwarmDataBucketRemovePointAtIndex(swarm->db,idx);
939:   PetscLogEventEnd(DMSWARM_RemovePoints,0,0,0,0);
940:   return(0);
941: }

943: /*@C
944:    DMSwarmCopyPoint - Copy point pj to point pi in the DMSwarm

946:    Not collective

948:    Input parameters:
949: +  dm - a DMSwarm
950: .  pi - the index of the point to copy
951: -  pj - the point index where the copy should be located

953:  Level: beginner

955: .seealso: DMSwarmRemovePoint()
956: @*/
957: PETSC_EXTERN PetscErrorCode DMSwarmCopyPoint(DM dm,PetscInt pi,PetscInt pj)
958: {
959:   DM_Swarm       *swarm = (DM_Swarm*)dm->data;

963:   if (!swarm->issetup) {DMSetUp(dm);}
964:   DMSwarmDataBucketCopyPoint(swarm->db,pi,swarm->db,pj);
965:   return(0);
966: }

968: PetscErrorCode DMSwarmMigrate_Basic(DM dm,PetscBool remove_sent_points)
969: {

973:   DMSwarmMigrate_Push_Basic(dm,remove_sent_points);
974:   return(0);
975: }

977: /*@C
978:    DMSwarmMigrate - Relocates points defined in the DMSwarm to other MPI-ranks

980:    Collective on dm

982:    Input parameters:
983: +  dm - the DMSwarm
984: -  remove_sent_points - flag indicating if sent points should be removed from the current MPI-rank

986:    Notes:
987:    The DM will be modified to accomodate received points.
988:    If remove_sent_points = PETSC_TRUE, any points that were sent will be removed from the DM.
989:    Different styles of migration are supported. See DMSwarmSetMigrateType().

991:    Level: advanced

993: .seealso: DMSwarmSetMigrateType()
994: @*/
995: PETSC_EXTERN PetscErrorCode DMSwarmMigrate(DM dm,PetscBool remove_sent_points)
996: {
997:   DM_Swarm       *swarm = (DM_Swarm*)dm->data;

1001:   PetscLogEventBegin(DMSWARM_Migrate,0,0,0,0);
1002:   switch (swarm->migrate_type) {
1003:     case DMSWARM_MIGRATE_BASIC:
1004:       DMSwarmMigrate_Basic(dm,remove_sent_points);
1005:       break;
1006:     case DMSWARM_MIGRATE_DMCELLNSCATTER:
1007:       DMSwarmMigrate_CellDMScatter(dm,remove_sent_points);
1008:       break;
1009:     case DMSWARM_MIGRATE_DMCELLEXACT:
1010:       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_MIGRATE_DMCELLEXACT not implemented");
1011:       break;
1012:     case DMSWARM_MIGRATE_USER:
1013:       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_MIGRATE_USER not implemented");
1014:       break;
1015:     default:
1016:       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_MIGRATE type unknown");
1017:       break;
1018:   }
1019:   PetscLogEventEnd(DMSWARM_Migrate,0,0,0,0);
1020:   return(0);
1021: }

1023: PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm,PetscInt *globalsize);

1025: /*
1026:  DMSwarmCollectViewCreate

1028:  * Applies a collection method and gathers point neighbour points into dm

1030:  Notes:
1031:  Users should call DMSwarmCollectViewDestroy() after
1032:  they have finished computations associated with the collected points
1033: */

1035: /*@C
1036:    DMSwarmCollectViewCreate - Applies a collection method and gathers points
1037:    in neighbour MPI-ranks into the DMSwarm

1039:    Collective on dm

1041:    Input parameter:
1042: .  dm - the DMSwarm

1044:    Notes:
1045:    Users should call DMSwarmCollectViewDestroy() after
1046:    they have finished computations associated with the collected points
1047:    Different collect methods are supported. See DMSwarmSetCollectType().

1049:    Level: advanced

1051: .seealso: DMSwarmCollectViewDestroy(), DMSwarmSetCollectType()
1052: @*/
1053: PETSC_EXTERN PetscErrorCode DMSwarmCollectViewCreate(DM dm)
1054: {
1056:   DM_Swarm *swarm = (DM_Swarm*)dm->data;
1057:   PetscInt ng;

1060:   if (swarm->collect_view_active) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"CollectView currently active");
1061:   DMSwarmGetLocalSize(dm,&ng);
1062:   switch (swarm->collect_type) {

1064:     case DMSWARM_COLLECT_BASIC:
1065:       DMSwarmMigrate_GlobalToLocal_Basic(dm,&ng);
1066:       break;
1067:     case DMSWARM_COLLECT_DMDABOUNDINGBOX:
1068:       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_COLLECT_DMDABOUNDINGBOX not implemented");
1069:       break;
1070:     case DMSWARM_COLLECT_GENERAL:
1071:       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_COLLECT_GENERAL not implemented");
1072:       break;
1073:     default:
1074:       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_COLLECT type unknown");
1075:       break;
1076:   }
1077:   swarm->collect_view_active = PETSC_TRUE;
1078:   swarm->collect_view_reset_nlocal = ng;
1079:   return(0);
1080: }

1082: /*@C
1083:    DMSwarmCollectViewDestroy - Resets the DMSwarm to the size prior to calling DMSwarmCollectViewCreate()

1085:    Collective on dm

1087:    Input parameters:
1088: .  dm - the DMSwarm

1090:    Notes:
1091:    Users should call DMSwarmCollectViewCreate() before this function is called.

1093:    Level: advanced

1095: .seealso: DMSwarmCollectViewCreate(), DMSwarmSetCollectType()
1096: @*/
1097: PETSC_EXTERN PetscErrorCode DMSwarmCollectViewDestroy(DM dm)
1098: {
1100:   DM_Swarm       *swarm = (DM_Swarm*)dm->data;

1103:   if (!swarm->collect_view_active) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"CollectView is currently not active");
1104:   DMSwarmSetLocalSizes(dm,swarm->collect_view_reset_nlocal,-1);
1105:   swarm->collect_view_active = PETSC_FALSE;
1106:   return(0);
1107: }

1109: PetscErrorCode DMSwarmSetUpPIC(DM dm)
1110: {
1111:   PetscInt       dim;

1115:   DMGetDimension(dm,&dim);
1116:   if (dim < 1) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Dimension must be 1,2,3 - found %D",dim);
1117:   if (dim > 3) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Dimension must be 1,2,3 - found %D",dim);
1118:   DMSwarmRegisterPetscDatatypeField(dm,DMSwarmPICField_coor,dim,PETSC_DOUBLE);
1119:   DMSwarmRegisterPetscDatatypeField(dm,DMSwarmPICField_cellid,1,PETSC_INT);
1120:   return(0);
1121: }

1123: /*@C
1124:    DMSwarmSetType - Set particular flavor of DMSwarm

1126:    Collective on dm

1128:    Input parameters:
1129: +  dm - the DMSwarm
1130: -  stype - the DMSwarm type (e.g. DMSWARM_PIC)

1132:    Level: advanced

1134: .seealso: DMSwarmSetMigrateType(), DMSwarmSetCollectType()
1135: @*/
1136: PETSC_EXTERN PetscErrorCode DMSwarmSetType(DM dm,DMSwarmType stype)
1137: {
1138:   DM_Swarm       *swarm = (DM_Swarm*)dm->data;

1142:   swarm->swarm_type = stype;
1143:   if (swarm->swarm_type == DMSWARM_PIC) {
1144:     DMSwarmSetUpPIC(dm);
1145:   }
1146:   return(0);
1147: }

1149: PetscErrorCode DMSetup_Swarm(DM dm)
1150: {
1151:   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1153:   PetscMPIInt    rank;
1154:   PetscInt       p,npoints,*rankval;

1157:   if (swarm->issetup) return(0);

1159:   swarm->issetup = PETSC_TRUE;

1161:   if (swarm->swarm_type == DMSWARM_PIC) {
1162:     /* check dmcell exists */
1163:     if (!swarm->dmcell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSWARM_PIC requires you call DMSwarmSetCellDM");

1165:     if (swarm->dmcell->ops->locatepointssubdomain) {
1166:       /* check methods exists for exact ownership identificiation */
1167:       PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->ops->LocatePointsSubdomain\n");
1168:       swarm->migrate_type = DMSWARM_MIGRATE_DMCELLEXACT;
1169:     } else {
1170:       /* check methods exist for point location AND rank neighbor identification */
1171:       if (swarm->dmcell->ops->locatepoints) {
1172:         PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->LocatePoints\n");
1173:       } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSWARM_PIC requires the method CellDM->ops->locatepoints be defined");

1175:       if (swarm->dmcell->ops->getneighbors) {
1176:         PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->GetNeigbors\n");
1177:       } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSWARM_PIC requires the method CellDM->ops->getneighbors be defined");

1179:       swarm->migrate_type = DMSWARM_MIGRATE_DMCELLNSCATTER;
1180:     }
1181:   }

1183:   DMSwarmFinalizeFieldRegister(dm);

1185:   /* check some fields were registered */
1186:   if (swarm->db->nfields <= 2) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"At least one field user must be registered via DMSwarmRegisterXXX()");

1188:   /* check local sizes were set */
1189:   if (swarm->db->L == -1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Local sizes must be set via DMSwarmSetLocalSizes()");

1191:   /* initialize values in pid and rank placeholders */
1192:   /* TODO: [pid - use MPI_Scan] */
1193:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);
1194:   DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);
1195:   DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);
1196:   for (p=0; p<npoints; p++) {
1197:     rankval[p] = (PetscInt)rank;
1198:   }
1199:   DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);
1200:   return(0);
1201: }

1203: extern PetscErrorCode DMSwarmSortDestroy(DMSwarmSort *_ctx);

1205: PetscErrorCode DMDestroy_Swarm(DM dm)
1206: {
1207:   DM_Swarm       *swarm = (DM_Swarm*)dm->data;

1211:   DMSwarmDataBucketDestroy(&swarm->db);
1212:   if (swarm->sort_context) {
1213:     DMSwarmSortDestroy(&swarm->sort_context);
1214:   }
1215:   PetscFree(swarm);
1216:   return(0);
1217: }

1219: PetscErrorCode DMSwarmView_Draw(DM dm, PetscViewer viewer)
1220: {
1221:   DM             cdm;
1222:   PetscDraw      draw;
1223:   PetscReal     *coords, oldPause;
1224:   PetscInt       Np, p, bs;

1228:   PetscViewerDrawGetDraw(viewer, 0, &draw);
1229:   DMSwarmGetCellDM(dm, &cdm);
1230:   PetscDrawGetPause(draw, &oldPause);
1231:   PetscDrawSetPause(draw, 0.0);
1232:   DMView(cdm, viewer);
1233:   PetscDrawSetPause(draw, oldPause);

1235:   DMSwarmGetLocalSize(dm, &Np);
1236:   DMSwarmGetField(dm, DMSwarmPICField_coor, &bs, NULL, (void **) &coords);
1237:   for (p = 0; p < Np; ++p) {
1238:     const PetscInt i = p*bs;

1240:     PetscDrawEllipse(draw, coords[i], coords[i+1], 0.01, 0.01, PETSC_DRAW_BLUE);
1241:   }
1242:   DMSwarmRestoreField(dm, DMSwarmPICField_coor, &bs, NULL, (void **) &coords);
1243:   PetscDrawFlush(draw);
1244:   PetscDrawPause(draw);
1245:   PetscDrawSave(draw);
1246:   return(0);
1247: }

1249: PetscErrorCode DMView_Swarm(DM dm, PetscViewer viewer)
1250: {
1251:   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
1252:   PetscBool      iascii,ibinary,ishdf5,isvtk,isdraw;

1258:   PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii);
1259:   PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY,&ibinary);
1260:   PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK,   &isvtk);
1261:   PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5,  &ishdf5);
1262:   PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW,  &isdraw);
1263:   if (iascii) {
1264:     DMSwarmDataBucketView(PetscObjectComm((PetscObject)dm),swarm->db,NULL,DATABUCKET_VIEW_STDOUT);
1265:   } else if (ibinary) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO Binary support");
1266: #if defined(PETSC_HAVE_HDF5)
1267:   else if (ishdf5) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO HDF5 support");
1268: #else
1269:   else if (ishdf5) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"HDF5 not supported. Please reconfigure using --download-hdf5");
1270: #endif
1271:   else if (isvtk) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO VTK support");
1272:   else if (isdraw) {
1273:     DMSwarmView_Draw(dm, viewer);
1274:   }
1275:   return(0);
1276: }

1278: /*MC

1280:  DMSWARM = "swarm" - A DM object used to represent arrays of data (fields) of arbitrary data type.
1281:  This implementation was designed for particle methods in which the underlying
1282:  data required to be represented is both (i) dynamic in length, (ii) and of arbitrary data type.

1284:  User data can be represented by DMSwarm through a registering "fields".
1285:  To register a field, the user must provide;
1286:  (a) a unique name;
1287:  (b) the data type (or size in bytes);
1288:  (c) the block size of the data.

1290:  For example, suppose the application requires a unique id, energy, momentum and density to be stored
1291:  on a set of particles. Then the following code could be used

1293: $    DMSwarmInitializeFieldRegister(dm)
1294: $    DMSwarmRegisterPetscDatatypeField(dm,"uid",1,PETSC_LONG);
1295: $    DMSwarmRegisterPetscDatatypeField(dm,"energy",1,PETSC_REAL);
1296: $    DMSwarmRegisterPetscDatatypeField(dm,"momentum",3,PETSC_REAL);
1297: $    DMSwarmRegisterPetscDatatypeField(dm,"density",1,PETSC_FLOAT);
1298: $    DMSwarmFinalizeFieldRegister(dm)

1300:  The fields represented by DMSwarm are dynamic and can be re-sized at any time.
1301:  The only restriction imposed by DMSwarm is that all fields contain the same number of points.

1303:  To support particle methods, "migration" techniques are provided. These methods migrate data
1304:  between MPI-ranks.

1306:  DMSwarm supports the methods DMCreateGlobalVector() and DMCreateLocalVector().
1307:  As a DMSwarm may internally define and store values of different data types,
1308:  before calling DMCreateGlobalVector() or DMCreateLocalVector(), the user must inform DMSwarm which
1309:  fields should be used to define a Vec object via
1310:    DMSwarmVectorDefineField()
1311:  The specified field can be changed at any time - thereby permitting vectors
1312:  compatible with different fields to be created.

1314:  A dual representation of fields in the DMSwarm and a Vec object is permitted via
1315:    DMSwarmCreateGlobalVectorFromField()
1316:  Here the data defining the field in the DMSwarm is shared with a Vec.
1317:  This is inherently unsafe if you alter the size of the field at any time between
1318:  calls to DMSwarmCreateGlobalVectorFromField() and DMSwarmDestroyGlobalVectorFromField().
1319:  If the local size of the DMSwarm does not match the local size of the global vector
1320:  when DMSwarmDestroyGlobalVectorFromField() is called, an error is thrown.

1322:  Additional high-level support is provided for Particle-In-Cell methods.
1323:  Please refer to the man page for DMSwarmSetType().

1325:  Level: beginner

1327: .seealso: DMType, DMCreate(), DMSetType()
1328: M*/
1329: PETSC_EXTERN PetscErrorCode DMCreate_Swarm(DM dm)
1330: {
1331:   DM_Swarm      *swarm;

1336:   PetscNewLog(dm,&swarm);
1337:   dm->data = swarm;

1339:   DMSwarmDataBucketCreate(&swarm->db);
1340:   DMSwarmInitializeFieldRegister(dm);

1342:   swarm->vec_field_set = PETSC_FALSE;
1343:   swarm->issetup = PETSC_FALSE;
1344:   swarm->swarm_type = DMSWARM_BASIC;
1345:   swarm->migrate_type = DMSWARM_MIGRATE_BASIC;
1346:   swarm->collect_type = DMSWARM_COLLECT_BASIC;
1347:   swarm->migrate_error_on_missing_point = PETSC_FALSE;

1349:   swarm->dmcell = NULL;
1350:   swarm->collect_view_active = PETSC_FALSE;
1351:   swarm->collect_view_reset_nlocal = -1;

1353:   dm->dim  = 0;
1354:   dm->ops->view                            = DMView_Swarm;
1355:   dm->ops->load                            = NULL;
1356:   dm->ops->setfromoptions                  = NULL;
1357:   dm->ops->clone                           = NULL;
1358:   dm->ops->setup                           = DMSetup_Swarm;
1359:   dm->ops->createlocalsection              = NULL;
1360:   dm->ops->createdefaultconstraints        = NULL;
1361:   dm->ops->createglobalvector              = DMCreateGlobalVector_Swarm;
1362:   dm->ops->createlocalvector               = DMCreateLocalVector_Swarm;
1363:   dm->ops->getlocaltoglobalmapping         = NULL;
1364:   dm->ops->createfieldis                   = NULL;
1365:   dm->ops->createcoordinatedm              = NULL;
1366:   dm->ops->getcoloring                     = NULL;
1367:   dm->ops->creatematrix                    = NULL;
1368:   dm->ops->createinterpolation             = NULL;
1369:   dm->ops->createinjection                 = NULL;
1370:   dm->ops->createmassmatrix                = DMCreateMassMatrix_Swarm;
1371:   dm->ops->refine                          = NULL;
1372:   dm->ops->coarsen                         = NULL;
1373:   dm->ops->refinehierarchy                 = NULL;
1374:   dm->ops->coarsenhierarchy                = NULL;
1375:   dm->ops->globaltolocalbegin              = NULL;
1376:   dm->ops->globaltolocalend                = NULL;
1377:   dm->ops->localtoglobalbegin              = NULL;
1378:   dm->ops->localtoglobalend                = NULL;
1379:   dm->ops->destroy                         = DMDestroy_Swarm;
1380:   dm->ops->createsubdm                     = NULL;
1381:   dm->ops->getdimpoints                    = NULL;
1382:   dm->ops->locatepoints                    = NULL;
1383:   return(0);
1384: }