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 PetscErrorCode ierr; 19 PetscMPIInt rank; 20 21 ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 22 ierr = VecGetLocalSize(pos,&n);CHKERRQ(ierr); 23 ierr = VecGetBlockSize(pos,&bs);CHKERRQ(ierr); 24 npoints = n/bs; 25 26 ierr = PetscMalloc1(npoints,&cellidx);CHKERRQ(ierr); 27 ierr = DMGetApplicationContext(dm,(void**)&dmregular);CHKERRQ(ierr); 28 ierr = DMDAGetCorners(dmregular,&si,&sj,NULL,&milocal,&mjlocal,NULL);CHKERRQ(ierr); 29 ierr = DMDAGetInfo(dmregular,NULL,&mx,&my,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 30 31 dx = 2.0/((PetscReal)mx); 32 dy = 2.0/((PetscReal)my); 33 34 ierr = VecGetArrayRead(pos,&coor);CHKERRQ(ierr); 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 ierr = VecRestoreArrayRead(pos,&coor);CHKERRQ(ierr); 58 ierr = ISCreateGeneral(PETSC_COMM_SELF,npoints,cellidx,PETSC_OWN_POINTER,iscell);CHKERRQ(ierr); 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 PetscErrorCode ierr; 69 70 ierr = _DMLocatePoints_DMDARegular_IS(dm,pos,&iscell);CHKERRQ(ierr); 71 ierr = VecGetLocalSize(pos,&npoints);CHKERRQ(ierr); 72 ierr = VecGetBlockSize(pos,&bs);CHKERRQ(ierr); 73 npoints = npoints / bs; 74 75 ierr = PetscMalloc1(npoints, &cells);CHKERRQ(ierr); 76 ierr = ISGetIndices(iscell, &boxCells);CHKERRQ(ierr); 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 ierr = ISRestoreIndices(iscell, &boxCells);CHKERRQ(ierr); 84 ierr = ISDestroy(&iscell);CHKERRQ(ierr); 85 nfound = npoints; 86 ierr = PetscSFSetGraph(cellSF, npoints, nfound, NULL, PETSC_OWN_POINTER, cells, PETSC_OWN_POINTER);CHKERRQ(ierr); 87 ierr = ISDestroy(&iscell);CHKERRQ(ierr); 88 PetscFunctionReturn(0); 89 } 90 91 PetscErrorCode DMGetNeighbors_DMDARegular(DM dm,PetscInt *nneighbors,const PetscMPIInt **neighbors) 92 { 93 DM dmregular; 94 PetscErrorCode ierr; 95 96 ierr = DMGetApplicationContext(dm,(void**)&dmregular);CHKERRQ(ierr); 97 ierr = DMGetNeighbors(dmregular,nneighbors,neighbors);CHKERRQ(ierr); 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 PetscErrorCode ierr; 110 111 ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 112 ierr = PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"%s-rank%d.gp",prefix,rank);CHKERRQ(ierr); 113 fp = fopen(name,"w"); 114 if (!fp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file %s",name); 115 ierr = DMSwarmGetLocalSize(dms,&npoints);CHKERRQ(ierr); 116 ierr = DMSwarmGetField(dms,DMSwarmPICField_coor,&bs,NULL,(void**)&array);CHKERRQ(ierr); 117 ierr = DMSwarmGetField(dms,"itag",NULL,NULL,(void**)&iarray);CHKERRQ(ierr); 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 ierr = DMSwarmRestoreField(dms,"itag",NULL,NULL,(void**)&iarray);CHKERRQ(ierr); 122 ierr = DMSwarmRestoreField(dms,DMSwarmPICField_coor,&bs,NULL,(void**)&array);CHKERRQ(ierr); 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 PetscErrorCode ierr; 141 142 ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 143 144 /* Create a regularly spaced DMDA */ 145 mx = 40; 146 overlap = 0; 147 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); 148 ierr = DMSetFromOptions(dmregular);CHKERRQ(ierr); 149 ierr = DMSetUp(dmregular);CHKERRQ(ierr); 150 151 dx = 2.0/((PetscReal)mx); 152 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); 153 154 /* Create a DMShell for point location purposes */ 155 ierr = DMShellCreate(PETSC_COMM_WORLD,&dmcell);CHKERRQ(ierr); 156 ierr = DMSetApplicationContext(dmcell,(void*)dmregular);CHKERRQ(ierr); 157 dmcell->ops->locatepoints = DMLocatePoints_DMDARegular; 158 dmcell->ops->getneighbors = DMGetNeighbors_DMDARegular; 159 160 /* Create the swarm */ 161 ierr = DMCreate(PETSC_COMM_WORLD,&dms);CHKERRQ(ierr); 162 ierr = DMSetType(dms,DMSWARM);CHKERRQ(ierr); 163 ierr = DMSetDimension(dms,2);CHKERRQ(ierr); 164 165 ierr = DMSwarmSetType(dms,DMSWARM_PIC);CHKERRQ(ierr); 166 ierr = DMSwarmSetCellDM(dms,dmcell);CHKERRQ(ierr); 167 168 ierr = DMSwarmRegisterPetscDatatypeField(dms,"itag",1,PETSC_INT);CHKERRQ(ierr); 169 ierr = DMSwarmFinalizeFieldRegister(dms);CHKERRQ(ierr); 170 { 171 PetscInt si,sj,milocal,mjlocal; 172 const PetscScalar *LA_coors; 173 Vec coors; 174 PetscInt cnt; 175 176 ierr = DMDAGetCorners(dmregular,&si,&sj,NULL,&milocal,&mjlocal,NULL);CHKERRQ(ierr); 177 ierr = DMGetCoordinates(dmregular,&coors);CHKERRQ(ierr); 178 ierr = VecGetArrayRead(coors,&LA_coors);CHKERRQ(ierr); 179 ierr = DMSwarmSetLocalSizes(dms,milocal*mjlocal,4);CHKERRQ(ierr); 180 ierr = DMSwarmGetLocalSize(dms,&nlocal);CHKERRQ(ierr); 181 ierr = DMSwarmGetField(dms,DMSwarmPICField_coor,&bs,NULL,(void**)&array);CHKERRQ(ierr); 182 cnt = 0; 183 ierr = PetscRandomCreate(PETSC_COMM_SELF,&rand);CHKERRQ(ierr); 184 ierr = PetscRandomSetInterval(rand,-1.0,1.0);CHKERRQ(ierr); 185 for (p=0; p<nlocal; p++) { 186 PetscReal px,py,rx,ry,r2; 187 188 ierr = PetscRandomGetValueReal(rand,&rx);CHKERRQ(ierr); 189 ierr = PetscRandomGetValueReal(rand,&ry);CHKERRQ(ierr); 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 ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr); 202 ierr = DMSwarmRestoreField(dms,DMSwarmPICField_coor,&bs,NULL,(void**)&array);CHKERRQ(ierr); 203 ierr = VecRestoreArrayRead(coors,&LA_coors);CHKERRQ(ierr); 204 ierr = DMSwarmSetLocalSizes(dms,cnt,4);CHKERRQ(ierr); 205 206 ierr = DMSwarmGetLocalSize(dms,&nlocal);CHKERRQ(ierr); 207 ierr = DMSwarmGetField(dms,"itag",&bs,NULL,(void**)&iarray);CHKERRQ(ierr); 208 for (p=0; p<nlocal; p++) { 209 iarray[p] = (PetscInt)rank; 210 } 211 ierr = DMSwarmRestoreField(dms,"itag",&bs,NULL,(void**)&iarray);CHKERRQ(ierr); 212 } 213 214 ierr = DMView(dms,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 215 ierr = SwarmViewGP(dms,"step0");CHKERRQ(ierr); 216 217 dt = 0.1; 218 for (tk=1; tk<20; tk++) { 219 char prefix[PETSC_MAX_PATH_LEN]; 220 ierr = PetscPrintf(PETSC_COMM_WORLD,"Step %D \n",tk);CHKERRQ(ierr); 221 /* push points */ 222 ierr = DMSwarmGetLocalSize(dms,&nlocal);CHKERRQ(ierr); 223 ierr = DMSwarmGetField(dms,DMSwarmPICField_coor,&bs,NULL,(void**)&array);CHKERRQ(ierr); 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 ierr = DMSwarmRestoreField(dms,DMSwarmPICField_coor,&bs,NULL,(void**)&array);CHKERRQ(ierr); 236 237 /* migrate points */ 238 ierr = DMSwarmMigrate(dms,PETSC_TRUE);CHKERRQ(ierr); 239 /* view points */ 240 ierr = PetscSNPrintf(prefix,PETSC_MAX_PATH_LEN-1,"step%d",tk);CHKERRQ(ierr); 241 /* should use the regular SwarmView() api, not one for a particular type */ 242 ierr = SwarmViewGP(dms,prefix);CHKERRQ(ierr); 243 } 244 ierr = DMDestroy(&dmregular);CHKERRQ(ierr); 245 ierr = DMDestroy(&dmcell);CHKERRQ(ierr); 246 ierr = DMDestroy(&dms);CHKERRQ(ierr); 247 PetscFunctionReturn(0); 248 } 249 250 int main(int argc,char **argv) 251 { 252 PetscErrorCode ierr; 253 254 ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 255 ierr = ex3_1();CHKERRQ(ierr); 256 ierr = PetscFinalize(); 257 return ierr; 258 } 259 260 /*TEST 261 262 build: 263 requires: double !complex 264 265 test: 266 filter: grep -v atomic 267 filter_output: grep -v atomic 268 TEST*/ 269