xref: /petsc/src/dm/tutorials/swarm_ex3.c (revision ad540459ab38c4a232710a68077e487eb6627fe2)
1c4762a1bSJed Brown 
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 
119371c9d4SSatish Balay PetscErrorCode _DMLocatePoints_DMDARegular_IS(DM dm, Vec pos, IS *iscell) {
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)) {
47*ad540459SPierre 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));
56c4762a1bSJed Brown   PetscFunctionReturn(0);
57c4762a1bSJed Brown }
58c4762a1bSJed Brown 
599371c9d4SSatish Balay PetscErrorCode DMLocatePoints_DMDARegular(DM dm, Vec pos, DMPointLocationType ltype, PetscSF cellSF) {
60c4762a1bSJed Brown   IS              iscell;
61c4762a1bSJed Brown   PetscSFNode    *cells;
62c4762a1bSJed Brown   PetscInt        p, bs, npoints, nfound;
63c4762a1bSJed Brown   const PetscInt *boxCells;
64c4762a1bSJed Brown 
65362febeeSStefano Zampini   PetscFunctionBegin;
669566063dSJacob Faibussowitsch   PetscCall(_DMLocatePoints_DMDARegular_IS(dm, pos, &iscell));
679566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(pos, &npoints));
689566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(pos, &bs));
69c4762a1bSJed Brown   npoints = npoints / bs;
70c4762a1bSJed Brown 
719566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(npoints, &cells));
729566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(iscell, &boxCells));
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   for (p = 0; p < npoints; p++) {
75c4762a1bSJed Brown     cells[p].rank  = 0;
76c4762a1bSJed Brown     cells[p].index = DMLOCATEPOINT_POINT_NOT_FOUND;
77c4762a1bSJed Brown     cells[p].index = boxCells[p];
78c4762a1bSJed Brown   }
799566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(iscell, &boxCells));
809566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iscell));
81c4762a1bSJed Brown   nfound = npoints;
829566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(cellSF, npoints, nfound, NULL, PETSC_OWN_POINTER, cells, PETSC_OWN_POINTER));
839566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iscell));
84c4762a1bSJed Brown   PetscFunctionReturn(0);
85c4762a1bSJed Brown }
86c4762a1bSJed Brown 
879371c9d4SSatish Balay PetscErrorCode DMGetNeighbors_DMDARegular(DM dm, PetscInt *nneighbors, const PetscMPIInt **neighbors) {
88c4762a1bSJed Brown   DM dmregular;
89c4762a1bSJed Brown 
90362febeeSStefano Zampini   PetscFunctionBegin;
919566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(dm, &dmregular));
929566063dSJacob Faibussowitsch   PetscCall(DMGetNeighbors(dmregular, nneighbors, neighbors));
93c4762a1bSJed Brown   PetscFunctionReturn(0);
94c4762a1bSJed Brown }
95c4762a1bSJed Brown 
969371c9d4SSatish Balay PetscErrorCode SwarmViewGP(DM dms, const char prefix[]) {
97c4762a1bSJed Brown   PetscReal  *array;
98c4762a1bSJed Brown   PetscInt   *iarray;
99c4762a1bSJed Brown   PetscInt    npoints, p, bs;
100c4762a1bSJed Brown   FILE       *fp;
101c4762a1bSJed Brown   char        name[PETSC_MAX_PATH_LEN];
102c4762a1bSJed Brown   PetscMPIInt rank;
103c4762a1bSJed Brown 
104362febeeSStefano Zampini   PetscFunctionBegin;
1059566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1069566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "%s-rank%d.gp", prefix, rank));
107c4762a1bSJed Brown   fp = fopen(name, "w");
10828b400f6SJacob Faibussowitsch   PetscCheck(fp, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open file %s", name);
1099566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(dms, &npoints));
1109566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dms, DMSwarmPICField_coor, &bs, NULL, (void **)&array));
1119566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dms, "itag", NULL, NULL, (void **)&iarray));
112*ad540459SPierre 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]);
1139566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dms, "itag", NULL, NULL, (void **)&iarray));
1149566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dms, DMSwarmPICField_coor, &bs, NULL, (void **)&array));
115c4762a1bSJed Brown   fclose(fp);
116c4762a1bSJed Brown   PetscFunctionReturn(0);
117c4762a1bSJed Brown }
118c4762a1bSJed Brown 
119c4762a1bSJed Brown /*
120c4762a1bSJed Brown  Create a DMShell and attach a regularly spaced DMDA for point location
121c4762a1bSJed Brown  Override methods for point location
122c4762a1bSJed Brown */
1239371c9d4SSatish Balay PetscErrorCode ex3_1(void) {
124c4762a1bSJed Brown   DM          dms, dmcell, dmregular;
125c4762a1bSJed Brown   PetscMPIInt rank;
126c4762a1bSJed Brown   PetscInt    p, bs, nlocal, overlap, mx, tk;
127c4762a1bSJed Brown   PetscReal   dx;
128c4762a1bSJed Brown   PetscReal  *array, dt;
129c4762a1bSJed Brown   PetscInt   *iarray;
130c4762a1bSJed Brown   PetscRandom rand;
131c4762a1bSJed Brown 
132362febeeSStefano Zampini   PetscFunctionBegin;
1339566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
134c4762a1bSJed Brown 
135c4762a1bSJed Brown   /* Create a regularly spaced DMDA */
136c4762a1bSJed Brown   mx      = 40;
137c4762a1bSJed Brown   overlap = 0;
1389566063dSJacob 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));
1399566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dmregular));
1409566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dmregular));
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   dx = 2.0 / ((PetscReal)mx);
1439566063dSJacob 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));
144c4762a1bSJed Brown 
145c4762a1bSJed Brown   /* Create a DMShell for point location purposes */
1469566063dSJacob Faibussowitsch   PetscCall(DMShellCreate(PETSC_COMM_WORLD, &dmcell));
1479566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(dmcell, dmregular));
148c4762a1bSJed Brown   dmcell->ops->locatepoints = DMLocatePoints_DMDARegular;
149c4762a1bSJed Brown   dmcell->ops->getneighbors = DMGetNeighbors_DMDARegular;
150c4762a1bSJed Brown 
151c4762a1bSJed Brown   /* Create the swarm */
1529566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dms));
1539566063dSJacob Faibussowitsch   PetscCall(DMSetType(dms, DMSWARM));
1549566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(dms, 2));
1556a5217c0SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)dms, "Particles"));
156c4762a1bSJed Brown 
1579566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetType(dms, DMSWARM_PIC));
1589566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(dms, dmcell));
159c4762a1bSJed Brown 
1609566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "itag", 1, PETSC_INT));
1619566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(dms));
162c4762a1bSJed Brown   {
163c4762a1bSJed Brown     PetscInt           si, sj, milocal, mjlocal;
164c4762a1bSJed Brown     const PetscScalar *LA_coors;
165c4762a1bSJed Brown     Vec                coors;
166c4762a1bSJed Brown     PetscInt           cnt;
167c4762a1bSJed Brown 
1689566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(dmregular, &si, &sj, NULL, &milocal, &mjlocal, NULL));
1699566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(dmregular, &coors));
1709566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(coors, &LA_coors));
1719566063dSJacob Faibussowitsch     PetscCall(DMSwarmSetLocalSizes(dms, milocal * mjlocal, 4));
1729566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &nlocal));
1739566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, DMSwarmPICField_coor, &bs, NULL, (void **)&array));
174c4762a1bSJed Brown     cnt = 0;
1759566063dSJacob Faibussowitsch     PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rand));
1769566063dSJacob Faibussowitsch     PetscCall(PetscRandomSetInterval(rand, -1.0, 1.0));
177c4762a1bSJed Brown     for (p = 0; p < nlocal; p++) {
178c4762a1bSJed Brown       PetscReal px, py, rx, ry, r2;
179c4762a1bSJed Brown 
1809566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValueReal(rand, &rx));
1819566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValueReal(rand, &ry));
182c4762a1bSJed Brown 
183c4762a1bSJed Brown       px = PetscRealPart(LA_coors[2 * p + 0]) + 0.1 * rx * dx;
184c4762a1bSJed Brown       py = PetscRealPart(LA_coors[2 * p + 1]) + 0.1 * ry * dx;
185c4762a1bSJed Brown 
186c4762a1bSJed Brown       r2 = px * px + py * py;
187c4762a1bSJed Brown       if (r2 < 0.75 * 0.75) {
188c4762a1bSJed Brown         array[bs * cnt + 0] = px;
189c4762a1bSJed Brown         array[bs * cnt + 1] = py;
190c4762a1bSJed Brown         cnt++;
191c4762a1bSJed Brown       }
192c4762a1bSJed Brown     }
1939566063dSJacob Faibussowitsch     PetscCall(PetscRandomDestroy(&rand));
1949566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, DMSwarmPICField_coor, &bs, NULL, (void **)&array));
1959566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(coors, &LA_coors));
1969566063dSJacob Faibussowitsch     PetscCall(DMSwarmSetLocalSizes(dms, cnt, 4));
197c4762a1bSJed Brown 
1989566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &nlocal));
1999566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "itag", &bs, NULL, (void **)&iarray));
200*ad540459SPierre Jolivet     for (p = 0; p < nlocal; p++) iarray[p] = (PetscInt)rank;
2019566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "itag", &bs, NULL, (void **)&iarray));
202c4762a1bSJed Brown   }
203c4762a1bSJed Brown 
2049566063dSJacob Faibussowitsch   PetscCall(DMView(dms, PETSC_VIEWER_STDOUT_WORLD));
2059566063dSJacob Faibussowitsch   PetscCall(SwarmViewGP(dms, "step0"));
206c4762a1bSJed Brown 
207c4762a1bSJed Brown   dt = 0.1;
208c4762a1bSJed Brown   for (tk = 1; tk < 20; tk++) {
209c4762a1bSJed Brown     char prefix[PETSC_MAX_PATH_LEN];
21063a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Step %" PetscInt_FMT " \n", tk));
211c4762a1bSJed Brown     /* push points */
2129566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &nlocal));
2139566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, DMSwarmPICField_coor, &bs, NULL, (void **)&array));
214c4762a1bSJed Brown     for (p = 0; p < nlocal; p++) {
215c4762a1bSJed Brown       PetscReal cx, cy, vx, vy;
216c4762a1bSJed Brown 
217c4762a1bSJed Brown       cx = array[2 * p];
218c4762a1bSJed Brown       cy = array[2 * p + 1];
219c4762a1bSJed Brown       vx = cy;
220c4762a1bSJed Brown       vy = -cx;
221c4762a1bSJed Brown 
222c4762a1bSJed Brown       array[2 * p] += dt * vx;
223c4762a1bSJed Brown       array[2 * p + 1] += dt * vy;
224c4762a1bSJed Brown     }
2259566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, DMSwarmPICField_coor, &bs, NULL, (void **)&array));
226c4762a1bSJed Brown 
227c4762a1bSJed Brown     /* migrate points */
2289566063dSJacob Faibussowitsch     PetscCall(DMSwarmMigrate(dms, PETSC_TRUE));
229c4762a1bSJed Brown     /* view points */
23063a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN - 1, "step%" PetscInt_FMT, tk));
231c4762a1bSJed Brown     /* should use the regular SwarmView() api, not one for a particular type */
2329566063dSJacob Faibussowitsch     PetscCall(SwarmViewGP(dms, prefix));
233c4762a1bSJed Brown   }
2349566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmregular));
2359566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmcell));
2369566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dms));
237c4762a1bSJed Brown   PetscFunctionReturn(0);
238c4762a1bSJed Brown }
239c4762a1bSJed Brown 
2409371c9d4SSatish Balay int main(int argc, char **argv) {
241327415f7SBarry Smith   PetscFunctionBeginUser;
2429566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
2439566063dSJacob Faibussowitsch   PetscCall(ex3_1());
2449566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
245b122ec5aSJacob Faibussowitsch   return 0;
246c4762a1bSJed Brown }
247c4762a1bSJed Brown 
248c4762a1bSJed Brown /*TEST
249c4762a1bSJed Brown 
250c4762a1bSJed Brown    build:
251c4762a1bSJed Brown       requires: double !complex
252c4762a1bSJed Brown 
253c4762a1bSJed Brown    test:
254c4762a1bSJed Brown       filter: grep -v atomic
255c4762a1bSJed Brown       filter_output: grep -v atomic
256c4762a1bSJed Brown TEST*/
257