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