1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests DMCreateDomainDecomposition.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /* 5c4762a1bSJed Brown Use the options 6c4762a1bSJed Brown -da_grid_x <nx> - number of grid points in x direction, if M < 0 7c4762a1bSJed Brown -da_grid_y <ny> - number of grid points in y direction, if N < 0 8c4762a1bSJed Brown -da_processors_x <MX> number of processors in x directio 9c4762a1bSJed Brown -da_processors_y <MY> number of processors in x direction 10c4762a1bSJed Brown */ 11c4762a1bSJed Brown 12c4762a1bSJed Brown #include <petscdm.h> 13c4762a1bSJed Brown #include <petscdmda.h> 14c4762a1bSJed Brown 15c4762a1bSJed Brown PetscErrorCode FillLocalSubdomain(DM da, Vec gvec) 16c4762a1bSJed Brown { 17c4762a1bSJed Brown DMDALocalInfo info; 18c4762a1bSJed Brown PetscMPIInt rank; 19c4762a1bSJed Brown PetscInt i,j,k,l; 20c4762a1bSJed Brown PetscErrorCode ierr; 21c4762a1bSJed Brown 22c4762a1bSJed Brown PetscFunctionBeginUser; 235f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 245f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetLocalInfo(da,&info)); 25c4762a1bSJed Brown 26c4762a1bSJed Brown if (info.dim == 3) { 27c4762a1bSJed Brown PetscScalar ***g; 285f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,gvec,&g)); 29c4762a1bSJed Brown /* loop over ghosts */ 30c4762a1bSJed Brown for (k=info.zs; k<info.zs+info.zm; k++) { 31c4762a1bSJed Brown for (j=info.ys; j<info.ys+info.ym; j++) { 32c4762a1bSJed Brown for (i=info.xs; i<info.xs+info.xm; i++) { 33c4762a1bSJed Brown g[k][j][info.dof*i+0] = i; 34c4762a1bSJed Brown g[k][j][info.dof*i+1] = j; 35c4762a1bSJed Brown g[k][j][info.dof*i+2] = k; 36c4762a1bSJed Brown } 37c4762a1bSJed Brown } 38c4762a1bSJed Brown } 395f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,gvec,&g)); 40c4762a1bSJed Brown } 41c4762a1bSJed Brown if (info.dim == 2) { 42c4762a1bSJed Brown PetscScalar **g; 435f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,gvec,&g)); 44c4762a1bSJed Brown /* loop over ghosts */ 45c4762a1bSJed Brown for (j=info.ys; j<info.ys+info.ym; j++) { 46c4762a1bSJed Brown for (i=info.xs; i<info.xs+info.xm; i++) { 47c4762a1bSJed Brown for (l = 0;l<info.dof;l++) { 48c4762a1bSJed Brown g[j][info.dof*i+0] = i; 49c4762a1bSJed Brown g[j][info.dof*i+1] = j; 50c4762a1bSJed Brown g[j][info.dof*i+2] = rank; 51c4762a1bSJed Brown } 52c4762a1bSJed Brown } 53c4762a1bSJed Brown } 545f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,gvec,&g)); 55c4762a1bSJed Brown } 56c4762a1bSJed Brown PetscFunctionReturn(0); 57c4762a1bSJed Brown } 58c4762a1bSJed Brown 59c4762a1bSJed Brown int main(int argc,char **argv) 60c4762a1bSJed Brown { 61c4762a1bSJed Brown PetscErrorCode ierr; 62c4762a1bSJed Brown DM da,*subda; 63c4762a1bSJed Brown PetscInt i,dim = 3; 643b5e53cdSSajid Ali PetscInt M = 25, N = 25, P = 25; 65c4762a1bSJed Brown PetscMPIInt size,rank; 66c4762a1bSJed Brown Vec v; 67c4762a1bSJed Brown Vec slvec,sgvec; 68c4762a1bSJed Brown IS *ois,*iis; 69c4762a1bSJed Brown VecScatter oscata; 70c4762a1bSJed Brown VecScatter *iscat,*oscat,*gscat; 71c4762a1bSJed Brown DMDALocalInfo info; 723b5e53cdSSajid Ali PetscBool patchis_offproc = PETSC_TRUE; 73c4762a1bSJed Brown 74*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL)); 76c4762a1bSJed Brown 77c4762a1bSJed Brown /* Create distributed array and get vectors */ 785f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 795f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 80c4762a1bSJed Brown if (dim == 2) { 815f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,M,N,PETSC_DECIDE,PETSC_DECIDE,3,1,NULL,NULL,&da)); 82c4762a1bSJed Brown } else if (dim == 3) { 835f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate3d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,M,N,P,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,3,1,NULL,NULL,NULL,&da)); 84c4762a1bSJed Brown } 855f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(da)); 865f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da)); 875f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetLocalInfo(da,&info)); 88c4762a1bSJed Brown 895f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateDomainDecomposition(da,NULL,NULL,&iis,&ois,&subda)); 905f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateDomainDecompositionScatters(da,1,subda,&iscat,&oscat,&gscat)); 91c4762a1bSJed Brown 92c4762a1bSJed Brown { 93c4762a1bSJed Brown DMDALocalInfo subinfo; 94c4762a1bSJed Brown MatStencil lower,upper; 953b5e53cdSSajid Ali IS patchis; 96c4762a1bSJed Brown Vec smallvec; 97c4762a1bSJed Brown Vec largevec; 98c4762a1bSJed Brown VecScatter patchscat; 99c4762a1bSJed Brown 1005f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetLocalInfo(subda[0],&subinfo)); 101c4762a1bSJed Brown 102c4762a1bSJed Brown lower.i = info.xs; 103c4762a1bSJed Brown lower.j = info.ys; 104c4762a1bSJed Brown lower.k = info.zs; 105c4762a1bSJed Brown upper.i = info.xs+info.xm; 106c4762a1bSJed Brown upper.j = info.ys+info.ym; 107c4762a1bSJed Brown upper.k = info.zs+info.zm; 108c4762a1bSJed Brown 109c4762a1bSJed Brown /* test the patch IS as a thing to scatter to/from */ 1105f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreatePatchIS(da,&lower,&upper,&patchis,patchis_offproc)); 1115f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(da,&largevec)); 112c4762a1bSJed Brown 1135f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_SELF,&smallvec)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(smallvec,info.dof*(upper.i - lower.i)*(upper.j - lower.j)*(upper.k - lower.k),PETSC_DECIDE)); 1155f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(smallvec)); 1165f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterCreate(smallvec,NULL,largevec,patchis,&patchscat)); 117c4762a1bSJed Brown 1185f80ce2aSJacob Faibussowitsch CHKERRQ(FillLocalSubdomain(subda[0],smallvec)); 1195f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(largevec,0)); 120c4762a1bSJed Brown 1215f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(patchscat,smallvec,largevec,ADD_VALUES,SCATTER_FORWARD)); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(patchscat,smallvec,largevec,ADD_VALUES,SCATTER_FORWARD)); 1235f80ce2aSJacob Faibussowitsch CHKERRQ(ISView(patchis,PETSC_VIEWER_STDOUT_WORLD)); 1245f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterView(patchscat,PETSC_VIEWER_STDOUT_WORLD)); 125c4762a1bSJed Brown 126c4762a1bSJed Brown for (i = 0; i < size; i++) { 127c4762a1bSJed Brown if (i == rank) { 1285f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(smallvec,PETSC_VIEWER_STDOUT_SELF)); 129c4762a1bSJed Brown } 1305f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Barrier(PETSC_COMM_WORLD)); 131c4762a1bSJed Brown } 132c4762a1bSJed Brown 1335f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Barrier(PETSC_COMM_WORLD)); 1345f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(largevec,PETSC_VIEWER_STDOUT_WORLD)); 135c4762a1bSJed Brown 1365f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&smallvec)); 1375f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(da,&largevec)); 1385f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&patchis)); 1395f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&patchscat)); 140c4762a1bSJed Brown } 141c4762a1bSJed Brown 142c4762a1bSJed Brown /* view the various parts */ 143c4762a1bSJed Brown { 144c4762a1bSJed Brown for (i = 0; i < size; i++) { 145c4762a1bSJed Brown if (i == rank) { 1465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Processor %d: \n",i)); 1475f80ce2aSJacob Faibussowitsch CHKERRQ(DMView(subda[0],PETSC_VIEWER_STDOUT_SELF)); 148c4762a1bSJed Brown } 1495f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Barrier(PETSC_COMM_WORLD)); 150c4762a1bSJed Brown } 151c4762a1bSJed Brown 1525f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(subda[0],&slvec)); 1535f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(subda[0],&sgvec)); 1545f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(da,&v)); 155c4762a1bSJed Brown 156c4762a1bSJed Brown /* test filling outer between the big DM and the small ones with the IS scatter*/ 1575f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterCreate(v,ois[0],sgvec,NULL,&oscata)); 158c4762a1bSJed Brown 1595f80ce2aSJacob Faibussowitsch CHKERRQ(FillLocalSubdomain(subda[0],sgvec)); 160c4762a1bSJed Brown 1615f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(oscata,sgvec,v,ADD_VALUES,SCATTER_REVERSE)); 1625f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(oscata,sgvec,v,ADD_VALUES,SCATTER_REVERSE)); 163c4762a1bSJed Brown 164c4762a1bSJed Brown /* test the local-to-local scatter */ 165c4762a1bSJed Brown 166c4762a1bSJed Brown /* fill up the local subdomain and then add them together */ 1675f80ce2aSJacob Faibussowitsch CHKERRQ(FillLocalSubdomain(da,v)); 168c4762a1bSJed Brown 1695f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(gscat[0],v,slvec,ADD_VALUES,SCATTER_FORWARD)); 1705f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(gscat[0],v,slvec,ADD_VALUES,SCATTER_FORWARD)); 171c4762a1bSJed Brown 1725f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(v,PETSC_VIEWER_STDOUT_WORLD)); 173c4762a1bSJed Brown 174c4762a1bSJed Brown /* test ghost scattering backwards */ 175c4762a1bSJed Brown 1765f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(v,0)); 177c4762a1bSJed Brown 1785f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(gscat[0],slvec,v,ADD_VALUES,SCATTER_REVERSE)); 1795f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(gscat[0],slvec,v,ADD_VALUES,SCATTER_REVERSE)); 180c4762a1bSJed Brown 1815f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(v,PETSC_VIEWER_STDOUT_WORLD)); 182c4762a1bSJed Brown 183c4762a1bSJed Brown /* test overlap scattering backwards */ 184c4762a1bSJed Brown 1855f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToGlobalBegin(subda[0],slvec,ADD_VALUES,sgvec)); 1865f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToGlobalEnd(subda[0],slvec,ADD_VALUES,sgvec)); 187c4762a1bSJed Brown 1885f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(v,0)); 189c4762a1bSJed Brown 1905f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(oscat[0],sgvec,v,ADD_VALUES,SCATTER_REVERSE)); 1915f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(oscat[0],sgvec,v,ADD_VALUES,SCATTER_REVERSE)); 192c4762a1bSJed Brown 1935f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(v,PETSC_VIEWER_STDOUT_WORLD)); 194c4762a1bSJed Brown 195c4762a1bSJed Brown /* test interior scattering backwards */ 196c4762a1bSJed Brown 1975f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(v,0)); 198c4762a1bSJed Brown 1995f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(iscat[0],sgvec,v,ADD_VALUES,SCATTER_REVERSE)); 2005f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(iscat[0],sgvec,v,ADD_VALUES,SCATTER_REVERSE)); 201c4762a1bSJed Brown 2025f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(v,PETSC_VIEWER_STDOUT_WORLD)); 203c4762a1bSJed Brown 204c4762a1bSJed Brown /* test matrix allocation */ 205c4762a1bSJed Brown for (i = 0; i < size; i++) { 206c4762a1bSJed Brown if (i == rank) { 207c4762a1bSJed Brown Mat m; 2085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Processor %d: \n",i)); 2095f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetMatType(subda[0],MATAIJ)); 2105f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(subda[0],&m)); 2115f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(m,PETSC_VIEWER_STDOUT_SELF)); 2125f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&m)); 213c4762a1bSJed Brown } 2145f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Barrier(PETSC_COMM_WORLD)); 215c4762a1bSJed Brown } 2165f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(subda[0],&slvec)); 2175f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(subda[0],&sgvec)); 2185f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(da,&v)); 219c4762a1bSJed Brown } 220c4762a1bSJed Brown 2215f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&subda[0])); 2225f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&ois[0])); 2235f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iis[0])); 224c4762a1bSJed Brown 2255f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&iscat[0])); 2265f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&oscat[0])); 2275f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&gscat[0])); 2285f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&oscata)); 229c4762a1bSJed Brown 2305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(iscat)); 2315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(oscat)); 2325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(gscat)); 2335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(oscata)); 234c4762a1bSJed Brown 2355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(subda)); 2365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ois)); 2375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(iis)); 238c4762a1bSJed Brown 2395f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&da)); 240*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 241*b122ec5aSJacob Faibussowitsch return 0; 242c4762a1bSJed Brown } 243