xref: /petsc/src/dm/tutorials/swarm_ex3.c (revision c8025a5415d73fd1c6005393f2b0e60677bf5915)
1c4762a1bSJed Brown static char help[] = "Tests DMSwarm with DMShell\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscsf.h>
4c4762a1bSJed Brown #include <petscdm.h>
5c4762a1bSJed Brown #include <petscdmshell.h>
6c4762a1bSJed Brown #include <petscdmda.h>
7c4762a1bSJed Brown #include <petscdmswarm.h>
8c4762a1bSJed Brown #include <petsc/private/dmimpl.h>
9c4762a1bSJed Brown 
10d71ae5a4SJacob Faibussowitsch PetscErrorCode _DMLocatePoints_DMDARegular_IS(DM dm, Vec pos, IS *iscell)
11d71ae5a4SJacob Faibussowitsch {
12c4762a1bSJed Brown   PetscInt           p, n, bs, npoints, si, sj, milocal, mjlocal, mx, my;
13c4762a1bSJed Brown   DM                 dmregular;
14c4762a1bSJed Brown   PetscInt          *cellidx;
15c4762a1bSJed Brown   const PetscScalar *coor;
16c4762a1bSJed Brown   PetscReal          dx, dy;
17c4762a1bSJed Brown   PetscMPIInt        rank;
18c4762a1bSJed Brown 
19362febeeSStefano Zampini   PetscFunctionBegin;
209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
219566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(pos, &n));
229566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(pos, &bs));
23c4762a1bSJed Brown   npoints = n / bs;
24c4762a1bSJed Brown 
259566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(npoints, &cellidx));
269566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(dm, &dmregular));
279566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dmregular, &si, &sj, NULL, &milocal, &mjlocal, NULL));
289566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dmregular, NULL, &mx, &my, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
29c4762a1bSJed Brown 
30c4762a1bSJed Brown   dx = 2.0 / ((PetscReal)mx);
31c4762a1bSJed Brown   dy = 2.0 / ((PetscReal)my);
32c4762a1bSJed Brown 
339566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(pos, &coor));
34c4762a1bSJed Brown   for (p = 0; p < npoints; p++) {
35c4762a1bSJed Brown     PetscReal coorx, coory;
36c4762a1bSJed Brown     PetscInt  mi, mj;
37c4762a1bSJed Brown 
38c4762a1bSJed Brown     coorx = PetscRealPart(coor[2 * p]);
39c4762a1bSJed Brown     coory = PetscRealPart(coor[2 * p + 1]);
40c4762a1bSJed Brown 
41c4762a1bSJed Brown     mi = (PetscInt)((coorx - (-1.0)) / dx);
42c4762a1bSJed Brown     mj = (PetscInt)((coory - (-1.0)) / dy);
43c4762a1bSJed Brown 
44c4762a1bSJed Brown     cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;
45c4762a1bSJed Brown 
46c4762a1bSJed Brown     if ((mj >= sj) && (mj < sj + mjlocal)) {
47ad540459SPierre Jolivet       if ((mi >= si) && (mi < si + milocal)) cellidx[p] = (mi - si) + (mj - sj) * milocal;
48c4762a1bSJed Brown     }
49c4762a1bSJed Brown     if (coorx < -1.0) cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;
50c4762a1bSJed Brown     if (coorx > 1.0) cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;
51c4762a1bSJed Brown     if (coory < -1.0) cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;
52c4762a1bSJed Brown     if (coory > 1.0) cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;
53c4762a1bSJed Brown   }
549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(pos, &coor));
559566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, npoints, cellidx, PETSC_OWN_POINTER, iscell));
563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
57c4762a1bSJed Brown }
58c4762a1bSJed Brown 
59d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocatePoints_DMDARegular(DM dm, Vec pos, DMPointLocationType ltype, PetscSF cellSF)
60d71ae5a4SJacob Faibussowitsch {
61c4762a1bSJed Brown   IS              iscell;
62c4762a1bSJed Brown   PetscSFNode    *cells;
63c4762a1bSJed Brown   PetscInt        p, bs, npoints, nfound;
64c4762a1bSJed Brown   const PetscInt *boxCells;
65c4762a1bSJed Brown 
66362febeeSStefano Zampini   PetscFunctionBegin;
679566063dSJacob Faibussowitsch   PetscCall(_DMLocatePoints_DMDARegular_IS(dm, pos, &iscell));
689566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(pos, &npoints));
699566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(pos, &bs));
70c4762a1bSJed Brown   npoints = npoints / bs;
71c4762a1bSJed Brown 
729566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(npoints, &cells));
739566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(iscell, &boxCells));
74c4762a1bSJed Brown 
75c4762a1bSJed Brown   for (p = 0; p < npoints; p++) {
76c4762a1bSJed Brown     cells[p].rank  = 0;
77c4762a1bSJed Brown     cells[p].index = DMLOCATEPOINT_POINT_NOT_FOUND;
78c4762a1bSJed Brown     cells[p].index = boxCells[p];
79c4762a1bSJed Brown   }
809566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(iscell, &boxCells));
819566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iscell));
82c4762a1bSJed Brown   nfound = npoints;
839566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(cellSF, npoints, nfound, NULL, PETSC_OWN_POINTER, cells, PETSC_OWN_POINTER));
849566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iscell));
853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
86c4762a1bSJed Brown }
87c4762a1bSJed Brown 
88d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetNeighbors_DMDARegular(DM dm, PetscInt *nneighbors, const PetscMPIInt **neighbors)
89d71ae5a4SJacob Faibussowitsch {
90c4762a1bSJed Brown   DM dmregular;
91c4762a1bSJed Brown 
92362febeeSStefano Zampini   PetscFunctionBegin;
939566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(dm, &dmregular));
949566063dSJacob Faibussowitsch   PetscCall(DMGetNeighbors(dmregular, nneighbors, neighbors));
953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
96c4762a1bSJed Brown }
97c4762a1bSJed Brown 
98d71ae5a4SJacob Faibussowitsch PetscErrorCode SwarmViewGP(DM dms, const char prefix[])
99d71ae5a4SJacob Faibussowitsch {
100c4762a1bSJed Brown   PetscReal  *array;
101c4762a1bSJed Brown   PetscInt   *iarray;
102c4762a1bSJed Brown   PetscInt    npoints, p, bs;
103c4762a1bSJed Brown   FILE       *fp;
104c4762a1bSJed Brown   char        name[PETSC_MAX_PATH_LEN];
105c4762a1bSJed Brown   PetscMPIInt rank;
106c4762a1bSJed Brown 
107362febeeSStefano Zampini   PetscFunctionBegin;
1089566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1099566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "%s-rank%d.gp", prefix, rank));
110c4762a1bSJed Brown   fp = fopen(name, "w");
11128b400f6SJacob Faibussowitsch   PetscCheck(fp, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open file %s", name);
1129566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(dms, &npoints));
1139566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dms, DMSwarmPICField_coor, &bs, NULL, (void **)&array));
1149566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dms, "itag", NULL, NULL, (void **)&iarray));
115ad540459SPierre 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]);
1169566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dms, "itag", NULL, NULL, (void **)&iarray));
1179566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dms, DMSwarmPICField_coor, &bs, NULL, (void **)&array));
118c4762a1bSJed Brown   fclose(fp);
1193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
120c4762a1bSJed Brown }
121c4762a1bSJed Brown 
122c4762a1bSJed Brown /*
123c4762a1bSJed Brown  Create a DMShell and attach a regularly spaced DMDA for point location
124c4762a1bSJed Brown  Override methods for point location
125c4762a1bSJed Brown */
126d71ae5a4SJacob Faibussowitsch PetscErrorCode ex3_1(void)
127d71ae5a4SJacob Faibussowitsch {
128c4762a1bSJed Brown   DM          dms, dmcell, dmregular;
129c4762a1bSJed Brown   PetscMPIInt rank;
130c4762a1bSJed Brown   PetscInt    p, bs, nlocal, overlap, mx, tk;
131c4762a1bSJed Brown   PetscReal   dx;
132c4762a1bSJed Brown   PetscReal  *array, dt;
133c4762a1bSJed Brown   PetscInt   *iarray;
134c4762a1bSJed Brown   PetscRandom rand;
135c4762a1bSJed Brown 
136362febeeSStefano Zampini   PetscFunctionBegin;
1379566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
138c4762a1bSJed Brown 
139c4762a1bSJed Brown   /* Create a regularly spaced DMDA */
140c4762a1bSJed Brown   mx      = 40;
141c4762a1bSJed Brown   overlap = 0;
1429566063dSJacob 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));
1439566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dmregular));
1449566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dmregular));
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   dx = 2.0 / ((PetscReal)mx);
1479566063dSJacob 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));
148c4762a1bSJed Brown 
149c4762a1bSJed Brown   /* Create a DMShell for point location purposes */
1509566063dSJacob Faibussowitsch   PetscCall(DMShellCreate(PETSC_COMM_WORLD, &dmcell));
1519566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(dmcell, dmregular));
152c4762a1bSJed Brown   dmcell->ops->locatepoints = DMLocatePoints_DMDARegular;
153c4762a1bSJed Brown   dmcell->ops->getneighbors = DMGetNeighbors_DMDARegular;
154c4762a1bSJed Brown 
155c4762a1bSJed Brown   /* Create the swarm */
1569566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dms));
1579566063dSJacob Faibussowitsch   PetscCall(DMSetType(dms, DMSWARM));
1589566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(dms, 2));
1596a5217c0SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)dms, "Particles"));
160c4762a1bSJed Brown 
1619566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetType(dms, DMSWARM_PIC));
1629566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(dms, dmcell));
163c4762a1bSJed Brown 
1649566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "itag", 1, PETSC_INT));
1659566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(dms));
166c4762a1bSJed Brown   {
167c4762a1bSJed Brown     PetscInt           si, sj, milocal, mjlocal;
168c4762a1bSJed Brown     const PetscScalar *LA_coors;
169c4762a1bSJed Brown     Vec                coors;
170c4762a1bSJed Brown     PetscInt           cnt;
171c4762a1bSJed Brown 
1729566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(dmregular, &si, &sj, NULL, &milocal, &mjlocal, NULL));
1739566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(dmregular, &coors));
1749566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(coors, &LA_coors));
1759566063dSJacob Faibussowitsch     PetscCall(DMSwarmSetLocalSizes(dms, milocal * mjlocal, 4));
1769566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &nlocal));
1779566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, DMSwarmPICField_coor, &bs, NULL, (void **)&array));
178c4762a1bSJed Brown     cnt = 0;
1799566063dSJacob Faibussowitsch     PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rand));
1809566063dSJacob Faibussowitsch     PetscCall(PetscRandomSetInterval(rand, -1.0, 1.0));
181c4762a1bSJed Brown     for (p = 0; p < nlocal; p++) {
182c4762a1bSJed Brown       PetscReal px, py, rx, ry, r2;
183c4762a1bSJed Brown 
1849566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValueReal(rand, &rx));
1859566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValueReal(rand, &ry));
186c4762a1bSJed Brown 
187c4762a1bSJed Brown       px = PetscRealPart(LA_coors[2 * p + 0]) + 0.1 * rx * dx;
188c4762a1bSJed Brown       py = PetscRealPart(LA_coors[2 * p + 1]) + 0.1 * ry * dx;
189c4762a1bSJed Brown 
190c4762a1bSJed Brown       r2 = px * px + py * py;
191c4762a1bSJed Brown       if (r2 < 0.75 * 0.75) {
192c4762a1bSJed Brown         array[bs * cnt + 0] = px;
193c4762a1bSJed Brown         array[bs * cnt + 1] = py;
194c4762a1bSJed Brown         cnt++;
195c4762a1bSJed Brown       }
196c4762a1bSJed Brown     }
1979566063dSJacob Faibussowitsch     PetscCall(PetscRandomDestroy(&rand));
1989566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, DMSwarmPICField_coor, &bs, NULL, (void **)&array));
1999566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(coors, &LA_coors));
2009566063dSJacob Faibussowitsch     PetscCall(DMSwarmSetLocalSizes(dms, cnt, 4));
201c4762a1bSJed Brown 
2029566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &nlocal));
2039566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "itag", &bs, NULL, (void **)&iarray));
204ad540459SPierre Jolivet     for (p = 0; p < nlocal; p++) iarray[p] = (PetscInt)rank;
2059566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "itag", &bs, NULL, (void **)&iarray));
206c4762a1bSJed Brown   }
207c4762a1bSJed Brown 
2089566063dSJacob Faibussowitsch   PetscCall(DMView(dms, PETSC_VIEWER_STDOUT_WORLD));
2099566063dSJacob Faibussowitsch   PetscCall(SwarmViewGP(dms, "step0"));
210c4762a1bSJed Brown 
211c4762a1bSJed Brown   dt = 0.1;
212c4762a1bSJed Brown   for (tk = 1; tk < 20; tk++) {
213c4762a1bSJed Brown     char prefix[PETSC_MAX_PATH_LEN];
21463a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Step %" PetscInt_FMT " \n", tk));
215c4762a1bSJed Brown     /* push points */
2169566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &nlocal));
2179566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, DMSwarmPICField_coor, &bs, NULL, (void **)&array));
218c4762a1bSJed Brown     for (p = 0; p < nlocal; p++) {
219c4762a1bSJed Brown       PetscReal cx, cy, vx, vy;
220c4762a1bSJed Brown 
221c4762a1bSJed Brown       cx = array[2 * p];
222c4762a1bSJed Brown       cy = array[2 * p + 1];
223c4762a1bSJed Brown       vx = cy;
224c4762a1bSJed Brown       vy = -cx;
225c4762a1bSJed Brown 
226c4762a1bSJed Brown       array[2 * p] += dt * vx;
227c4762a1bSJed Brown       array[2 * p + 1] += dt * vy;
228c4762a1bSJed Brown     }
2299566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, DMSwarmPICField_coor, &bs, NULL, (void **)&array));
230c4762a1bSJed Brown 
231c4762a1bSJed Brown     /* migrate points */
2329566063dSJacob Faibussowitsch     PetscCall(DMSwarmMigrate(dms, PETSC_TRUE));
233c4762a1bSJed Brown     /* view points */
23463a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN - 1, "step%" PetscInt_FMT, tk));
235c4762a1bSJed Brown     /* should use the regular SwarmView() api, not one for a particular type */
2369566063dSJacob Faibussowitsch     PetscCall(SwarmViewGP(dms, prefix));
237c4762a1bSJed Brown   }
2389566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmregular));
2399566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmcell));
2409566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dms));
2413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
242c4762a1bSJed Brown }
243c4762a1bSJed Brown 
244d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
245d71ae5a4SJacob Faibussowitsch {
246327415f7SBarry Smith   PetscFunctionBeginUser;
247*c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
2489566063dSJacob Faibussowitsch   PetscCall(ex3_1());
2499566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
250b122ec5aSJacob Faibussowitsch   return 0;
251c4762a1bSJed Brown }
252c4762a1bSJed Brown 
253c4762a1bSJed Brown /*TEST
254c4762a1bSJed Brown 
255c4762a1bSJed Brown    build:
256c4762a1bSJed Brown       requires: double !complex
257c4762a1bSJed Brown 
258c4762a1bSJed Brown    test:
259c4762a1bSJed Brown       filter: grep -v atomic
260c4762a1bSJed Brown       filter_output: grep -v atomic
261c4762a1bSJed Brown TEST*/
262