xref: /petsc/src/dm/tutorials/swarm_ex3.c (revision ffa8c5705e8ab2cf85ee1d14dbe507a6e2eb5283)
1 
2 static char help[] = "Tests DMSwarm with DMShell\n\n";
3 
4 #include <petscsf.h>
5 #include <petscdm.h>
6 #include <petscdmshell.h>
7 #include <petscdmda.h>
8 #include <petscdmswarm.h>
9 #include <petsc/private/dmimpl.h>
10 
11 PetscErrorCode _DMLocatePoints_DMDARegular_IS(DM dm,Vec pos,IS *iscell)
12 {
13   PetscInt       p,n,bs,npoints,si,sj,milocal,mjlocal,mx,my;
14   DM             dmregular;
15   PetscInt       *cellidx;
16   const PetscScalar *coor;
17   PetscReal      dx,dy;
18   PetscMPIInt    rank;
19 
20   PetscFunctionBegin;
21   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
22   PetscCall(VecGetLocalSize(pos,&n));
23   PetscCall(VecGetBlockSize(pos,&bs));
24   npoints = n/bs;
25 
26   PetscCall(PetscMalloc1(npoints,&cellidx));
27   PetscCall(DMGetApplicationContext(dm,&dmregular));
28   PetscCall(DMDAGetCorners(dmregular,&si,&sj,NULL,&milocal,&mjlocal,NULL));
29   PetscCall(DMDAGetInfo(dmregular,NULL,&mx,&my,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL));
30 
31   dx = 2.0/((PetscReal)mx);
32   dy = 2.0/((PetscReal)my);
33 
34   PetscCall(VecGetArrayRead(pos,&coor));
35   for (p=0; p<npoints; p++) {
36     PetscReal coorx,coory;
37     PetscInt  mi,mj;
38 
39     coorx = PetscRealPart(coor[2*p]);
40     coory = PetscRealPart(coor[2*p+1]);
41 
42     mi = (PetscInt)( (coorx - (-1.0))/dx);
43     mj = (PetscInt)( (coory - (-1.0))/dy);
44 
45     cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;
46 
47     if ((mj >= sj) && (mj < sj + mjlocal)) {
48       if ((mi >= si) && (mi < si + milocal)) {
49         cellidx[p] = (mi-si) + (mj-sj) * milocal;
50       }
51     }
52     if (coorx < -1.0) cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;
53     if (coorx >  1.0) cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;
54     if (coory < -1.0) cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;
55     if (coory >  1.0) cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;
56   }
57   PetscCall(VecRestoreArrayRead(pos,&coor));
58   PetscCall(ISCreateGeneral(PETSC_COMM_SELF,npoints,cellidx,PETSC_OWN_POINTER,iscell));
59   PetscFunctionReturn(0);
60 }
61 
62 PetscErrorCode DMLocatePoints_DMDARegular(DM dm,Vec pos,DMPointLocationType ltype, PetscSF cellSF)
63 {
64   IS             iscell;
65   PetscSFNode    *cells;
66   PetscInt       p,bs,npoints,nfound;
67   const PetscInt *boxCells;
68 
69   PetscFunctionBegin;
70   PetscCall(_DMLocatePoints_DMDARegular_IS(dm,pos,&iscell));
71   PetscCall(VecGetLocalSize(pos,&npoints));
72   PetscCall(VecGetBlockSize(pos,&bs));
73   npoints = npoints / bs;
74 
75   PetscCall(PetscMalloc1(npoints, &cells));
76   PetscCall(ISGetIndices(iscell, &boxCells));
77 
78   for (p=0; p<npoints; p++) {
79     cells[p].rank  = 0;
80     cells[p].index = DMLOCATEPOINT_POINT_NOT_FOUND;
81     cells[p].index = boxCells[p];
82   }
83   PetscCall(ISRestoreIndices(iscell, &boxCells));
84   PetscCall(ISDestroy(&iscell));
85   nfound = npoints;
86   PetscCall(PetscSFSetGraph(cellSF, npoints, nfound, NULL, PETSC_OWN_POINTER, cells, PETSC_OWN_POINTER));
87   PetscCall(ISDestroy(&iscell));
88   PetscFunctionReturn(0);
89 }
90 
91 PetscErrorCode DMGetNeighbors_DMDARegular(DM dm,PetscInt *nneighbors,const PetscMPIInt **neighbors)
92 {
93   DM             dmregular;
94 
95   PetscFunctionBegin;
96   PetscCall(DMGetApplicationContext(dm,&dmregular));
97   PetscCall(DMGetNeighbors(dmregular,nneighbors,neighbors));
98   PetscFunctionReturn(0);
99 }
100 
101 PetscErrorCode SwarmViewGP(DM dms,const char prefix[])
102 {
103   PetscReal      *array;
104   PetscInt       *iarray;
105   PetscInt       npoints,p,bs;
106   FILE           *fp;
107   char           name[PETSC_MAX_PATH_LEN];
108   PetscMPIInt    rank;
109 
110   PetscFunctionBegin;
111   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
112   PetscCall(PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"%s-rank%d.gp",prefix,rank));
113   fp = fopen(name,"w");
114   PetscCheck(fp,PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file %s",name);
115   PetscCall(DMSwarmGetLocalSize(dms,&npoints));
116   PetscCall(DMSwarmGetField(dms,DMSwarmPICField_coor,&bs,NULL,(void**)&array));
117   PetscCall(DMSwarmGetField(dms,"itag",NULL,NULL,(void**)&iarray));
118   for (p=0; p<npoints; p++) {
119     fprintf(fp,"%+1.4e %+1.4e %1.4e\n",array[2*p],array[2*p+1],(double)iarray[p]);
120   }
121   PetscCall(DMSwarmRestoreField(dms,"itag",NULL,NULL,(void**)&iarray));
122   PetscCall(DMSwarmRestoreField(dms,DMSwarmPICField_coor,&bs,NULL,(void**)&array));
123   fclose(fp);
124   PetscFunctionReturn(0);
125 }
126 
127 /*
128  Create a DMShell and attach a regularly spaced DMDA for point location
129  Override methods for point location
130 */
131 PetscErrorCode ex3_1(void)
132 {
133   DM             dms,dmcell,dmregular;
134   PetscMPIInt    rank;
135   PetscInt       p,bs,nlocal,overlap,mx,tk;
136   PetscReal      dx;
137   PetscReal      *array,dt;
138   PetscInt       *iarray;
139   PetscRandom    rand;
140 
141   PetscFunctionBegin;
142   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
143 
144   /* Create a regularly spaced DMDA */
145   mx = 40;
146   overlap = 0;
147   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));
148   PetscCall(DMSetFromOptions(dmregular));
149   PetscCall(DMSetUp(dmregular));
150 
151   dx = 2.0/((PetscReal)mx);
152   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));
153 
154   /* Create a DMShell for point location purposes */
155   PetscCall(DMShellCreate(PETSC_COMM_WORLD,&dmcell));
156   PetscCall(DMSetApplicationContext(dmcell,dmregular));
157   dmcell->ops->locatepoints = DMLocatePoints_DMDARegular;
158   dmcell->ops->getneighbors = DMGetNeighbors_DMDARegular;
159 
160   /* Create the swarm */
161   PetscCall(DMCreate(PETSC_COMM_WORLD,&dms));
162   PetscCall(DMSetType(dms,DMSWARM));
163   PetscCall(DMSetDimension(dms,2));
164 
165   PetscCall(DMSwarmSetType(dms,DMSWARM_PIC));
166   PetscCall(DMSwarmSetCellDM(dms,dmcell));
167 
168   PetscCall(DMSwarmRegisterPetscDatatypeField(dms,"itag",1,PETSC_INT));
169   PetscCall(DMSwarmFinalizeFieldRegister(dms));
170   {
171     PetscInt  si,sj,milocal,mjlocal;
172     const PetscScalar *LA_coors;
173     Vec       coors;
174     PetscInt  cnt;
175 
176     PetscCall(DMDAGetCorners(dmregular,&si,&sj,NULL,&milocal,&mjlocal,NULL));
177     PetscCall(DMGetCoordinates(dmregular,&coors));
178     PetscCall(VecGetArrayRead(coors,&LA_coors));
179     PetscCall(DMSwarmSetLocalSizes(dms,milocal*mjlocal,4));
180     PetscCall(DMSwarmGetLocalSize(dms,&nlocal));
181     PetscCall(DMSwarmGetField(dms,DMSwarmPICField_coor,&bs,NULL,(void**)&array));
182     cnt = 0;
183     PetscCall(PetscRandomCreate(PETSC_COMM_SELF,&rand));
184     PetscCall(PetscRandomSetInterval(rand,-1.0,1.0));
185     for (p=0; p<nlocal; p++) {
186       PetscReal px,py,rx,ry,r2;
187 
188       PetscCall(PetscRandomGetValueReal(rand,&rx));
189       PetscCall(PetscRandomGetValueReal(rand,&ry));
190 
191       px = PetscRealPart(LA_coors[2*p+0]) + 0.1*rx*dx;
192       py = PetscRealPart(LA_coors[2*p+1]) + 0.1*ry*dx;
193 
194       r2 = px*px + py*py;
195       if (r2 < 0.75*0.75) {
196         array[bs*cnt+0] = px;
197         array[bs*cnt+1] = py;
198         cnt++;
199       }
200     }
201     PetscCall(PetscRandomDestroy(&rand));
202     PetscCall(DMSwarmRestoreField(dms,DMSwarmPICField_coor,&bs,NULL,(void**)&array));
203     PetscCall(VecRestoreArrayRead(coors,&LA_coors));
204     PetscCall(DMSwarmSetLocalSizes(dms,cnt,4));
205 
206     PetscCall(DMSwarmGetLocalSize(dms,&nlocal));
207     PetscCall(DMSwarmGetField(dms,"itag",&bs,NULL,(void**)&iarray));
208     for (p=0; p<nlocal; p++) {
209       iarray[p] = (PetscInt)rank;
210     }
211     PetscCall(DMSwarmRestoreField(dms,"itag",&bs,NULL,(void**)&iarray));
212   }
213 
214   PetscCall(DMView(dms,PETSC_VIEWER_STDOUT_WORLD));
215   PetscCall(SwarmViewGP(dms,"step0"));
216 
217   dt = 0.1;
218   for (tk=1; tk<20; tk++) {
219     char prefix[PETSC_MAX_PATH_LEN];
220     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Step %D \n",tk));
221     /* push points */
222     PetscCall(DMSwarmGetLocalSize(dms,&nlocal));
223     PetscCall(DMSwarmGetField(dms,DMSwarmPICField_coor,&bs,NULL,(void**)&array));
224     for (p=0; p<nlocal; p++) {
225       PetscReal cx,cy,vx,vy;
226 
227       cx = array[2*p];
228       cy = array[2*p+1];
229       vx =  cy;
230       vy = -cx;
231 
232       array[2*p  ] += dt * vx;
233       array[2*p+1] += dt * vy;
234     }
235     PetscCall(DMSwarmRestoreField(dms,DMSwarmPICField_coor,&bs,NULL,(void**)&array));
236 
237     /* migrate points */
238     PetscCall(DMSwarmMigrate(dms,PETSC_TRUE));
239     /* view points */
240     PetscCall(PetscSNPrintf(prefix,PETSC_MAX_PATH_LEN-1,"step%d",tk));
241     /* should use the regular SwarmView() api, not one for a particular type */
242     PetscCall(SwarmViewGP(dms,prefix));
243   }
244   PetscCall(DMDestroy(&dmregular));
245   PetscCall(DMDestroy(&dmcell));
246   PetscCall(DMDestroy(&dms));
247   PetscFunctionReturn(0);
248 }
249 
250 int main(int argc,char **argv)
251 {
252 
253   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
254   PetscCall(ex3_1());
255   PetscCall(PetscFinalize());
256   return 0;
257 }
258 
259 /*TEST
260 
261    build:
262       requires: double !complex
263 
264    test:
265       filter: grep -v atomic
266       filter_output: grep -v atomic
267 TEST*/
268