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