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; 23*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 24*9566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da,&info)); 25c4762a1bSJed Brown 26c4762a1bSJed Brown if (info.dim == 3) { 27c4762a1bSJed Brown PetscScalar ***g; 28*9566063dSJacob Faibussowitsch PetscCall(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 } 39*9566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,gvec,&g)); 40c4762a1bSJed Brown } 41c4762a1bSJed Brown if (info.dim == 2) { 42c4762a1bSJed Brown PetscScalar **g; 43*9566063dSJacob Faibussowitsch PetscCall(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 } 54*9566063dSJacob Faibussowitsch PetscCall(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*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 75*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL)); 76c4762a1bSJed Brown 77c4762a1bSJed Brown /* Create distributed array and get vectors */ 78*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 79*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 80c4762a1bSJed Brown if (dim == 2) { 81*9566063dSJacob Faibussowitsch PetscCall(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) { 83*9566063dSJacob Faibussowitsch PetscCall(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 } 85*9566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 86*9566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 87*9566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da,&info)); 88c4762a1bSJed Brown 89*9566063dSJacob Faibussowitsch PetscCall(DMCreateDomainDecomposition(da,NULL,NULL,&iis,&ois,&subda)); 90*9566063dSJacob Faibussowitsch PetscCall(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 100*9566063dSJacob Faibussowitsch PetscCall(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 */ 110*9566063dSJacob Faibussowitsch PetscCall(DMDACreatePatchIS(da,&lower,&upper,&patchis,patchis_offproc)); 111*9566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da,&largevec)); 112c4762a1bSJed Brown 113*9566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF,&smallvec)); 114*9566063dSJacob Faibussowitsch PetscCall(VecSetSizes(smallvec,info.dof*(upper.i - lower.i)*(upper.j - lower.j)*(upper.k - lower.k),PETSC_DECIDE)); 115*9566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(smallvec)); 116*9566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(smallvec,NULL,largevec,patchis,&patchscat)); 117c4762a1bSJed Brown 118*9566063dSJacob Faibussowitsch PetscCall(FillLocalSubdomain(subda[0],smallvec)); 119*9566063dSJacob Faibussowitsch PetscCall(VecSet(largevec,0)); 120c4762a1bSJed Brown 121*9566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(patchscat,smallvec,largevec,ADD_VALUES,SCATTER_FORWARD)); 122*9566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(patchscat,smallvec,largevec,ADD_VALUES,SCATTER_FORWARD)); 123*9566063dSJacob Faibussowitsch PetscCall(ISView(patchis,PETSC_VIEWER_STDOUT_WORLD)); 124*9566063dSJacob Faibussowitsch PetscCall(VecScatterView(patchscat,PETSC_VIEWER_STDOUT_WORLD)); 125c4762a1bSJed Brown 126c4762a1bSJed Brown for (i = 0; i < size; i++) { 127c4762a1bSJed Brown if (i == rank) { 128*9566063dSJacob Faibussowitsch PetscCall(VecView(smallvec,PETSC_VIEWER_STDOUT_SELF)); 129c4762a1bSJed Brown } 130*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD)); 131c4762a1bSJed Brown } 132c4762a1bSJed Brown 133*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD)); 134*9566063dSJacob Faibussowitsch PetscCall(VecView(largevec,PETSC_VIEWER_STDOUT_WORLD)); 135c4762a1bSJed Brown 136*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&smallvec)); 137*9566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(da,&largevec)); 138*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&patchis)); 139*9566063dSJacob Faibussowitsch PetscCall(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) { 146*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Processor %d: \n",i)); 147*9566063dSJacob Faibussowitsch PetscCall(DMView(subda[0],PETSC_VIEWER_STDOUT_SELF)); 148c4762a1bSJed Brown } 149*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD)); 150c4762a1bSJed Brown } 151c4762a1bSJed Brown 152*9566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(subda[0],&slvec)); 153*9566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(subda[0],&sgvec)); 154*9566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da,&v)); 155c4762a1bSJed Brown 156c4762a1bSJed Brown /* test filling outer between the big DM and the small ones with the IS scatter*/ 157*9566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(v,ois[0],sgvec,NULL,&oscata)); 158c4762a1bSJed Brown 159*9566063dSJacob Faibussowitsch PetscCall(FillLocalSubdomain(subda[0],sgvec)); 160c4762a1bSJed Brown 161*9566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(oscata,sgvec,v,ADD_VALUES,SCATTER_REVERSE)); 162*9566063dSJacob Faibussowitsch PetscCall(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 */ 167*9566063dSJacob Faibussowitsch PetscCall(FillLocalSubdomain(da,v)); 168c4762a1bSJed Brown 169*9566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat[0],v,slvec,ADD_VALUES,SCATTER_FORWARD)); 170*9566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat[0],v,slvec,ADD_VALUES,SCATTER_FORWARD)); 171c4762a1bSJed Brown 172*9566063dSJacob Faibussowitsch PetscCall(VecView(v,PETSC_VIEWER_STDOUT_WORLD)); 173c4762a1bSJed Brown 174c4762a1bSJed Brown /* test ghost scattering backwards */ 175c4762a1bSJed Brown 176*9566063dSJacob Faibussowitsch PetscCall(VecSet(v,0)); 177c4762a1bSJed Brown 178*9566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat[0],slvec,v,ADD_VALUES,SCATTER_REVERSE)); 179*9566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat[0],slvec,v,ADD_VALUES,SCATTER_REVERSE)); 180c4762a1bSJed Brown 181*9566063dSJacob Faibussowitsch PetscCall(VecView(v,PETSC_VIEWER_STDOUT_WORLD)); 182c4762a1bSJed Brown 183c4762a1bSJed Brown /* test overlap scattering backwards */ 184c4762a1bSJed Brown 185*9566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(subda[0],slvec,ADD_VALUES,sgvec)); 186*9566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(subda[0],slvec,ADD_VALUES,sgvec)); 187c4762a1bSJed Brown 188*9566063dSJacob Faibussowitsch PetscCall(VecSet(v,0)); 189c4762a1bSJed Brown 190*9566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(oscat[0],sgvec,v,ADD_VALUES,SCATTER_REVERSE)); 191*9566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(oscat[0],sgvec,v,ADD_VALUES,SCATTER_REVERSE)); 192c4762a1bSJed Brown 193*9566063dSJacob Faibussowitsch PetscCall(VecView(v,PETSC_VIEWER_STDOUT_WORLD)); 194c4762a1bSJed Brown 195c4762a1bSJed Brown /* test interior scattering backwards */ 196c4762a1bSJed Brown 197*9566063dSJacob Faibussowitsch PetscCall(VecSet(v,0)); 198c4762a1bSJed Brown 199*9566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(iscat[0],sgvec,v,ADD_VALUES,SCATTER_REVERSE)); 200*9566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(iscat[0],sgvec,v,ADD_VALUES,SCATTER_REVERSE)); 201c4762a1bSJed Brown 202*9566063dSJacob Faibussowitsch PetscCall(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; 208*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Processor %d: \n",i)); 209*9566063dSJacob Faibussowitsch PetscCall(DMSetMatType(subda[0],MATAIJ)); 210*9566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(subda[0],&m)); 211*9566063dSJacob Faibussowitsch PetscCall(MatView(m,PETSC_VIEWER_STDOUT_SELF)); 212*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&m)); 213c4762a1bSJed Brown } 214*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD)); 215c4762a1bSJed Brown } 216*9566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(subda[0],&slvec)); 217*9566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(subda[0],&sgvec)); 218*9566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(da,&v)); 219c4762a1bSJed Brown } 220c4762a1bSJed Brown 221*9566063dSJacob Faibussowitsch PetscCall(DMDestroy(&subda[0])); 222*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ois[0])); 223*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iis[0])); 224c4762a1bSJed Brown 225*9566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&iscat[0])); 226*9566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&oscat[0])); 227*9566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&gscat[0])); 228*9566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&oscata)); 229c4762a1bSJed Brown 230*9566063dSJacob Faibussowitsch PetscCall(PetscFree(iscat)); 231*9566063dSJacob Faibussowitsch PetscCall(PetscFree(oscat)); 232*9566063dSJacob Faibussowitsch PetscCall(PetscFree(gscat)); 233*9566063dSJacob Faibussowitsch PetscCall(PetscFree(oscata)); 234c4762a1bSJed Brown 235*9566063dSJacob Faibussowitsch PetscCall(PetscFree(subda)); 236*9566063dSJacob Faibussowitsch PetscCall(PetscFree(ois)); 237*9566063dSJacob Faibussowitsch PetscCall(PetscFree(iis)); 238c4762a1bSJed Brown 239*9566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 240*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 241b122ec5aSJacob Faibussowitsch return 0; 242c4762a1bSJed Brown } 243