1*19307e5cSMatthew G. Knepley #include "petscsys.h"
2c4762a1bSJed Brown static char help[] = "Tests DMSwarm with DMShell\n\n";
3c4762a1bSJed Brown
4c4762a1bSJed Brown #include <petscsf.h>
5c4762a1bSJed Brown #include <petscdm.h>
6c4762a1bSJed Brown #include <petscdmshell.h>
7c4762a1bSJed Brown #include <petscdmda.h>
8c4762a1bSJed Brown #include <petscdmswarm.h>
9c4762a1bSJed Brown #include <petsc/private/dmimpl.h>
10c4762a1bSJed Brown
_DMLocatePoints_DMDARegular_IS(DM dm,Vec pos,IS * iscell)11d71ae5a4SJacob Faibussowitsch PetscErrorCode _DMLocatePoints_DMDARegular_IS(DM dm, Vec pos, IS *iscell)
12d71ae5a4SJacob Faibussowitsch {
13c4762a1bSJed Brown PetscInt p, n, bs, npoints, si, sj, milocal, mjlocal, mx, my;
14c4762a1bSJed Brown DM dmregular;
15c4762a1bSJed Brown PetscInt *cellidx;
16c4762a1bSJed Brown const PetscScalar *coor;
17c4762a1bSJed Brown PetscReal dx, dy;
18c4762a1bSJed Brown PetscMPIInt rank;
19c4762a1bSJed Brown
20362febeeSStefano Zampini PetscFunctionBegin;
219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
229566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(pos, &n));
239566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(pos, &bs));
24c4762a1bSJed Brown npoints = n / bs;
25c4762a1bSJed Brown
269566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npoints, &cellidx));
279566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(dm, &dmregular));
289566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dmregular, &si, &sj, NULL, &milocal, &mjlocal, NULL));
299566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dmregular, NULL, &mx, &my, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
30c4762a1bSJed Brown
31c4762a1bSJed Brown dx = 2.0 / ((PetscReal)mx);
32c4762a1bSJed Brown dy = 2.0 / ((PetscReal)my);
33c4762a1bSJed Brown
349566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pos, &coor));
35c4762a1bSJed Brown for (p = 0; p < npoints; p++) {
36c4762a1bSJed Brown PetscReal coorx, coory;
37c4762a1bSJed Brown PetscInt mi, mj;
38c4762a1bSJed Brown
39c4762a1bSJed Brown coorx = PetscRealPart(coor[2 * p]);
40c4762a1bSJed Brown coory = PetscRealPart(coor[2 * p + 1]);
41c4762a1bSJed Brown
42c4762a1bSJed Brown mi = (PetscInt)((coorx - (-1.0)) / dx);
43c4762a1bSJed Brown mj = (PetscInt)((coory - (-1.0)) / dy);
44c4762a1bSJed Brown
45c4762a1bSJed Brown cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;
46c4762a1bSJed Brown
47c4762a1bSJed Brown if ((mj >= sj) && (mj < sj + mjlocal)) {
48ad540459SPierre Jolivet if ((mi >= si) && (mi < si + milocal)) cellidx[p] = (mi - si) + (mj - sj) * milocal;
49c4762a1bSJed Brown }
50c4762a1bSJed Brown if (coorx < -1.0) cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;
51c4762a1bSJed Brown if (coorx > 1.0) cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;
52c4762a1bSJed Brown if (coory < -1.0) cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;
53c4762a1bSJed Brown if (coory > 1.0) cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;
54c4762a1bSJed Brown }
559566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pos, &coor));
569566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, npoints, cellidx, PETSC_OWN_POINTER, iscell));
573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
58c4762a1bSJed Brown }
59c4762a1bSJed Brown
DMLocatePoints_DMDARegular(DM dm,Vec pos,DMPointLocationType ltype,PetscSF cellSF)60d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocatePoints_DMDARegular(DM dm, Vec pos, DMPointLocationType ltype, PetscSF cellSF)
61d71ae5a4SJacob Faibussowitsch {
62c4762a1bSJed Brown IS iscell;
63c4762a1bSJed Brown PetscSFNode *cells;
64c4762a1bSJed Brown PetscInt p, bs, npoints, nfound;
65c4762a1bSJed Brown const PetscInt *boxCells;
66c4762a1bSJed Brown
67362febeeSStefano Zampini PetscFunctionBegin;
689566063dSJacob Faibussowitsch PetscCall(_DMLocatePoints_DMDARegular_IS(dm, pos, &iscell));
699566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(pos, &npoints));
709566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(pos, &bs));
71c4762a1bSJed Brown npoints = npoints / bs;
72c4762a1bSJed Brown
739566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npoints, &cells));
749566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iscell, &boxCells));
75c4762a1bSJed Brown
76c4762a1bSJed Brown for (p = 0; p < npoints; p++) {
77c4762a1bSJed Brown cells[p].rank = 0;
78c4762a1bSJed Brown cells[p].index = DMLOCATEPOINT_POINT_NOT_FOUND;
79c4762a1bSJed Brown cells[p].index = boxCells[p];
80c4762a1bSJed Brown }
819566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscell, &boxCells));
829566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscell));
83c4762a1bSJed Brown nfound = npoints;
849566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(cellSF, npoints, nfound, NULL, PETSC_OWN_POINTER, cells, PETSC_OWN_POINTER));
859566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscell));
863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
87c4762a1bSJed Brown }
88c4762a1bSJed Brown
DMGetNeighbors_DMDARegular(DM dm,PetscInt * nneighbors,const PetscMPIInt ** neighbors)89d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetNeighbors_DMDARegular(DM dm, PetscInt *nneighbors, const PetscMPIInt **neighbors)
90d71ae5a4SJacob Faibussowitsch {
91c4762a1bSJed Brown DM dmregular;
92c4762a1bSJed Brown
93362febeeSStefano Zampini PetscFunctionBegin;
949566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(dm, &dmregular));
959566063dSJacob Faibussowitsch PetscCall(DMGetNeighbors(dmregular, nneighbors, neighbors));
963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
97c4762a1bSJed Brown }
98c4762a1bSJed Brown
SwarmViewGP(DM dms,const char prefix[])99d71ae5a4SJacob Faibussowitsch PetscErrorCode SwarmViewGP(DM dms, const char prefix[])
100d71ae5a4SJacob Faibussowitsch {
101c4762a1bSJed Brown PetscReal *array;
102c4762a1bSJed Brown PetscInt *iarray;
103c4762a1bSJed Brown PetscInt npoints, p, bs;
104c4762a1bSJed Brown FILE *fp;
105c4762a1bSJed Brown char name[PETSC_MAX_PATH_LEN];
106c4762a1bSJed Brown PetscMPIInt rank;
107c4762a1bSJed Brown
108362febeeSStefano Zampini PetscFunctionBegin;
1099566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1109566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "%s-rank%d.gp", prefix, rank));
111c4762a1bSJed Brown fp = fopen(name, "w");
11228b400f6SJacob Faibussowitsch PetscCheck(fp, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open file %s", name);
1139566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dms, &npoints));
1149566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dms, DMSwarmPICField_coor, &bs, NULL, (void **)&array));
1159566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dms, "itag", NULL, NULL, (void **)&iarray));
116ad540459SPierre Jolivet for (p = 0; p < npoints; p++) fprintf(fp, "%+1.4e %+1.4e %1.4e\n", array[2 * p], array[2 * p + 1], (double)iarray[p]);
1179566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dms, "itag", NULL, NULL, (void **)&iarray));
1189566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dms, DMSwarmPICField_coor, &bs, NULL, (void **)&array));
119c4762a1bSJed Brown fclose(fp);
1203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
121c4762a1bSJed Brown }
122c4762a1bSJed Brown
123c4762a1bSJed Brown /*
124c4762a1bSJed Brown Create a DMShell and attach a regularly spaced DMDA for point location
125c4762a1bSJed Brown Override methods for point location
126c4762a1bSJed Brown */
ex3_1(void)127d71ae5a4SJacob Faibussowitsch PetscErrorCode ex3_1(void)
128d71ae5a4SJacob Faibussowitsch {
129c4762a1bSJed Brown DM dms, dmcell, dmregular;
130c4762a1bSJed Brown PetscMPIInt rank;
131c4762a1bSJed Brown PetscInt p, bs, nlocal, overlap, mx, tk;
132c4762a1bSJed Brown PetscReal dx;
133c4762a1bSJed Brown PetscReal *array, dt;
134c4762a1bSJed Brown PetscInt *iarray;
135c4762a1bSJed Brown PetscRandom rand;
136c4762a1bSJed Brown
137362febeeSStefano Zampini PetscFunctionBegin;
1389566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
139c4762a1bSJed Brown
140c4762a1bSJed Brown /* Create a regularly spaced DMDA */
141c4762a1bSJed Brown mx = 40;
142c4762a1bSJed Brown overlap = 0;
1439566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, mx, mx, PETSC_DECIDE, PETSC_DECIDE, 1, overlap, NULL, NULL, &dmregular));
1449566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dmregular));
1459566063dSJacob Faibussowitsch PetscCall(DMSetUp(dmregular));
146c4762a1bSJed Brown
147c4762a1bSJed Brown dx = 2.0 / ((PetscReal)mx);
1489566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(dmregular, -1.0 + 0.5 * dx, 1.0 - 0.5 * dx, -1.0 + 0.5 * dx, 1.0 - 0.5 * dx, -1.0, 1.0));
149c4762a1bSJed Brown
150c4762a1bSJed Brown /* Create a DMShell for point location purposes */
1519566063dSJacob Faibussowitsch PetscCall(DMShellCreate(PETSC_COMM_WORLD, &dmcell));
152*19307e5cSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)dmcell, "celldm"));
1539566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(dmcell, dmregular));
154c4762a1bSJed Brown dmcell->ops->locatepoints = DMLocatePoints_DMDARegular;
155c4762a1bSJed Brown dmcell->ops->getneighbors = DMGetNeighbors_DMDARegular;
156c4762a1bSJed Brown
157c4762a1bSJed Brown /* Create the swarm */
1589566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD, &dms));
1599566063dSJacob Faibussowitsch PetscCall(DMSetType(dms, DMSWARM));
1609566063dSJacob Faibussowitsch PetscCall(DMSetDimension(dms, 2));
1616a5217c0SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)dms, "Particles"));
162c4762a1bSJed Brown
1639566063dSJacob Faibussowitsch PetscCall(DMSwarmSetType(dms, DMSWARM_PIC));
1649566063dSJacob Faibussowitsch PetscCall(DMSwarmSetCellDM(dms, dmcell));
165c4762a1bSJed Brown
1669566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "itag", 1, PETSC_INT));
1679566063dSJacob Faibussowitsch PetscCall(DMSwarmFinalizeFieldRegister(dms));
168c4762a1bSJed Brown {
169c4762a1bSJed Brown PetscInt si, sj, milocal, mjlocal;
170c4762a1bSJed Brown const PetscScalar *LA_coors;
171c4762a1bSJed Brown Vec coors;
172c4762a1bSJed Brown PetscInt cnt;
173c4762a1bSJed Brown
1749566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dmregular, &si, &sj, NULL, &milocal, &mjlocal, NULL));
1759566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(dmregular, &coors));
1769566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coors, &LA_coors));
1779566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dms, milocal * mjlocal, 4));
1789566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dms, &nlocal));
1799566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dms, DMSwarmPICField_coor, &bs, NULL, (void **)&array));
180c4762a1bSJed Brown cnt = 0;
1819566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rand));
1829566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rand, -1.0, 1.0));
183c4762a1bSJed Brown for (p = 0; p < nlocal; p++) {
184c4762a1bSJed Brown PetscReal px, py, rx, ry, r2;
185c4762a1bSJed Brown
1869566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(rand, &rx));
1879566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(rand, &ry));
188c4762a1bSJed Brown
189c4762a1bSJed Brown px = PetscRealPart(LA_coors[2 * p + 0]) + 0.1 * rx * dx;
190c4762a1bSJed Brown py = PetscRealPart(LA_coors[2 * p + 1]) + 0.1 * ry * dx;
191c4762a1bSJed Brown
192c4762a1bSJed Brown r2 = px * px + py * py;
193c4762a1bSJed Brown if (r2 < 0.75 * 0.75) {
194c4762a1bSJed Brown array[bs * cnt + 0] = px;
195c4762a1bSJed Brown array[bs * cnt + 1] = py;
196c4762a1bSJed Brown cnt++;
197c4762a1bSJed Brown }
198c4762a1bSJed Brown }
1999566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand));
2009566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dms, DMSwarmPICField_coor, &bs, NULL, (void **)&array));
2019566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coors, &LA_coors));
2029566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dms, cnt, 4));
203c4762a1bSJed Brown
2049566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dms, &nlocal));
2059566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dms, "itag", &bs, NULL, (void **)&iarray));
206ad540459SPierre Jolivet for (p = 0; p < nlocal; p++) iarray[p] = (PetscInt)rank;
2079566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dms, "itag", &bs, NULL, (void **)&iarray));
208c4762a1bSJed Brown }
209c4762a1bSJed Brown
2109566063dSJacob Faibussowitsch PetscCall(DMView(dms, PETSC_VIEWER_STDOUT_WORLD));
2119566063dSJacob Faibussowitsch PetscCall(SwarmViewGP(dms, "step0"));
212c4762a1bSJed Brown
213c4762a1bSJed Brown dt = 0.1;
214c4762a1bSJed Brown for (tk = 1; tk < 20; tk++) {
215c4762a1bSJed Brown char prefix[PETSC_MAX_PATH_LEN];
21663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Step %" PetscInt_FMT " \n", tk));
217c4762a1bSJed Brown /* push points */
2189566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dms, &nlocal));
2199566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dms, DMSwarmPICField_coor, &bs, NULL, (void **)&array));
220c4762a1bSJed Brown for (p = 0; p < nlocal; p++) {
221c4762a1bSJed Brown PetscReal cx, cy, vx, vy;
222c4762a1bSJed Brown
223c4762a1bSJed Brown cx = array[2 * p];
224c4762a1bSJed Brown cy = array[2 * p + 1];
225c4762a1bSJed Brown vx = cy;
226c4762a1bSJed Brown vy = -cx;
227c4762a1bSJed Brown
228c4762a1bSJed Brown array[2 * p] += dt * vx;
229c4762a1bSJed Brown array[2 * p + 1] += dt * vy;
230c4762a1bSJed Brown }
2319566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dms, DMSwarmPICField_coor, &bs, NULL, (void **)&array));
232c4762a1bSJed Brown
233c4762a1bSJed Brown /* migrate points */
2349566063dSJacob Faibussowitsch PetscCall(DMSwarmMigrate(dms, PETSC_TRUE));
235c4762a1bSJed Brown /* view points */
23663a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN - 1, "step%" PetscInt_FMT, tk));
237c4762a1bSJed Brown /* should use the regular SwarmView() api, not one for a particular type */
2389566063dSJacob Faibussowitsch PetscCall(SwarmViewGP(dms, prefix));
239c4762a1bSJed Brown }
2409566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmregular));
2419566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmcell));
2429566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dms));
2433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
244c4762a1bSJed Brown }
245c4762a1bSJed Brown
main(int argc,char ** argv)246d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
247d71ae5a4SJacob Faibussowitsch {
248327415f7SBarry Smith PetscFunctionBeginUser;
249c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
2509566063dSJacob Faibussowitsch PetscCall(ex3_1());
2519566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
252b122ec5aSJacob Faibussowitsch return 0;
253c4762a1bSJed Brown }
254c4762a1bSJed Brown
255c4762a1bSJed Brown /*TEST
256c4762a1bSJed Brown
257c4762a1bSJed Brown build:
258c4762a1bSJed Brown requires: double !complex
259c4762a1bSJed Brown
260c4762a1bSJed Brown test:
261c4762a1bSJed Brown filter: grep -v atomic
262c4762a1bSJed Brown filter_output: grep -v atomic
263c4762a1bSJed Brown TEST*/
264