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