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