Actual source code: swarm.c
petsc-master 2020-08-25
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", NULL };
15: const char* DMSwarmMigrateTypeNames[] = { "basic", "dmcellnscatter", "dmcellexact", "user", NULL };
16: const char* DMSwarmCollectTypeNames[] = { "basic", "boundingbox", "general", "user", NULL };
17: const char* DMSwarmPICLayoutTypeNames[] = { "regular", "gauss", "subdivision", NULL };
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: PETSC_EXTERN PetscErrorCode VecView_Seq(Vec, PetscViewer);
25: PETSC_EXTERN PetscErrorCode VecView_MPI(Vec, PetscViewer);
27: #if defined(PETSC_HAVE_HDF5)
28: #include <petscviewerhdf5.h>
30: PetscErrorCode VecView_Swarm_HDF5_Internal(Vec v, PetscViewer viewer)
31: {
32: DM dm;
33: PetscReal seqval;
34: PetscInt seqnum, bs;
35: PetscBool isseq;
39: VecGetDM(v, &dm);
40: VecGetBlockSize(v, &bs);
41: PetscViewerHDF5PushGroup(viewer, "/particle_fields");
42: PetscObjectTypeCompare((PetscObject) v, VECSEQ, &isseq);
43: DMGetOutputSequenceNumber(dm, &seqnum, &seqval);
44: PetscViewerHDF5SetTimestep(viewer, seqnum);
45: /* DMSequenceView_HDF5(dm, "time", seqnum, (PetscScalar) seqval, viewer); */
46: if (isseq) {VecView_Seq(v, viewer);}
47: else {VecView_MPI(v, viewer);}
48: PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject) v, "Nc", PETSC_INT, (void *) &bs);
49: PetscViewerHDF5PopGroup(viewer);
50: return(0);
51: }
53: PetscErrorCode DMSwarmView_HDF5(DM dm, PetscViewer viewer)
54: {
55: Vec coordinates;
56: PetscInt Np;
57: PetscBool isseq;
61: DMSwarmGetSize(dm, &Np);
62: DMSwarmCreateGlobalVectorFromField(dm, DMSwarmPICField_coor, &coordinates);
63: PetscObjectSetName((PetscObject) coordinates, "coordinates");
64: PetscViewerHDF5PushGroup(viewer, "/particles");
65: PetscObjectTypeCompare((PetscObject) coordinates, VECSEQ, &isseq);
66: if (isseq) {VecView_Seq(coordinates, viewer);}
67: else {VecView_MPI(coordinates, viewer);}
68: PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject) coordinates, "Np", PETSC_INT, (void *) &Np);
69: PetscViewerHDF5PopGroup(viewer);
70: DMSwarmDestroyGlobalVectorFromField(dm, DMSwarmPICField_coor, &coordinates);
71: return(0);
72: }
73: #endif
75: PetscErrorCode VecView_Swarm(Vec v, PetscViewer viewer)
76: {
77: DM dm;
78: PetscBool ishdf5;
82: VecGetDM(v, &dm);
83: if (!dm) SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
84: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5);
85: if (ishdf5) {
86: #if defined(PETSC_HAVE_HDF5)
87: VecView_Swarm_HDF5_Internal(v, viewer);
88: #else
89: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
90: #endif
91: } else {
92: PetscBool isseq;
94: PetscObjectTypeCompare((PetscObject) v, VECSEQ, &isseq);
95: if (isseq) {VecView_Seq(v, viewer);}
96: else {VecView_MPI(v, viewer);}
97: }
98: return(0);
99: }
101: /*@C
102: DMSwarmVectorDefineField - Sets the field from which to define a Vec object
103: when DMCreateLocalVector(), or DMCreateGlobalVector() is called
105: Collective on dm
107: Input parameters:
108: + dm - a DMSwarm
109: - fieldname - the textual name given to a registered field
111: Level: beginner
113: Notes:
115: The field with name fieldname must be defined as having a data type of PetscScalar.
117: This function must be called prior to calling DMCreateLocalVector(), DMCreateGlobalVector().
118: Mutiple calls to DMSwarmVectorDefineField() are permitted.
120: .seealso: DMSwarmRegisterPetscDatatypeField(), DMCreateGlobalVector(), DMCreateLocalVector()
121: @*/
122: PetscErrorCode DMSwarmVectorDefineField(DM dm,const char fieldname[])
123: {
124: DM_Swarm *swarm = (DM_Swarm*)dm->data;
126: PetscInt bs,n;
127: PetscScalar *array;
128: PetscDataType type;
131: if (!swarm->issetup) { DMSetUp(dm); }
132: DMSwarmDataBucketGetSizes(swarm->db,&n,NULL,NULL);
133: DMSwarmGetField(dm,fieldname,&bs,&type,(void**)&array);
135: /* Check all fields are of type PETSC_REAL or PETSC_SCALAR */
136: if (type != PETSC_REAL) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid for PETSC_REAL");
137: PetscSNPrintf(swarm->vec_field_name,PETSC_MAX_PATH_LEN-1,"%s",fieldname);
138: swarm->vec_field_set = PETSC_TRUE;
139: swarm->vec_field_bs = bs;
140: swarm->vec_field_nlocal = n;
141: DMSwarmRestoreField(dm,fieldname,&bs,&type,(void**)&array);
142: return(0);
143: }
145: /* requires DMSwarmDefineFieldVector has been called */
146: PetscErrorCode DMCreateGlobalVector_Swarm(DM dm,Vec *vec)
147: {
148: DM_Swarm *swarm = (DM_Swarm*)dm->data;
150: Vec x;
151: char name[PETSC_MAX_PATH_LEN];
154: if (!swarm->issetup) { DMSetUp(dm); }
155: if (!swarm->vec_field_set) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmVectorDefineField first");
156: 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 */
158: PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarmField_%s",swarm->vec_field_name);
159: VecCreate(PetscObjectComm((PetscObject)dm),&x);
160: PetscObjectSetName((PetscObject)x,name);
161: VecSetSizes(x,swarm->db->L*swarm->vec_field_bs,PETSC_DETERMINE);
162: VecSetBlockSize(x,swarm->vec_field_bs);
163: VecSetDM(x,dm);
164: VecSetFromOptions(x);
165: VecSetDM(x, dm);
166: VecSetOperation(x, VECOP_VIEW, (void (*)(void)) VecView_Swarm);
167: *vec = x;
168: return(0);
169: }
171: /* requires DMSwarmDefineFieldVector has been called */
172: PetscErrorCode DMCreateLocalVector_Swarm(DM dm,Vec *vec)
173: {
174: DM_Swarm *swarm = (DM_Swarm*)dm->data;
176: Vec x;
177: char name[PETSC_MAX_PATH_LEN];
180: if (!swarm->issetup) { DMSetUp(dm); }
181: if (!swarm->vec_field_set) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmVectorDefineField first");
182: 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 */
184: PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarmField_%s",swarm->vec_field_name);
185: VecCreate(PETSC_COMM_SELF,&x);
186: PetscObjectSetName((PetscObject)x,name);
187: VecSetSizes(x,swarm->db->L*swarm->vec_field_bs,PETSC_DETERMINE);
188: VecSetBlockSize(x,swarm->vec_field_bs);
189: VecSetDM(x,dm);
190: VecSetFromOptions(x);
191: *vec = x;
192: return(0);
193: }
195: static PetscErrorCode DMSwarmDestroyVectorFromField_Private(DM dm, const char fieldname[], Vec *vec)
196: {
197: DM_Swarm *swarm = (DM_Swarm *) dm->data;
198: DMSwarmDataField gfield;
199: void (*fptr)(void);
200: PetscInt bs, nlocal;
201: char name[PETSC_MAX_PATH_LEN];
205: VecGetLocalSize(*vec, &nlocal);
206: VecGetBlockSize(*vec, &bs);
207: 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 */
208: DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield);
209: /* check vector is an inplace array */
210: PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "DMSwarm_VecFieldInPlace_%s", fieldname);
211: PetscObjectQueryFunction((PetscObject) *vec, name, &fptr);
212: if (!fptr) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_USER, "Vector being destroyed was not created from DMSwarm field(%s)", fieldname);
213: DMSwarmDataFieldRestoreAccess(gfield);
214: VecDestroy(vec);
215: return(0);
216: }
218: static PetscErrorCode DMSwarmCreateVectorFromField_Private(DM dm, const char fieldname[], MPI_Comm comm, Vec *vec)
219: {
220: DM_Swarm *swarm = (DM_Swarm *) dm->data;
221: PetscDataType type;
222: PetscScalar *array;
223: PetscInt bs, n;
224: char name[PETSC_MAX_PATH_LEN];
225: PetscMPIInt size;
229: if (!swarm->issetup) {DMSetUp(dm);}
230: DMSwarmDataBucketGetSizes(swarm->db, &n, NULL, NULL);
231: DMSwarmGetField(dm, fieldname, &bs, &type, (void **) &array);
232: if (type != PETSC_REAL) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Only valid for PETSC_REAL");
234: MPI_Comm_size(comm, &size);
235: if (size == 1) {
236: VecCreateSeqWithArray(comm, bs, n*bs, array, vec);
237: } else {
238: VecCreateMPIWithArray(comm, bs, n*bs, PETSC_DETERMINE, array, vec);
239: }
240: PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "DMSwarmSharedField_%s", fieldname);
241: PetscObjectSetName((PetscObject) *vec, name);
243: /* Set guard */
244: PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "DMSwarm_VecFieldInPlace_%s", fieldname);
245: PetscObjectComposeFunction((PetscObject) *vec, name, DMSwarmDestroyVectorFromField_Private);
247: VecSetDM(*vec, dm);
248: VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Swarm);
249: return(0);
250: }
252: /* This creates a "mass matrix" between a finite element and particle space. If a finite element interpolant is given by
254: \hat f = \sum_i f_i \phi_i
256: and a particle function is given by
258: f = \sum_i w_i \delta(x - x_i)
260: then we want to require that
262: M \hat f = M_p f
264: where the particle mass matrix is given by
266: (M_p)_{ij} = \int \phi_i \delta(x - x_j)
268: The way Dave May does particles, they amount to quadratue weights rather than delta functions, so he has |J| is in
269: his integral. We allow this with the boolean flag.
270: */
271: static PetscErrorCode DMSwarmComputeMassMatrix_Private(DM dmc, DM dmf, Mat mass, PetscBool useDeltaFunction, void *ctx)
272: {
273: const char *name = "Mass Matrix";
274: MPI_Comm comm;
275: PetscDS prob;
276: PetscSection fsection, globalFSection;
277: PetscHSetIJ ht;
278: PetscLayout rLayout, colLayout;
279: PetscInt *dnz, *onz;
280: PetscInt locRows, locCols, rStart, colStart, colEnd, *rowIDXs;
281: PetscReal *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0};
282: PetscScalar *elemMat;
283: PetscInt dim, Nf, field, cStart, cEnd, cell, totDim, maxC = 0;
287: PetscObjectGetComm((PetscObject) mass, &comm);
288: DMGetCoordinateDim(dmf, &dim);
289: DMGetDS(dmf, &prob);
290: PetscDSGetNumFields(prob, &Nf);
291: PetscDSGetTotalDimension(prob, &totDim);
292: PetscMalloc3(dim, &v0, dim*dim, &J, dim*dim,&invJ);
293: DMGetLocalSection(dmf, &fsection);
294: DMGetGlobalSection(dmf, &globalFSection);
295: DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd);
296: MatGetLocalSize(mass, &locRows, &locCols);
298: PetscLayoutCreate(comm, &colLayout);
299: PetscLayoutSetLocalSize(colLayout, locCols);
300: PetscLayoutSetBlockSize(colLayout, 1);
301: PetscLayoutSetUp(colLayout);
302: PetscLayoutGetRange(colLayout, &colStart, &colEnd);
303: PetscLayoutDestroy(&colLayout);
305: PetscLayoutCreate(comm, &rLayout);
306: PetscLayoutSetLocalSize(rLayout, locRows);
307: PetscLayoutSetBlockSize(rLayout, 1);
308: PetscLayoutSetUp(rLayout);
309: PetscLayoutGetRange(rLayout, &rStart, NULL);
310: PetscLayoutDestroy(&rLayout);
312: PetscCalloc2(locRows, &dnz, locRows, &onz);
313: PetscHSetIJCreate(&ht);
315: PetscSynchronizedFlush(comm, NULL);
316: /* count non-zeros */
317: DMSwarmSortGetAccess(dmc);
318: for (field = 0; field < Nf; ++field) {
319: for (cell = cStart; cell < cEnd; ++cell) {
320: PetscInt c, i;
321: PetscInt *findices, *cindices; /* fine is vertices, coarse is particles */
322: PetscInt numFIndices, numCIndices;
324: DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL);
325: DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices);
326: maxC = PetscMax(maxC, numCIndices);
327: {
328: PetscHashIJKey key;
329: PetscBool missing;
330: for (i = 0; i < numFIndices; ++i) {
331: key.j = findices[i]; /* global column (from Plex) */
332: if (key.j >= 0) {
333: /* Get indices for coarse elements */
334: for (c = 0; c < numCIndices; ++c) {
335: key.i = cindices[c] + rStart; /* global cols (from Swarm) */
336: if (key.i < 0) continue;
337: PetscHSetIJQueryAdd(ht, key, &missing);
338: if (missing) {
339: if ((key.j >= colStart) && (key.j < colEnd)) ++dnz[key.i - rStart];
340: else ++onz[key.i - rStart];
341: } else SETERRQ2(PetscObjectComm((PetscObject) dmf), PETSC_ERR_SUP, "Set new value at %D,%D", key.i, key.j);
342: }
343: }
344: }
345: PetscFree(cindices);
346: }
347: DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL);
348: }
349: }
350: PetscHSetIJDestroy(&ht);
351: MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL);
352: MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE);
353: PetscFree2(dnz, onz);
354: PetscMalloc3(maxC*totDim, &elemMat, maxC, &rowIDXs, maxC*dim, &xi);
355: for (field = 0; field < Nf; ++field) {
356: PetscTabulation Tcoarse;
357: PetscObject obj;
358: PetscReal *coords;
359: PetscInt Nc, i;
361: PetscDSGetDiscretization(prob, field, &obj);
362: PetscFEGetNumComponents((PetscFE) obj, &Nc);
363: if (Nc != 1) SETERRQ1(PetscObjectComm((PetscObject) dmf), PETSC_ERR_SUP, "Can only interpolate a scalar field from particles, Nc = %D", Nc);
364: DMSwarmGetField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);
365: for (cell = cStart; cell < cEnd; ++cell) {
366: PetscInt *findices , *cindices;
367: PetscInt numFIndices, numCIndices;
368: PetscInt p, c;
370: /* TODO: Use DMField instead of assuming affine */
371: DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);
372: DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL);
373: DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices);
374: for (p = 0; p < numCIndices; ++p) {
375: CoordinatesRealToRef(dim, dim, v0ref, v0, invJ, &coords[cindices[p]*dim], &xi[p*dim]);
376: }
377: PetscFECreateTabulation((PetscFE) obj, 1, numCIndices, xi, 0, &Tcoarse);
378: /* Get elemMat entries by multiplying by weight */
379: PetscArrayzero(elemMat, numCIndices*totDim);
380: for (i = 0; i < numFIndices; ++i) {
381: for (p = 0; p < numCIndices; ++p) {
382: for (c = 0; c < Nc; ++c) {
383: /* B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c */
384: elemMat[p*numFIndices+i] += Tcoarse->T[0][(p*numFIndices + i)*Nc + c]*(useDeltaFunction ? 1.0 : detJ);
385: }
386: }
387: }
388: for (p = 0; p < numCIndices; ++p) rowIDXs[p] = cindices[p] + rStart;
389: if (0) {DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat);}
390: MatSetValues(mass, numCIndices, rowIDXs, numFIndices, findices, elemMat, ADD_VALUES);
391: PetscFree(cindices);
392: DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL);
393: PetscTabulationDestroy(&Tcoarse);
394: }
395: DMSwarmRestoreField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);
396: }
397: PetscFree3(elemMat, rowIDXs, xi);
398: DMSwarmSortRestoreAccess(dmc);
399: PetscFree3(v0, J, invJ);
400: MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY);
401: MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY);
402: return(0);
403: }
405: /* FEM cols, Particle rows */
406: static PetscErrorCode DMCreateMassMatrix_Swarm(DM dmCoarse, DM dmFine, Mat *mass)
407: {
408: PetscSection gsf;
409: PetscInt m, n;
410: void *ctx;
414: DMGetGlobalSection(dmFine, &gsf);
415: PetscSectionGetConstrainedStorageSize(gsf, &m);
416: DMSwarmGetLocalSize(dmCoarse, &n);
417: MatCreate(PetscObjectComm((PetscObject) dmCoarse), mass);
418: MatSetSizes(*mass, n, m, PETSC_DETERMINE, PETSC_DETERMINE);
419: MatSetType(*mass, dmCoarse->mattype);
420: DMGetApplicationContext(dmFine, &ctx);
422: DMSwarmComputeMassMatrix_Private(dmCoarse, dmFine, *mass, PETSC_TRUE, ctx);
423: MatViewFromOptions(*mass, NULL, "-mass_mat_view");
424: return(0);
425: }
427: /*@C
428: DMSwarmCreateGlobalVectorFromField - Creates a Vec object sharing the array associated with a given field
430: Collective on dm
432: Input parameters:
433: + dm - a DMSwarm
434: - fieldname - the textual name given to a registered field
436: Output parameter:
437: . vec - the vector
439: Level: beginner
441: Notes:
442: The vector must be returned using a matching call to DMSwarmDestroyGlobalVectorFromField().
444: .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmDestroyGlobalVectorFromField()
445: @*/
446: PetscErrorCode DMSwarmCreateGlobalVectorFromField(DM dm,const char fieldname[],Vec *vec)
447: {
448: MPI_Comm comm = PetscObjectComm((PetscObject) dm);
452: DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec);
453: return(0);
454: }
456: /*@C
457: DMSwarmDestroyGlobalVectorFromField - Destroys the Vec object which share the array associated with a given field
459: Collective on dm
461: Input parameters:
462: + dm - a DMSwarm
463: - fieldname - the textual name given to a registered field
465: Output parameter:
466: . vec - the vector
468: Level: beginner
470: .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmCreateGlobalVectorFromField()
471: @*/
472: PetscErrorCode DMSwarmDestroyGlobalVectorFromField(DM dm,const char fieldname[],Vec *vec)
473: {
477: DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec);
478: return(0);
479: }
481: /*@C
482: DMSwarmCreateLocalVectorFromField - Creates a Vec object sharing the array associated with a given field
484: Collective on dm
486: Input parameters:
487: + dm - a DMSwarm
488: - fieldname - the textual name given to a registered field
490: Output parameter:
491: . vec - the vector
493: Level: beginner
495: Notes:
496: The vector must be returned using a matching call to DMSwarmDestroyLocalVectorFromField().
498: .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmDestroyLocalVectorFromField()
499: @*/
500: PetscErrorCode DMSwarmCreateLocalVectorFromField(DM dm,const char fieldname[],Vec *vec)
501: {
502: MPI_Comm comm = PETSC_COMM_SELF;
506: DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec);
507: return(0);
508: }
510: /*@C
511: DMSwarmDestroyLocalVectorFromField - Destroys the Vec object which share the array associated with a given field
513: Collective on dm
515: Input parameters:
516: + dm - a DMSwarm
517: - fieldname - the textual name given to a registered field
519: Output parameter:
520: . vec - the vector
522: Level: beginner
524: .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmCreateLocalVectorFromField()
525: @*/
526: PetscErrorCode DMSwarmDestroyLocalVectorFromField(DM dm,const char fieldname[],Vec *vec)
527: {
531: DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec);
532: return(0);
533: }
535: /*
536: PetscErrorCode DMSwarmCreateGlobalVectorFromFields(DM dm,const PetscInt nf,const char *fieldnames[],Vec *vec)
537: {
538: return(0);
539: }
541: PetscErrorCode DMSwarmRestoreGlobalVectorFromFields(DM dm,Vec *vec)
542: {
543: return(0);
544: }
545: */
547: /*@C
548: DMSwarmInitializeFieldRegister - Initiates the registration of fields to a DMSwarm
550: Collective on dm
552: Input parameter:
553: . dm - a DMSwarm
555: Level: beginner
557: Notes:
558: After all fields have been registered, you must call DMSwarmFinalizeFieldRegister().
560: .seealso: DMSwarmFinalizeFieldRegister(), DMSwarmRegisterPetscDatatypeField(),
561: DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField()
562: @*/
563: PetscErrorCode DMSwarmInitializeFieldRegister(DM dm)
564: {
565: DM_Swarm *swarm = (DM_Swarm *) dm->data;
569: if (!swarm->field_registration_initialized) {
570: swarm->field_registration_initialized = PETSC_TRUE;
571: DMSwarmRegisterPetscDatatypeField(dm,DMSwarmField_pid,1,PETSC_INT64); /* unique identifer */
572: DMSwarmRegisterPetscDatatypeField(dm,DMSwarmField_rank,1,PETSC_INT); /* used for communication */
573: }
574: return(0);
575: }
577: /*@
578: DMSwarmFinalizeFieldRegister - Finalizes the registration of fields to a DMSwarm
580: Collective on dm
582: Input parameter:
583: . dm - a DMSwarm
585: Level: beginner
587: Notes:
588: After DMSwarmFinalizeFieldRegister() has been called, no new fields can be defined on the DMSwarm.
590: .seealso: DMSwarmInitializeFieldRegister(), DMSwarmRegisterPetscDatatypeField(),
591: DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField()
592: @*/
593: PetscErrorCode DMSwarmFinalizeFieldRegister(DM dm)
594: {
595: DM_Swarm *swarm = (DM_Swarm*)dm->data;
599: if (!swarm->field_registration_finalized) {
600: DMSwarmDataBucketFinalize(swarm->db);
601: }
602: swarm->field_registration_finalized = PETSC_TRUE;
603: return(0);
604: }
606: /*@
607: DMSwarmSetLocalSizes - Sets the length of all registered fields on the DMSwarm
609: Not collective
611: Input parameters:
612: + dm - a DMSwarm
613: . nlocal - the length of each registered field
614: - buffer - the length of the buffer used to efficient dynamic re-sizing
616: Level: beginner
618: .seealso: DMSwarmGetLocalSize()
619: @*/
620: PetscErrorCode DMSwarmSetLocalSizes(DM dm,PetscInt nlocal,PetscInt buffer)
621: {
622: DM_Swarm *swarm = (DM_Swarm*)dm->data;
626: PetscLogEventBegin(DMSWARM_SetSizes,0,0,0,0);
627: DMSwarmDataBucketSetSizes(swarm->db,nlocal,buffer);
628: PetscLogEventEnd(DMSWARM_SetSizes,0,0,0,0);
629: return(0);
630: }
632: /*@
633: DMSwarmSetCellDM - Attachs a DM to a DMSwarm
635: Collective on dm
637: Input parameters:
638: + dm - a DMSwarm
639: - dmcell - the DM to attach to the DMSwarm
641: Level: beginner
643: Notes:
644: The attached DM (dmcell) will be queried for point location and
645: neighbor MPI-rank information if DMSwarmMigrate() is called.
647: .seealso: DMSwarmSetType(), DMSwarmGetCellDM(), DMSwarmMigrate()
648: @*/
649: PetscErrorCode DMSwarmSetCellDM(DM dm,DM dmcell)
650: {
651: DM_Swarm *swarm = (DM_Swarm*)dm->data;
654: swarm->dmcell = dmcell;
655: return(0);
656: }
658: /*@
659: DMSwarmGetCellDM - Fetches the attached cell DM
661: Collective on dm
663: Input parameter:
664: . dm - a DMSwarm
666: Output parameter:
667: . dmcell - the DM which was attached to the DMSwarm
669: Level: beginner
671: .seealso: DMSwarmSetCellDM()
672: @*/
673: PetscErrorCode DMSwarmGetCellDM(DM dm,DM *dmcell)
674: {
675: DM_Swarm *swarm = (DM_Swarm*)dm->data;
678: *dmcell = swarm->dmcell;
679: return(0);
680: }
682: /*@
683: DMSwarmGetLocalSize - Retrives the local length of fields registered
685: Not collective
687: Input parameter:
688: . dm - a DMSwarm
690: Output parameter:
691: . nlocal - the length of each registered field
693: Level: beginner
695: .seealso: DMSwarmGetSize(), DMSwarmSetLocalSizes()
696: @*/
697: PetscErrorCode DMSwarmGetLocalSize(DM dm,PetscInt *nlocal)
698: {
699: DM_Swarm *swarm = (DM_Swarm*)dm->data;
703: if (nlocal) {DMSwarmDataBucketGetSizes(swarm->db,nlocal,NULL,NULL);}
704: return(0);
705: }
707: /*@
708: DMSwarmGetSize - Retrives the total length of fields registered
710: Collective on dm
712: Input parameter:
713: . dm - a DMSwarm
715: Output parameter:
716: . n - the total length of each registered field
718: Level: beginner
720: Note:
721: This calls MPI_Allreduce upon each call (inefficient but safe)
723: .seealso: DMSwarmGetLocalSize(), DMSwarmSetLocalSizes()
724: @*/
725: PetscErrorCode DMSwarmGetSize(DM dm,PetscInt *n)
726: {
727: DM_Swarm *swarm = (DM_Swarm*)dm->data;
729: PetscInt nlocal,ng;
732: DMSwarmDataBucketGetSizes(swarm->db,&nlocal,NULL,NULL);
733: MPI_Allreduce(&nlocal,&ng,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));
734: if (n) { *n = ng; }
735: return(0);
736: }
738: /*@C
739: DMSwarmRegisterPetscDatatypeField - Register a field to a DMSwarm with a native PETSc data type
741: Collective on dm
743: Input parameters:
744: + dm - a DMSwarm
745: . fieldname - the textual name to identify this field
746: . blocksize - the number of each data type
747: - type - a valid PETSc data type (PETSC_CHAR, PETSC_SHORT, PETSC_INT, PETSC_FLOAT, PETSC_REAL, PETSC_LONG)
749: Level: beginner
751: Notes:
752: The textual name for each registered field must be unique.
754: .seealso: DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField()
755: @*/
756: PetscErrorCode DMSwarmRegisterPetscDatatypeField(DM dm,const char fieldname[],PetscInt blocksize,PetscDataType type)
757: {
759: DM_Swarm *swarm = (DM_Swarm*)dm->data;
760: size_t size;
763: if (!swarm->field_registration_initialized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmInitializeFieldRegister() first");
764: if (swarm->field_registration_finalized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Cannot register additional fields after calling DMSwarmFinalizeFieldRegister() first");
766: if (type == PETSC_OBJECT) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
767: if (type == PETSC_FUNCTION) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
768: if (type == PETSC_STRING) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
769: if (type == PETSC_STRUCT) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
770: if (type == PETSC_DATATYPE_UNKNOWN) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
772: PetscDataTypeGetSize(type, &size);
773: /* Load a specific data type into data bucket, specifying textual name and its size in bytes */
774: DMSwarmDataBucketRegisterField(swarm->db,"DMSwarmRegisterPetscDatatypeField",fieldname,blocksize*size,NULL);
775: {
776: DMSwarmDataField gfield;
778: DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield);
779: DMSwarmDataFieldSetBlockSize(gfield,blocksize);
780: }
781: swarm->db->field[swarm->db->nfields-1]->petsc_type = type;
782: return(0);
783: }
785: /*@C
786: DMSwarmRegisterUserStructField - Register a user defined struct to a DMSwarm
788: Collective on dm
790: Input parameters:
791: + dm - a DMSwarm
792: . fieldname - the textual name to identify this field
793: - size - the size in bytes of the user struct of each data type
795: Level: beginner
797: Notes:
798: The textual name for each registered field must be unique.
800: .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmRegisterUserDatatypeField()
801: @*/
802: PetscErrorCode DMSwarmRegisterUserStructField(DM dm,const char fieldname[],size_t size)
803: {
805: DM_Swarm *swarm = (DM_Swarm*)dm->data;
808: DMSwarmDataBucketRegisterField(swarm->db,"DMSwarmRegisterUserStructField",fieldname,size,NULL);
809: swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_STRUCT ;
810: return(0);
811: }
813: /*@C
814: DMSwarmRegisterUserDatatypeField - Register a user defined data type to a DMSwarm
816: Collective on dm
818: Input parameters:
819: + dm - a DMSwarm
820: . fieldname - the textual name to identify this field
821: . size - the size in bytes of the user data type
822: - blocksize - the number of each data type
824: Level: beginner
826: Notes:
827: The textual name for each registered field must be unique.
829: .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField()
830: @*/
831: PetscErrorCode DMSwarmRegisterUserDatatypeField(DM dm,const char fieldname[],size_t size,PetscInt blocksize)
832: {
833: DM_Swarm *swarm = (DM_Swarm*)dm->data;
837: DMSwarmDataBucketRegisterField(swarm->db,"DMSwarmRegisterUserDatatypeField",fieldname,blocksize*size,NULL);
838: {
839: DMSwarmDataField gfield;
841: DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield);
842: DMSwarmDataFieldSetBlockSize(gfield,blocksize);
843: }
844: swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_DATATYPE_UNKNOWN;
845: return(0);
846: }
848: /*@C
849: DMSwarmGetField - Get access to the underlying array storing all entries associated with a registered field
851: Not collective
853: Input parameters:
854: + dm - a DMSwarm
855: - fieldname - the textual name to identify this field
857: Output parameters:
858: + blocksize - the number of each data type
859: . type - the data type
860: - data - pointer to raw array
862: Level: beginner
864: Notes:
865: The array must be returned using a matching call to DMSwarmRestoreField().
867: .seealso: DMSwarmRestoreField()
868: @*/
869: PetscErrorCode DMSwarmGetField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data)
870: {
871: DM_Swarm *swarm = (DM_Swarm*)dm->data;
872: DMSwarmDataField gfield;
876: if (!swarm->issetup) { DMSetUp(dm); }
877: DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield);
878: DMSwarmDataFieldGetAccess(gfield);
879: DMSwarmDataFieldGetEntries(gfield,data);
880: if (blocksize) {*blocksize = gfield->bs; }
881: if (type) { *type = gfield->petsc_type; }
882: return(0);
883: }
885: /*@C
886: DMSwarmRestoreField - Restore access to the underlying array storing all entries associated with a registered field
888: Not collective
890: Input parameters:
891: + dm - a DMSwarm
892: - fieldname - the textual name to identify this field
894: Output parameters:
895: + blocksize - the number of each data type
896: . type - the data type
897: - data - pointer to raw array
899: Level: beginner
901: Notes:
902: The user must call DMSwarmGetField() prior to calling DMSwarmRestoreField().
904: .seealso: DMSwarmGetField()
905: @*/
906: PetscErrorCode DMSwarmRestoreField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data)
907: {
908: DM_Swarm *swarm = (DM_Swarm*)dm->data;
909: DMSwarmDataField gfield;
913: DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield);
914: DMSwarmDataFieldRestoreAccess(gfield);
915: if (data) *data = NULL;
916: return(0);
917: }
919: /*@
920: DMSwarmAddPoint - Add space for one new point in the DMSwarm
922: Not collective
924: Input parameter:
925: . dm - a DMSwarm
927: Level: beginner
929: Notes:
930: The new point will have all fields initialized to zero.
932: .seealso: DMSwarmAddNPoints()
933: @*/
934: PetscErrorCode DMSwarmAddPoint(DM dm)
935: {
936: DM_Swarm *swarm = (DM_Swarm*)dm->data;
940: if (!swarm->issetup) {DMSetUp(dm);}
941: PetscLogEventBegin(DMSWARM_AddPoints,0,0,0,0);
942: DMSwarmDataBucketAddPoint(swarm->db);
943: PetscLogEventEnd(DMSWARM_AddPoints,0,0,0,0);
944: return(0);
945: }
947: /*@
948: DMSwarmAddNPoints - Add space for a number of new points in the DMSwarm
950: Not collective
952: Input parameters:
953: + dm - a DMSwarm
954: - npoints - the number of new points to add
956: Level: beginner
958: Notes:
959: The new point will have all fields initialized to zero.
961: .seealso: DMSwarmAddPoint()
962: @*/
963: PetscErrorCode DMSwarmAddNPoints(DM dm,PetscInt npoints)
964: {
965: DM_Swarm *swarm = (DM_Swarm*)dm->data;
967: PetscInt nlocal;
970: PetscLogEventBegin(DMSWARM_AddPoints,0,0,0,0);
971: DMSwarmDataBucketGetSizes(swarm->db,&nlocal,NULL,NULL);
972: nlocal = nlocal + npoints;
973: DMSwarmDataBucketSetSizes(swarm->db,nlocal,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT);
974: PetscLogEventEnd(DMSWARM_AddPoints,0,0,0,0);
975: return(0);
976: }
978: /*@
979: DMSwarmRemovePoint - Remove the last point from the DMSwarm
981: Not collective
983: Input parameter:
984: . dm - a DMSwarm
986: Level: beginner
988: .seealso: DMSwarmRemovePointAtIndex()
989: @*/
990: PetscErrorCode DMSwarmRemovePoint(DM dm)
991: {
992: DM_Swarm *swarm = (DM_Swarm*)dm->data;
996: PetscLogEventBegin(DMSWARM_RemovePoints,0,0,0,0);
997: DMSwarmDataBucketRemovePoint(swarm->db);
998: PetscLogEventEnd(DMSWARM_RemovePoints,0,0,0,0);
999: return(0);
1000: }
1002: /*@
1003: DMSwarmRemovePointAtIndex - Removes a specific point from the DMSwarm
1005: Not collective
1007: Input parameters:
1008: + dm - a DMSwarm
1009: - idx - index of point to remove
1011: Level: beginner
1013: .seealso: DMSwarmRemovePoint()
1014: @*/
1015: PetscErrorCode DMSwarmRemovePointAtIndex(DM dm,PetscInt idx)
1016: {
1017: DM_Swarm *swarm = (DM_Swarm*)dm->data;
1021: PetscLogEventBegin(DMSWARM_RemovePoints,0,0,0,0);
1022: DMSwarmDataBucketRemovePointAtIndex(swarm->db,idx);
1023: PetscLogEventEnd(DMSWARM_RemovePoints,0,0,0,0);
1024: return(0);
1025: }
1027: /*@
1028: DMSwarmCopyPoint - Copy point pj to point pi in the DMSwarm
1030: Not collective
1032: Input parameters:
1033: + dm - a DMSwarm
1034: . pi - the index of the point to copy
1035: - pj - the point index where the copy should be located
1037: Level: beginner
1039: .seealso: DMSwarmRemovePoint()
1040: @*/
1041: PetscErrorCode DMSwarmCopyPoint(DM dm,PetscInt pi,PetscInt pj)
1042: {
1043: DM_Swarm *swarm = (DM_Swarm*)dm->data;
1047: if (!swarm->issetup) {DMSetUp(dm);}
1048: DMSwarmDataBucketCopyPoint(swarm->db,pi,swarm->db,pj);
1049: return(0);
1050: }
1052: PetscErrorCode DMSwarmMigrate_Basic(DM dm,PetscBool remove_sent_points)
1053: {
1057: DMSwarmMigrate_Push_Basic(dm,remove_sent_points);
1058: return(0);
1059: }
1061: /*@
1062: DMSwarmMigrate - Relocates points defined in the DMSwarm to other MPI-ranks
1064: Collective on dm
1066: Input parameters:
1067: + dm - the DMSwarm
1068: - remove_sent_points - flag indicating if sent points should be removed from the current MPI-rank
1070: Notes:
1071: The DM will be modified to accomodate received points.
1072: If remove_sent_points = PETSC_TRUE, any points that were sent will be removed from the DM.
1073: Different styles of migration are supported. See DMSwarmSetMigrateType().
1075: Level: advanced
1077: .seealso: DMSwarmSetMigrateType()
1078: @*/
1079: PetscErrorCode DMSwarmMigrate(DM dm,PetscBool remove_sent_points)
1080: {
1081: DM_Swarm *swarm = (DM_Swarm*)dm->data;
1085: PetscLogEventBegin(DMSWARM_Migrate,0,0,0,0);
1086: switch (swarm->migrate_type) {
1087: case DMSWARM_MIGRATE_BASIC:
1088: DMSwarmMigrate_Basic(dm,remove_sent_points);
1089: break;
1090: case DMSWARM_MIGRATE_DMCELLNSCATTER:
1091: DMSwarmMigrate_CellDMScatter(dm,remove_sent_points);
1092: break;
1093: case DMSWARM_MIGRATE_DMCELLEXACT:
1094: SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_MIGRATE_DMCELLEXACT not implemented");
1095: break;
1096: case DMSWARM_MIGRATE_USER:
1097: SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_MIGRATE_USER not implemented");
1098: break;
1099: default:
1100: SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_MIGRATE type unknown");
1101: break;
1102: }
1103: PetscLogEventEnd(DMSWARM_Migrate,0,0,0,0);
1104: return(0);
1105: }
1107: PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm,PetscInt *globalsize);
1109: /*
1110: DMSwarmCollectViewCreate
1112: * Applies a collection method and gathers point neighbour points into dm
1114: Notes:
1115: Users should call DMSwarmCollectViewDestroy() after
1116: they have finished computations associated with the collected points
1117: */
1119: /*@
1120: DMSwarmCollectViewCreate - Applies a collection method and gathers points
1121: in neighbour MPI-ranks into the DMSwarm
1123: Collective on dm
1125: Input parameter:
1126: . dm - the DMSwarm
1128: Notes:
1129: Users should call DMSwarmCollectViewDestroy() after
1130: they have finished computations associated with the collected points
1131: Different collect methods are supported. See DMSwarmSetCollectType().
1133: Level: advanced
1135: .seealso: DMSwarmCollectViewDestroy(), DMSwarmSetCollectType()
1136: @*/
1137: PetscErrorCode DMSwarmCollectViewCreate(DM dm)
1138: {
1140: DM_Swarm *swarm = (DM_Swarm*)dm->data;
1141: PetscInt ng;
1144: if (swarm->collect_view_active) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"CollectView currently active");
1145: DMSwarmGetLocalSize(dm,&ng);
1146: switch (swarm->collect_type) {
1148: case DMSWARM_COLLECT_BASIC:
1149: DMSwarmMigrate_GlobalToLocal_Basic(dm,&ng);
1150: break;
1151: case DMSWARM_COLLECT_DMDABOUNDINGBOX:
1152: SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_COLLECT_DMDABOUNDINGBOX not implemented");
1153: break;
1154: case DMSWARM_COLLECT_GENERAL:
1155: SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_COLLECT_GENERAL not implemented");
1156: break;
1157: default:
1158: SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_COLLECT type unknown");
1159: break;
1160: }
1161: swarm->collect_view_active = PETSC_TRUE;
1162: swarm->collect_view_reset_nlocal = ng;
1163: return(0);
1164: }
1166: /*@
1167: DMSwarmCollectViewDestroy - Resets the DMSwarm to the size prior to calling DMSwarmCollectViewCreate()
1169: Collective on dm
1171: Input parameters:
1172: . dm - the DMSwarm
1174: Notes:
1175: Users should call DMSwarmCollectViewCreate() before this function is called.
1177: Level: advanced
1179: .seealso: DMSwarmCollectViewCreate(), DMSwarmSetCollectType()
1180: @*/
1181: PetscErrorCode DMSwarmCollectViewDestroy(DM dm)
1182: {
1184: DM_Swarm *swarm = (DM_Swarm*)dm->data;
1187: if (!swarm->collect_view_active) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"CollectView is currently not active");
1188: DMSwarmSetLocalSizes(dm,swarm->collect_view_reset_nlocal,-1);
1189: swarm->collect_view_active = PETSC_FALSE;
1190: return(0);
1191: }
1193: PetscErrorCode DMSwarmSetUpPIC(DM dm)
1194: {
1195: PetscInt dim;
1199: DMGetDimension(dm,&dim);
1200: if (dim < 1) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Dimension must be 1,2,3 - found %D",dim);
1201: if (dim > 3) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Dimension must be 1,2,3 - found %D",dim);
1202: DMSwarmRegisterPetscDatatypeField(dm,DMSwarmPICField_coor,dim,PETSC_DOUBLE);
1203: DMSwarmRegisterPetscDatatypeField(dm,DMSwarmPICField_cellid,1,PETSC_INT);
1204: return(0);
1205: }
1207: /*@
1208: DMSwarmSetType - Set particular flavor of DMSwarm
1210: Collective on dm
1212: Input parameters:
1213: + dm - the DMSwarm
1214: - stype - the DMSwarm type (e.g. DMSWARM_PIC)
1216: Level: advanced
1218: .seealso: DMSwarmSetMigrateType(), DMSwarmSetCollectType()
1219: @*/
1220: PetscErrorCode DMSwarmSetType(DM dm,DMSwarmType stype)
1221: {
1222: DM_Swarm *swarm = (DM_Swarm*)dm->data;
1226: swarm->swarm_type = stype;
1227: if (swarm->swarm_type == DMSWARM_PIC) {
1228: DMSwarmSetUpPIC(dm);
1229: }
1230: return(0);
1231: }
1233: PetscErrorCode DMSetup_Swarm(DM dm)
1234: {
1235: DM_Swarm *swarm = (DM_Swarm*)dm->data;
1237: PetscMPIInt rank;
1238: PetscInt p,npoints,*rankval;
1241: if (swarm->issetup) return(0);
1243: swarm->issetup = PETSC_TRUE;
1245: if (swarm->swarm_type == DMSWARM_PIC) {
1246: /* check dmcell exists */
1247: if (!swarm->dmcell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSWARM_PIC requires you call DMSwarmSetCellDM");
1249: if (swarm->dmcell->ops->locatepointssubdomain) {
1250: /* check methods exists for exact ownership identificiation */
1251: PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->ops->LocatePointsSubdomain\n");
1252: swarm->migrate_type = DMSWARM_MIGRATE_DMCELLEXACT;
1253: } else {
1254: /* check methods exist for point location AND rank neighbor identification */
1255: if (swarm->dmcell->ops->locatepoints) {
1256: PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->LocatePoints\n");
1257: } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSWARM_PIC requires the method CellDM->ops->locatepoints be defined");
1259: if (swarm->dmcell->ops->getneighbors) {
1260: PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->GetNeigbors\n");
1261: } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSWARM_PIC requires the method CellDM->ops->getneighbors be defined");
1263: swarm->migrate_type = DMSWARM_MIGRATE_DMCELLNSCATTER;
1264: }
1265: }
1267: DMSwarmFinalizeFieldRegister(dm);
1269: /* check some fields were registered */
1270: if (swarm->db->nfields <= 2) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"At least one field user must be registered via DMSwarmRegisterXXX()");
1272: /* check local sizes were set */
1273: if (swarm->db->L == -1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Local sizes must be set via DMSwarmSetLocalSizes()");
1275: /* initialize values in pid and rank placeholders */
1276: /* TODO: [pid - use MPI_Scan] */
1277: MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);
1278: DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);
1279: DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);
1280: for (p=0; p<npoints; p++) {
1281: rankval[p] = (PetscInt)rank;
1282: }
1283: DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);
1284: return(0);
1285: }
1287: extern PetscErrorCode DMSwarmSortDestroy(DMSwarmSort *_ctx);
1289: PetscErrorCode DMDestroy_Swarm(DM dm)
1290: {
1291: DM_Swarm *swarm = (DM_Swarm*)dm->data;
1295: DMSwarmDataBucketDestroy(&swarm->db);
1296: if (swarm->sort_context) {
1297: DMSwarmSortDestroy(&swarm->sort_context);
1298: }
1299: PetscFree(swarm);
1300: return(0);
1301: }
1303: PetscErrorCode DMSwarmView_Draw(DM dm, PetscViewer viewer)
1304: {
1305: DM cdm;
1306: PetscDraw draw;
1307: PetscReal *coords, oldPause;
1308: PetscInt Np, p, bs;
1312: PetscViewerDrawGetDraw(viewer, 0, &draw);
1313: DMSwarmGetCellDM(dm, &cdm);
1314: PetscDrawGetPause(draw, &oldPause);
1315: PetscDrawSetPause(draw, 0.0);
1316: DMView(cdm, viewer);
1317: PetscDrawSetPause(draw, oldPause);
1319: DMSwarmGetLocalSize(dm, &Np);
1320: DMSwarmGetField(dm, DMSwarmPICField_coor, &bs, NULL, (void **) &coords);
1321: for (p = 0; p < Np; ++p) {
1322: const PetscInt i = p*bs;
1324: PetscDrawEllipse(draw, coords[i], coords[i+1], 0.01, 0.01, PETSC_DRAW_BLUE);
1325: }
1326: DMSwarmRestoreField(dm, DMSwarmPICField_coor, &bs, NULL, (void **) &coords);
1327: PetscDrawFlush(draw);
1328: PetscDrawPause(draw);
1329: PetscDrawSave(draw);
1330: return(0);
1331: }
1333: PetscErrorCode DMView_Swarm(DM dm, PetscViewer viewer)
1334: {
1335: DM_Swarm *swarm = (DM_Swarm*)dm->data;
1336: PetscBool iascii,ibinary,ishdf5,isvtk,isdraw;
1342: PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii);
1343: PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY,&ibinary);
1344: PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk);
1345: PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5);
1346: PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw);
1347: if (iascii) {
1348: DMSwarmDataBucketView(PetscObjectComm((PetscObject)dm),swarm->db,NULL,DATABUCKET_VIEW_STDOUT);
1349: } else if (ibinary) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO Binary support");
1350: #if defined(PETSC_HAVE_HDF5)
1351: else if (ishdf5) {
1352: DMSwarmView_HDF5(dm, viewer);
1353: }
1354: #else
1355: else if (ishdf5) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"HDF5 not supported. Please reconfigure using --download-hdf5");
1356: #endif
1357: else if (isvtk) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO VTK support");
1358: else if (isdraw) {
1359: DMSwarmView_Draw(dm, viewer);
1360: }
1361: return(0);
1362: }
1364: /*MC
1366: DMSWARM = "swarm" - A DM object used to represent arrays of data (fields) of arbitrary data type.
1367: This implementation was designed for particle methods in which the underlying
1368: data required to be represented is both (i) dynamic in length, (ii) and of arbitrary data type.
1370: User data can be represented by DMSwarm through a registering "fields".
1371: To register a field, the user must provide;
1372: (a) a unique name;
1373: (b) the data type (or size in bytes);
1374: (c) the block size of the data.
1376: For example, suppose the Section 1.5 Writing Application Codes with PETSc requires a unique id, energy, momentum and density to be stored
1377: on a set of particles. Then the following code could be used
1379: $ DMSwarmInitializeFieldRegister(dm)
1380: $ DMSwarmRegisterPetscDatatypeField(dm,"uid",1,PETSC_LONG);
1381: $ DMSwarmRegisterPetscDatatypeField(dm,"energy",1,PETSC_REAL);
1382: $ DMSwarmRegisterPetscDatatypeField(dm,"momentum",3,PETSC_REAL);
1383: $ DMSwarmRegisterPetscDatatypeField(dm,"density",1,PETSC_FLOAT);
1384: $ DMSwarmFinalizeFieldRegister(dm)
1386: The fields represented by DMSwarm are dynamic and can be re-sized at any time.
1387: The only restriction imposed by DMSwarm is that all fields contain the same number of points.
1389: To support particle methods, "migration" techniques are provided. These methods migrate data
1390: between MPI-ranks.
1392: DMSwarm supports the methods DMCreateGlobalVector() and DMCreateLocalVector().
1393: As a DMSwarm may internally define and store values of different data types,
1394: before calling DMCreateGlobalVector() or DMCreateLocalVector(), the user must inform DMSwarm which
1395: fields should be used to define a Vec object via
1396: DMSwarmVectorDefineField()
1397: The specified field can be changed at any time - thereby permitting vectors
1398: compatible with different fields to be created.
1400: A dual representation of fields in the DMSwarm and a Vec object is permitted via
1401: DMSwarmCreateGlobalVectorFromField()
1402: Here the data defining the field in the DMSwarm is shared with a Vec.
1403: This is inherently unsafe if you alter the size of the field at any time between
1404: calls to DMSwarmCreateGlobalVectorFromField() and DMSwarmDestroyGlobalVectorFromField().
1405: If the local size of the DMSwarm does not match the local size of the global vector
1406: when DMSwarmDestroyGlobalVectorFromField() is called, an error is thrown.
1408: Additional high-level support is provided for Particle-In-Cell methods.
1409: Please refer to the man page for DMSwarmSetType().
1411: Level: beginner
1413: .seealso: DMType, DMCreate(), DMSetType()
1414: M*/
1415: PETSC_EXTERN PetscErrorCode DMCreate_Swarm(DM dm)
1416: {
1417: DM_Swarm *swarm;
1422: PetscNewLog(dm,&swarm);
1423: dm->data = swarm;
1425: DMSwarmDataBucketCreate(&swarm->db);
1426: DMSwarmInitializeFieldRegister(dm);
1428: swarm->vec_field_set = PETSC_FALSE;
1429: swarm->issetup = PETSC_FALSE;
1430: swarm->swarm_type = DMSWARM_BASIC;
1431: swarm->migrate_type = DMSWARM_MIGRATE_BASIC;
1432: swarm->collect_type = DMSWARM_COLLECT_BASIC;
1433: swarm->migrate_error_on_missing_point = PETSC_FALSE;
1435: swarm->dmcell = NULL;
1436: swarm->collect_view_active = PETSC_FALSE;
1437: swarm->collect_view_reset_nlocal = -1;
1439: dm->dim = 0;
1440: dm->ops->view = DMView_Swarm;
1441: dm->ops->load = NULL;
1442: dm->ops->setfromoptions = NULL;
1443: dm->ops->clone = NULL;
1444: dm->ops->setup = DMSetup_Swarm;
1445: dm->ops->createlocalsection = NULL;
1446: dm->ops->createdefaultconstraints = NULL;
1447: dm->ops->createglobalvector = DMCreateGlobalVector_Swarm;
1448: dm->ops->createlocalvector = DMCreateLocalVector_Swarm;
1449: dm->ops->getlocaltoglobalmapping = NULL;
1450: dm->ops->createfieldis = NULL;
1451: dm->ops->createcoordinatedm = NULL;
1452: dm->ops->getcoloring = NULL;
1453: dm->ops->creatematrix = NULL;
1454: dm->ops->createinterpolation = NULL;
1455: dm->ops->createinjection = NULL;
1456: dm->ops->createmassmatrix = DMCreateMassMatrix_Swarm;
1457: dm->ops->refine = NULL;
1458: dm->ops->coarsen = NULL;
1459: dm->ops->refinehierarchy = NULL;
1460: dm->ops->coarsenhierarchy = NULL;
1461: dm->ops->globaltolocalbegin = NULL;
1462: dm->ops->globaltolocalend = NULL;
1463: dm->ops->localtoglobalbegin = NULL;
1464: dm->ops->localtoglobalend = NULL;
1465: dm->ops->destroy = DMDestroy_Swarm;
1466: dm->ops->createsubdm = NULL;
1467: dm->ops->getdimpoints = NULL;
1468: dm->ops->locatepoints = NULL;
1469: return(0);
1470: }