xref: /petsc/src/dm/tutorials/swarm_ex3.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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 
11c4762a1bSJed Brown PetscErrorCode _DMLocatePoints_DMDARegular_IS(DM dm,Vec pos,IS *iscell)
12c4762a1bSJed Brown {
13c4762a1bSJed Brown   PetscInt       p,n,bs,npoints,si,sj,milocal,mjlocal,mx,my;
14c4762a1bSJed Brown   DM             dmregular;
15c4762a1bSJed Brown   PetscInt       *cellidx;
16c4762a1bSJed Brown   const PetscScalar *coor;
17c4762a1bSJed Brown   PetscReal      dx,dy;
18c4762a1bSJed Brown   PetscMPIInt    rank;
19c4762a1bSJed Brown 
20362febeeSStefano Zampini   PetscFunctionBegin;
215f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
225f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(pos,&n));
235f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetBlockSize(pos,&bs));
24c4762a1bSJed Brown   npoints = n/bs;
25c4762a1bSJed Brown 
265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(npoints,&cellidx));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetApplicationContext(dm,&dmregular));
285f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(dmregular,&si,&sj,NULL,&milocal,&mjlocal,NULL));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(dmregular,NULL,&mx,&my,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL));
30c4762a1bSJed Brown 
31c4762a1bSJed Brown   dx = 2.0/((PetscReal)mx);
32c4762a1bSJed Brown   dy = 2.0/((PetscReal)my);
33c4762a1bSJed Brown 
345f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(pos,&coor));
35c4762a1bSJed Brown   for (p=0; p<npoints; p++) {
36c4762a1bSJed Brown     PetscReal coorx,coory;
37c4762a1bSJed Brown     PetscInt  mi,mj;
38c4762a1bSJed Brown 
39c4762a1bSJed Brown     coorx = PetscRealPart(coor[2*p]);
40c4762a1bSJed Brown     coory = PetscRealPart(coor[2*p+1]);
41c4762a1bSJed Brown 
42c4762a1bSJed Brown     mi = (PetscInt)( (coorx - (-1.0))/dx);
43c4762a1bSJed Brown     mj = (PetscInt)( (coory - (-1.0))/dy);
44c4762a1bSJed Brown 
45c4762a1bSJed Brown     cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;
46c4762a1bSJed Brown 
47c4762a1bSJed Brown     if ((mj >= sj) && (mj < sj + mjlocal)) {
48c4762a1bSJed Brown       if ((mi >= si) && (mi < si + milocal)) {
49c4762a1bSJed Brown         cellidx[p] = (mi-si) + (mj-sj) * milocal;
50c4762a1bSJed Brown       }
51c4762a1bSJed Brown     }
52c4762a1bSJed Brown     if (coorx < -1.0) cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;
53c4762a1bSJed Brown     if (coorx >  1.0) cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;
54c4762a1bSJed Brown     if (coory < -1.0) cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;
55c4762a1bSJed Brown     if (coory >  1.0) cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;
56c4762a1bSJed Brown   }
575f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(pos,&coor));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,npoints,cellidx,PETSC_OWN_POINTER,iscell));
59c4762a1bSJed Brown   PetscFunctionReturn(0);
60c4762a1bSJed Brown }
61c4762a1bSJed Brown 
62c4762a1bSJed Brown PetscErrorCode DMLocatePoints_DMDARegular(DM dm,Vec pos,DMPointLocationType ltype, PetscSF cellSF)
63c4762a1bSJed Brown {
64c4762a1bSJed Brown   IS             iscell;
65c4762a1bSJed Brown   PetscSFNode    *cells;
66c4762a1bSJed Brown   PetscInt       p,bs,npoints,nfound;
67c4762a1bSJed Brown   const PetscInt *boxCells;
68c4762a1bSJed Brown 
69362febeeSStefano Zampini   PetscFunctionBegin;
705f80ce2aSJacob Faibussowitsch   CHKERRQ(_DMLocatePoints_DMDARegular_IS(dm,pos,&iscell));
715f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(pos,&npoints));
725f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetBlockSize(pos,&bs));
73c4762a1bSJed Brown   npoints = npoints / bs;
74c4762a1bSJed Brown 
755f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(npoints, &cells));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetIndices(iscell, &boxCells));
77c4762a1bSJed Brown 
78c4762a1bSJed Brown   for (p=0; p<npoints; p++) {
79c4762a1bSJed Brown     cells[p].rank  = 0;
80c4762a1bSJed Brown     cells[p].index = DMLOCATEPOINT_POINT_NOT_FOUND;
81c4762a1bSJed Brown     cells[p].index = boxCells[p];
82c4762a1bSJed Brown   }
835f80ce2aSJacob Faibussowitsch   CHKERRQ(ISRestoreIndices(iscell, &boxCells));
845f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&iscell));
85c4762a1bSJed Brown   nfound = npoints;
865f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetGraph(cellSF, npoints, nfound, NULL, PETSC_OWN_POINTER, cells, PETSC_OWN_POINTER));
875f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&iscell));
88c4762a1bSJed Brown   PetscFunctionReturn(0);
89c4762a1bSJed Brown }
90c4762a1bSJed Brown 
91c4762a1bSJed Brown PetscErrorCode DMGetNeighbors_DMDARegular(DM dm,PetscInt *nneighbors,const PetscMPIInt **neighbors)
92c4762a1bSJed Brown {
93c4762a1bSJed Brown   DM             dmregular;
94c4762a1bSJed Brown 
95362febeeSStefano Zampini   PetscFunctionBegin;
965f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetApplicationContext(dm,&dmregular));
975f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetNeighbors(dmregular,nneighbors,neighbors));
98c4762a1bSJed Brown   PetscFunctionReturn(0);
99c4762a1bSJed Brown }
100c4762a1bSJed Brown 
101c4762a1bSJed Brown PetscErrorCode SwarmViewGP(DM dms,const char prefix[])
102c4762a1bSJed Brown {
103c4762a1bSJed Brown   PetscReal      *array;
104c4762a1bSJed Brown   PetscInt       *iarray;
105c4762a1bSJed Brown   PetscInt       npoints,p,bs;
106c4762a1bSJed Brown   FILE           *fp;
107c4762a1bSJed Brown   char           name[PETSC_MAX_PATH_LEN];
108c4762a1bSJed Brown   PetscMPIInt    rank;
109c4762a1bSJed Brown 
110362febeeSStefano Zampini   PetscFunctionBegin;
1115f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
1125f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"%s-rank%d.gp",prefix,rank));
113c4762a1bSJed Brown   fp = fopen(name,"w");
11428b400f6SJacob Faibussowitsch   PetscCheck(fp,PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file %s",name);
1155f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetLocalSize(dms,&npoints));
1165f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(dms,DMSwarmPICField_coor,&bs,NULL,(void**)&array));
1175f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(dms,"itag",NULL,NULL,(void**)&iarray));
118c4762a1bSJed Brown   for (p=0; p<npoints; p++) {
119c4762a1bSJed Brown     fprintf(fp,"%+1.4e %+1.4e %1.4e\n",array[2*p],array[2*p+1],(double)iarray[p]);
120c4762a1bSJed Brown   }
1215f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(dms,"itag",NULL,NULL,(void**)&iarray));
1225f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(dms,DMSwarmPICField_coor,&bs,NULL,(void**)&array));
123c4762a1bSJed Brown   fclose(fp);
124c4762a1bSJed Brown   PetscFunctionReturn(0);
125c4762a1bSJed Brown }
126c4762a1bSJed Brown 
127c4762a1bSJed Brown /*
128c4762a1bSJed Brown  Create a DMShell and attach a regularly spaced DMDA for point location
129c4762a1bSJed Brown  Override methods for point location
130c4762a1bSJed Brown */
131c4762a1bSJed Brown PetscErrorCode ex3_1(void)
132c4762a1bSJed Brown {
133c4762a1bSJed Brown   DM             dms,dmcell,dmregular;
134c4762a1bSJed Brown   PetscMPIInt    rank;
135c4762a1bSJed Brown   PetscInt       p,bs,nlocal,overlap,mx,tk;
136c4762a1bSJed Brown   PetscReal      dx;
137c4762a1bSJed Brown   PetscReal      *array,dt;
138c4762a1bSJed Brown   PetscInt       *iarray;
139c4762a1bSJed Brown   PetscRandom    rand;
140c4762a1bSJed Brown 
141362febeeSStefano Zampini   PetscFunctionBegin;
1425f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   /* Create a regularly spaced DMDA */
145c4762a1bSJed Brown   mx = 40;
146c4762a1bSJed Brown   overlap = 0;
1475f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,mx,mx,PETSC_DECIDE,PETSC_DECIDE,1,overlap,NULL,NULL,&dmregular));
1485f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(dmregular));
1495f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(dmregular));
150c4762a1bSJed Brown 
151c4762a1bSJed Brown   dx = 2.0/((PetscReal)mx);
1525f80ce2aSJacob Faibussowitsch   CHKERRQ(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));
153c4762a1bSJed Brown 
154c4762a1bSJed Brown   /* Create a DMShell for point location purposes */
1555f80ce2aSJacob Faibussowitsch   CHKERRQ(DMShellCreate(PETSC_COMM_WORLD,&dmcell));
1565f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetApplicationContext(dmcell,dmregular));
157c4762a1bSJed Brown   dmcell->ops->locatepoints = DMLocatePoints_DMDARegular;
158c4762a1bSJed Brown   dmcell->ops->getneighbors = DMGetNeighbors_DMDARegular;
159c4762a1bSJed Brown 
160c4762a1bSJed Brown   /* Create the swarm */
1615f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(PETSC_COMM_WORLD,&dms));
1625f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(dms,DMSWARM));
1635f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetDimension(dms,2));
164c4762a1bSJed Brown 
1655f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetType(dms,DMSWARM_PIC));
1665f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetCellDM(dms,dmcell));
167c4762a1bSJed Brown 
1685f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRegisterPetscDatatypeField(dms,"itag",1,PETSC_INT));
1695f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmFinalizeFieldRegister(dms));
170c4762a1bSJed Brown   {
171c4762a1bSJed Brown     PetscInt  si,sj,milocal,mjlocal;
172c4762a1bSJed Brown     const PetscScalar *LA_coors;
173c4762a1bSJed Brown     Vec       coors;
174c4762a1bSJed Brown     PetscInt  cnt;
175c4762a1bSJed Brown 
1765f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAGetCorners(dmregular,&si,&sj,NULL,&milocal,&mjlocal,NULL));
1775f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetCoordinates(dmregular,&coors));
1785f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArrayRead(coors,&LA_coors));
1795f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmSetLocalSizes(dms,milocal*mjlocal,4));
1805f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmGetLocalSize(dms,&nlocal));
1815f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmGetField(dms,DMSwarmPICField_coor,&bs,NULL,(void**)&array));
182c4762a1bSJed Brown     cnt = 0;
1835f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF,&rand));
1845f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomSetInterval(rand,-1.0,1.0));
185c4762a1bSJed Brown     for (p=0; p<nlocal; p++) {
186c4762a1bSJed Brown       PetscReal px,py,rx,ry,r2;
187c4762a1bSJed Brown 
1885f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscRandomGetValueReal(rand,&rx));
1895f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscRandomGetValueReal(rand,&ry));
190c4762a1bSJed Brown 
191c4762a1bSJed Brown       px = PetscRealPart(LA_coors[2*p+0]) + 0.1*rx*dx;
192c4762a1bSJed Brown       py = PetscRealPart(LA_coors[2*p+1]) + 0.1*ry*dx;
193c4762a1bSJed Brown 
194c4762a1bSJed Brown       r2 = px*px + py*py;
195c4762a1bSJed Brown       if (r2 < 0.75*0.75) {
196c4762a1bSJed Brown         array[bs*cnt+0] = px;
197c4762a1bSJed Brown         array[bs*cnt+1] = py;
198c4762a1bSJed Brown         cnt++;
199c4762a1bSJed Brown       }
200c4762a1bSJed Brown     }
2015f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomDestroy(&rand));
2025f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmRestoreField(dms,DMSwarmPICField_coor,&bs,NULL,(void**)&array));
2035f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArrayRead(coors,&LA_coors));
2045f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmSetLocalSizes(dms,cnt,4));
205c4762a1bSJed Brown 
2065f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmGetLocalSize(dms,&nlocal));
2075f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmGetField(dms,"itag",&bs,NULL,(void**)&iarray));
208c4762a1bSJed Brown     for (p=0; p<nlocal; p++) {
209c4762a1bSJed Brown       iarray[p] = (PetscInt)rank;
210c4762a1bSJed Brown     }
2115f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmRestoreField(dms,"itag",&bs,NULL,(void**)&iarray));
212c4762a1bSJed Brown   }
213c4762a1bSJed Brown 
2145f80ce2aSJacob Faibussowitsch   CHKERRQ(DMView(dms,PETSC_VIEWER_STDOUT_WORLD));
2155f80ce2aSJacob Faibussowitsch   CHKERRQ(SwarmViewGP(dms,"step0"));
216c4762a1bSJed Brown 
217c4762a1bSJed Brown   dt = 0.1;
218c4762a1bSJed Brown   for (tk=1; tk<20; tk++) {
219c4762a1bSJed Brown     char prefix[PETSC_MAX_PATH_LEN];
2205f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Step %D \n",tk));
221c4762a1bSJed Brown     /* push points */
2225f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmGetLocalSize(dms,&nlocal));
2235f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmGetField(dms,DMSwarmPICField_coor,&bs,NULL,(void**)&array));
224c4762a1bSJed Brown     for (p=0; p<nlocal; p++) {
225c4762a1bSJed Brown       PetscReal cx,cy,vx,vy;
226c4762a1bSJed Brown 
227c4762a1bSJed Brown       cx = array[2*p];
228c4762a1bSJed Brown       cy = array[2*p+1];
229c4762a1bSJed Brown       vx =  cy;
230c4762a1bSJed Brown       vy = -cx;
231c4762a1bSJed Brown 
232c4762a1bSJed Brown       array[2*p  ] += dt * vx;
233c4762a1bSJed Brown       array[2*p+1] += dt * vy;
234c4762a1bSJed Brown     }
2355f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmRestoreField(dms,DMSwarmPICField_coor,&bs,NULL,(void**)&array));
236c4762a1bSJed Brown 
237c4762a1bSJed Brown     /* migrate points */
2385f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmMigrate(dms,PETSC_TRUE));
239c4762a1bSJed Brown     /* view points */
2405f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSNPrintf(prefix,PETSC_MAX_PATH_LEN-1,"step%d",tk));
241c4762a1bSJed Brown     /* should use the regular SwarmView() api, not one for a particular type */
2425f80ce2aSJacob Faibussowitsch     CHKERRQ(SwarmViewGP(dms,prefix));
243c4762a1bSJed Brown   }
2445f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dmregular));
2455f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dmcell));
2465f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dms));
247c4762a1bSJed Brown   PetscFunctionReturn(0);
248c4762a1bSJed Brown }
249c4762a1bSJed Brown 
250c4762a1bSJed Brown int main(int argc,char **argv)
251c4762a1bSJed Brown {
252c4762a1bSJed Brown 
253*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
2545f80ce2aSJacob Faibussowitsch   CHKERRQ(ex3_1());
255*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
256*b122ec5aSJacob Faibussowitsch   return 0;
257c4762a1bSJed Brown }
258c4762a1bSJed Brown 
259c4762a1bSJed Brown /*TEST
260c4762a1bSJed Brown 
261c4762a1bSJed Brown    build:
262c4762a1bSJed Brown       requires: double !complex
263c4762a1bSJed Brown 
264c4762a1bSJed Brown    test:
265c4762a1bSJed Brown       filter: grep -v atomic
266c4762a1bSJed Brown       filter_output: grep -v atomic
267c4762a1bSJed Brown TEST*/
268