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