1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /* 5c4762a1bSJed Brown Demonstrates using DM_BOUNDARY_GHOSTED how to handle a rotated boundary conditions where one edge 6c4762a1bSJed Brown is connected to its immediate neighbor 7c4762a1bSJed Brown 8c4762a1bSJed Brown Consider the domain (with natural numbering) 9c4762a1bSJed Brown 10c4762a1bSJed Brown 6 7 8 11c4762a1bSJed Brown 3 4 5 12c4762a1bSJed Brown 0 1 2 13c4762a1bSJed Brown 14c4762a1bSJed Brown The ghost points along the bottom (directly below the three columns above) should be 0 3 and 6 15c4762a1bSJed Brown while the ghost points along the left side should be 0 1 2 16c4762a1bSJed Brown 17c4762a1bSJed Brown Note that the ghosted local vectors extend in both the x and y directions so, for example if we have a 18c4762a1bSJed Brown single MPI process the ghosted vector has (in the original natural numbering) 19c4762a1bSJed Brown 20c4762a1bSJed Brown x x x x x 21c4762a1bSJed Brown 2 6 7 8 x 22c4762a1bSJed Brown 1 3 4 5 x 23c4762a1bSJed Brown 0 0 1 2 x 24c4762a1bSJed Brown x 0 3 6 x 25c4762a1bSJed Brown 26c4762a1bSJed Brown where x indicates a location that is not updated by the communication and should be used. 27c4762a1bSJed Brown 28c4762a1bSJed Brown For this to make sense the number of grid points in the x and y directions must be the same 29c4762a1bSJed Brown 30c4762a1bSJed Brown This ghost point mapping was suggested by: Wenbo Zhao <zhaowenbo.npic@gmail.com> 31c4762a1bSJed Brown */ 32c4762a1bSJed Brown 33c4762a1bSJed Brown #include <petscdm.h> 34c4762a1bSJed Brown #include <petscdmda.h> 35c4762a1bSJed Brown 36c4762a1bSJed Brown int main(int argc,char **argv) 37c4762a1bSJed Brown { 38c4762a1bSJed Brown PetscInt M = 6; 39c4762a1bSJed Brown DM da; 40c4762a1bSJed Brown Vec local,global,natural; 41c4762a1bSJed Brown PetscInt i,start,end,*ifrom,x,y,xm,ym; 42c4762a1bSJed Brown PetscScalar *xnatural; 43c4762a1bSJed Brown IS from,to; 44c4762a1bSJed Brown AO ao; 45c4762a1bSJed Brown VecScatter scatter1,scatter2; 46c4762a1bSJed Brown PetscViewer subviewer; 47c4762a1bSJed Brown 48*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 49c4762a1bSJed Brown 50c4762a1bSJed Brown /* Create distributed array and get vectors */ 515f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_GHOSTED,DM_BOUNDARY_GHOSTED,DMDA_STENCIL_STAR,M,M,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da)); 525f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da)); 535f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(da,&global)); 545f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLocalVector(da,&local)); 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* construct global to local scatter for the left side of the domain to the ghost on the bottom */ 575f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&x,&y,NULL,&xm,&ym,NULL)); 58c4762a1bSJed Brown if (!y) { /* only processes on the bottom of the domain fill up the ghost locations */ 595f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PETSC_COMM_SELF,xm,1,1,&to)); 60c4762a1bSJed Brown } else { 615f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PETSC_COMM_SELF,0,0,0,&to)); 62c4762a1bSJed Brown } 635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(xm,&ifrom)); 64c4762a1bSJed Brown for (i=x;i<x+xm;i++) { 65c4762a1bSJed Brown ifrom[i-x] = M*i; 66c4762a1bSJed Brown } 675f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetAO(da,&ao)); 685f80ce2aSJacob Faibussowitsch CHKERRQ(AOApplicationToPetsc(ao,xm,ifrom)); 69c4762a1bSJed Brown if (!y) { 705f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,xm,ifrom,PETSC_OWN_POINTER,&from)); 71c4762a1bSJed Brown } else { 725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ifrom)); 735f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,0,NULL,PETSC_COPY_VALUES,&from)); 74c4762a1bSJed Brown } 755f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterCreate(global,from,local,to,&scatter1)); 765f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&to)); 775f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&from)); 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* construct global to local scatter for the bottom side of the domain to the ghost on the right */ 80c4762a1bSJed Brown if (!x) { /* only processes on the left side of the domain fill up the ghost locations */ 815f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PETSC_COMM_SELF,ym,xm+2,xm+2,&to)); 82c4762a1bSJed Brown } else { 835f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PETSC_COMM_SELF,0,0,0,&to)); 84c4762a1bSJed Brown } 855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(ym,&ifrom)); 86c4762a1bSJed Brown for (i=y;i<y+ym;i++) { 87c4762a1bSJed Brown ifrom[i-y] = i; 88c4762a1bSJed Brown } 895f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetAO(da,&ao)); 905f80ce2aSJacob Faibussowitsch CHKERRQ(AOApplicationToPetsc(ao,ym,ifrom)); 91c4762a1bSJed Brown if (!x) { 925f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,ym,ifrom,PETSC_OWN_POINTER,&from)); 93c4762a1bSJed Brown } else { 945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ifrom)); 955f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,0,NULL,PETSC_COPY_VALUES,&from)); 96c4762a1bSJed Brown } 975f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterCreate(global,from,local,to,&scatter2)); 985f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&to)); 995f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&from)); 100c4762a1bSJed Brown 101c4762a1bSJed Brown /* 102c4762a1bSJed Brown fill the global vector with the natural global numbering for each local entry 103c4762a1bSJed Brown this is only done for testing purposes since it is easy to see if the scatter worked correctly 104c4762a1bSJed Brown */ 1055f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreateNaturalVector(da,&natural)); 1065f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetOwnershipRange(natural,&start,&end)); 1075f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(natural,&xnatural)); 108c4762a1bSJed Brown for (i=start; i<end; i++) { 109c4762a1bSJed Brown xnatural[i-start] = i; 110c4762a1bSJed Brown } 1115f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(natural,&xnatural)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(DMDANaturalToGlobalBegin(da,natural,INSERT_VALUES,global)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(DMDANaturalToGlobalEnd(da,natural,INSERT_VALUES,global)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&natural)); 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* scatter from global to local */ 1175f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(scatter1,global,local,INSERT_VALUES,SCATTER_FORWARD)); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(scatter1,global,local,INSERT_VALUES,SCATTER_FORWARD)); 1195f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(scatter2,global,local,INSERT_VALUES,SCATTER_FORWARD)); 1205f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(scatter2,global,local,INSERT_VALUES,SCATTER_FORWARD)); 121c4762a1bSJed Brown /* 122c4762a1bSJed Brown normally here you would also call 1235f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(da,global,INSERT_VALUES,local)); 1245f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(da,global,INSERT_VALUES,local)); 125c4762a1bSJed Brown to update all the interior ghost cells between neighboring processes. 126c4762a1bSJed Brown We don't do it here since this is only a test of "special" ghost points. 127c4762a1bSJed Brown */ 128c4762a1bSJed Brown 129c4762a1bSJed Brown /* view each local ghosted vector */ 1305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&subviewer)); 1315f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(local,subviewer)); 1325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&subviewer)); 133c4762a1bSJed Brown 1345f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&scatter1)); 1355f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&scatter2)); 1365f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&local)); 1375f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&global)); 1385f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&da)); 139*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 140*b122ec5aSJacob Faibussowitsch return 0; 141c4762a1bSJed Brown } 142c4762a1bSJed Brown 143c4762a1bSJed Brown /*TEST 144c4762a1bSJed Brown 145c4762a1bSJed Brown test: 146c4762a1bSJed Brown 147c4762a1bSJed Brown test: 148c4762a1bSJed Brown suffix: 2 149c4762a1bSJed Brown nsize: 2 150c4762a1bSJed Brown 151c4762a1bSJed Brown test: 152c4762a1bSJed Brown suffix: 4 153c4762a1bSJed Brown nsize: 4 154c4762a1bSJed Brown 155c4762a1bSJed Brown test: 156c4762a1bSJed Brown suffix: 9 157c4762a1bSJed Brown nsize: 9 158c4762a1bSJed Brown 159c4762a1bSJed Brown TEST*/ 160