1 2 static char help[] = "Tests DMCreateDomainDecomposition.\n\n"; 3 4 /* 5 Use the options 6 -da_grid_x <nx> - number of grid points in x direction, if M < 0 7 -da_grid_y <ny> - number of grid points in y direction, if N < 0 8 -da_processors_x <MX> number of processors in x directio 9 -da_processors_y <MY> number of processors in x direction 10 */ 11 12 #include <petscdm.h> 13 #include <petscdmda.h> 14 15 PetscErrorCode FillLocalSubdomain(DM da, Vec gvec) 16 { 17 DMDALocalInfo info; 18 PetscMPIInt rank; 19 PetscInt i,j,k,l; 20 PetscErrorCode ierr; 21 22 PetscFunctionBeginUser; 23 CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 24 CHKERRQ(DMDAGetLocalInfo(da,&info)); 25 26 if (info.dim == 3) { 27 PetscScalar ***g; 28 CHKERRQ(DMDAVecGetArray(da,gvec,&g)); 29 /* loop over ghosts */ 30 for (k=info.zs; k<info.zs+info.zm; k++) { 31 for (j=info.ys; j<info.ys+info.ym; j++) { 32 for (i=info.xs; i<info.xs+info.xm; i++) { 33 g[k][j][info.dof*i+0] = i; 34 g[k][j][info.dof*i+1] = j; 35 g[k][j][info.dof*i+2] = k; 36 } 37 } 38 } 39 CHKERRQ(DMDAVecRestoreArray(da,gvec,&g)); 40 } 41 if (info.dim == 2) { 42 PetscScalar **g; 43 CHKERRQ(DMDAVecGetArray(da,gvec,&g)); 44 /* loop over ghosts */ 45 for (j=info.ys; j<info.ys+info.ym; j++) { 46 for (i=info.xs; i<info.xs+info.xm; i++) { 47 for (l = 0;l<info.dof;l++) { 48 g[j][info.dof*i+0] = i; 49 g[j][info.dof*i+1] = j; 50 g[j][info.dof*i+2] = rank; 51 } 52 } 53 } 54 CHKERRQ(DMDAVecRestoreArray(da,gvec,&g)); 55 } 56 PetscFunctionReturn(0); 57 } 58 59 int main(int argc,char **argv) 60 { 61 PetscErrorCode ierr; 62 DM da,*subda; 63 PetscInt i,dim = 3; 64 PetscInt M = 25, N = 25, P = 25; 65 PetscMPIInt size,rank; 66 Vec v; 67 Vec slvec,sgvec; 68 IS *ois,*iis; 69 VecScatter oscata; 70 VecScatter *iscat,*oscat,*gscat; 71 DMDALocalInfo info; 72 PetscBool patchis_offproc = PETSC_TRUE; 73 74 CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 75 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL)); 76 77 /* Create distributed array and get vectors */ 78 CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 79 CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 80 if (dim == 2) { 81 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)); 82 } else if (dim == 3) { 83 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)); 84 } 85 CHKERRQ(DMSetFromOptions(da)); 86 CHKERRQ(DMSetUp(da)); 87 CHKERRQ(DMDAGetLocalInfo(da,&info)); 88 89 CHKERRQ(DMCreateDomainDecomposition(da,NULL,NULL,&iis,&ois,&subda)); 90 CHKERRQ(DMCreateDomainDecompositionScatters(da,1,subda,&iscat,&oscat,&gscat)); 91 92 { 93 DMDALocalInfo subinfo; 94 MatStencil lower,upper; 95 IS patchis; 96 Vec smallvec; 97 Vec largevec; 98 VecScatter patchscat; 99 100 CHKERRQ(DMDAGetLocalInfo(subda[0],&subinfo)); 101 102 lower.i = info.xs; 103 lower.j = info.ys; 104 lower.k = info.zs; 105 upper.i = info.xs+info.xm; 106 upper.j = info.ys+info.ym; 107 upper.k = info.zs+info.zm; 108 109 /* test the patch IS as a thing to scatter to/from */ 110 CHKERRQ(DMDACreatePatchIS(da,&lower,&upper,&patchis,patchis_offproc)); 111 CHKERRQ(DMGetGlobalVector(da,&largevec)); 112 113 CHKERRQ(VecCreate(PETSC_COMM_SELF,&smallvec)); 114 CHKERRQ(VecSetSizes(smallvec,info.dof*(upper.i - lower.i)*(upper.j - lower.j)*(upper.k - lower.k),PETSC_DECIDE)); 115 CHKERRQ(VecSetFromOptions(smallvec)); 116 CHKERRQ(VecScatterCreate(smallvec,NULL,largevec,patchis,&patchscat)); 117 118 CHKERRQ(FillLocalSubdomain(subda[0],smallvec)); 119 CHKERRQ(VecSet(largevec,0)); 120 121 CHKERRQ(VecScatterBegin(patchscat,smallvec,largevec,ADD_VALUES,SCATTER_FORWARD)); 122 CHKERRQ(VecScatterEnd(patchscat,smallvec,largevec,ADD_VALUES,SCATTER_FORWARD)); 123 CHKERRQ(ISView(patchis,PETSC_VIEWER_STDOUT_WORLD)); 124 CHKERRQ(VecScatterView(patchscat,PETSC_VIEWER_STDOUT_WORLD)); 125 126 for (i = 0; i < size; i++) { 127 if (i == rank) { 128 CHKERRQ(VecView(smallvec,PETSC_VIEWER_STDOUT_SELF)); 129 } 130 CHKERRMPI(MPI_Barrier(PETSC_COMM_WORLD)); 131 } 132 133 CHKERRMPI(MPI_Barrier(PETSC_COMM_WORLD)); 134 CHKERRQ(VecView(largevec,PETSC_VIEWER_STDOUT_WORLD)); 135 136 CHKERRQ(VecDestroy(&smallvec)); 137 CHKERRQ(DMRestoreGlobalVector(da,&largevec)); 138 CHKERRQ(ISDestroy(&patchis)); 139 CHKERRQ(VecScatterDestroy(&patchscat)); 140 } 141 142 /* view the various parts */ 143 { 144 for (i = 0; i < size; i++) { 145 if (i == rank) { 146 CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Processor %d: \n",i)); 147 CHKERRQ(DMView(subda[0],PETSC_VIEWER_STDOUT_SELF)); 148 } 149 CHKERRMPI(MPI_Barrier(PETSC_COMM_WORLD)); 150 } 151 152 CHKERRQ(DMGetLocalVector(subda[0],&slvec)); 153 CHKERRQ(DMGetGlobalVector(subda[0],&sgvec)); 154 CHKERRQ(DMGetGlobalVector(da,&v)); 155 156 /* test filling outer between the big DM and the small ones with the IS scatter*/ 157 CHKERRQ(VecScatterCreate(v,ois[0],sgvec,NULL,&oscata)); 158 159 CHKERRQ(FillLocalSubdomain(subda[0],sgvec)); 160 161 CHKERRQ(VecScatterBegin(oscata,sgvec,v,ADD_VALUES,SCATTER_REVERSE)); 162 CHKERRQ(VecScatterEnd(oscata,sgvec,v,ADD_VALUES,SCATTER_REVERSE)); 163 164 /* test the local-to-local scatter */ 165 166 /* fill up the local subdomain and then add them together */ 167 CHKERRQ(FillLocalSubdomain(da,v)); 168 169 CHKERRQ(VecScatterBegin(gscat[0],v,slvec,ADD_VALUES,SCATTER_FORWARD)); 170 CHKERRQ(VecScatterEnd(gscat[0],v,slvec,ADD_VALUES,SCATTER_FORWARD)); 171 172 CHKERRQ(VecView(v,PETSC_VIEWER_STDOUT_WORLD)); 173 174 /* test ghost scattering backwards */ 175 176 CHKERRQ(VecSet(v,0)); 177 178 CHKERRQ(VecScatterBegin(gscat[0],slvec,v,ADD_VALUES,SCATTER_REVERSE)); 179 CHKERRQ(VecScatterEnd(gscat[0],slvec,v,ADD_VALUES,SCATTER_REVERSE)); 180 181 CHKERRQ(VecView(v,PETSC_VIEWER_STDOUT_WORLD)); 182 183 /* test overlap scattering backwards */ 184 185 CHKERRQ(DMLocalToGlobalBegin(subda[0],slvec,ADD_VALUES,sgvec)); 186 CHKERRQ(DMLocalToGlobalEnd(subda[0],slvec,ADD_VALUES,sgvec)); 187 188 CHKERRQ(VecSet(v,0)); 189 190 CHKERRQ(VecScatterBegin(oscat[0],sgvec,v,ADD_VALUES,SCATTER_REVERSE)); 191 CHKERRQ(VecScatterEnd(oscat[0],sgvec,v,ADD_VALUES,SCATTER_REVERSE)); 192 193 CHKERRQ(VecView(v,PETSC_VIEWER_STDOUT_WORLD)); 194 195 /* test interior scattering backwards */ 196 197 CHKERRQ(VecSet(v,0)); 198 199 CHKERRQ(VecScatterBegin(iscat[0],sgvec,v,ADD_VALUES,SCATTER_REVERSE)); 200 CHKERRQ(VecScatterEnd(iscat[0],sgvec,v,ADD_VALUES,SCATTER_REVERSE)); 201 202 CHKERRQ(VecView(v,PETSC_VIEWER_STDOUT_WORLD)); 203 204 /* test matrix allocation */ 205 for (i = 0; i < size; i++) { 206 if (i == rank) { 207 Mat m; 208 CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Processor %d: \n",i)); 209 CHKERRQ(DMSetMatType(subda[0],MATAIJ)); 210 CHKERRQ(DMCreateMatrix(subda[0],&m)); 211 CHKERRQ(MatView(m,PETSC_VIEWER_STDOUT_SELF)); 212 CHKERRQ(MatDestroy(&m)); 213 } 214 CHKERRMPI(MPI_Barrier(PETSC_COMM_WORLD)); 215 } 216 CHKERRQ(DMRestoreLocalVector(subda[0],&slvec)); 217 CHKERRQ(DMRestoreGlobalVector(subda[0],&sgvec)); 218 CHKERRQ(DMRestoreGlobalVector(da,&v)); 219 } 220 221 CHKERRQ(DMDestroy(&subda[0])); 222 CHKERRQ(ISDestroy(&ois[0])); 223 CHKERRQ(ISDestroy(&iis[0])); 224 225 CHKERRQ(VecScatterDestroy(&iscat[0])); 226 CHKERRQ(VecScatterDestroy(&oscat[0])); 227 CHKERRQ(VecScatterDestroy(&gscat[0])); 228 CHKERRQ(VecScatterDestroy(&oscata)); 229 230 CHKERRQ(PetscFree(iscat)); 231 CHKERRQ(PetscFree(oscat)); 232 CHKERRQ(PetscFree(gscat)); 233 CHKERRQ(PetscFree(oscata)); 234 235 CHKERRQ(PetscFree(subda)); 236 CHKERRQ(PetscFree(ois)); 237 CHKERRQ(PetscFree(iis)); 238 239 CHKERRQ(DMDestroy(&da)); 240 CHKERRQ(PetscFinalize()); 241 return 0; 242 } 243