1 2 #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 3 #include <petscdraw.h> 4 5 static PetscErrorCode DMView_DA_2d(DM da,PetscViewer viewer) 6 { 7 PetscMPIInt rank; 8 PetscBool iascii,isdraw,isglvis,isbinary; 9 DM_DA *dd = (DM_DA*)da->data; 10 #if defined(PETSC_HAVE_MATLAB_ENGINE) 11 PetscBool ismatlab; 12 #endif 13 14 PetscFunctionBegin; 15 CHKERRMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank)); 16 17 CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 18 CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw)); 19 CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis)); 20 CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary)); 21 #if defined(PETSC_HAVE_MATLAB_ENGINE) 22 CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab)); 23 #endif 24 if (iascii) { 25 PetscViewerFormat format; 26 27 CHKERRQ(PetscViewerGetFormat(viewer, &format)); 28 if (format == PETSC_VIEWER_LOAD_BALANCE) { 29 PetscInt i,nmax = 0,nmin = PETSC_MAX_INT,navg = 0,*nz,nzlocal; 30 DMDALocalInfo info; 31 PetscMPIInt size; 32 CHKERRMPI(MPI_Comm_size(PetscObjectComm((PetscObject)da),&size)); 33 CHKERRQ(DMDAGetLocalInfo(da,&info)); 34 nzlocal = info.xm*info.ym; 35 CHKERRQ(PetscMalloc1(size,&nz)); 36 CHKERRMPI(MPI_Allgather(&nzlocal,1,MPIU_INT,nz,1,MPIU_INT,PetscObjectComm((PetscObject)da))); 37 for (i=0; i<(PetscInt)size; i++) { 38 nmax = PetscMax(nmax,nz[i]); 39 nmin = PetscMin(nmin,nz[i]); 40 navg += nz[i]; 41 } 42 CHKERRQ(PetscFree(nz)); 43 navg = navg/size; 44 CHKERRQ(PetscViewerASCIIPrintf(viewer," Load Balance - Grid Points: Min %D avg %D max %D\n",nmin,navg,nmax)); 45 PetscFunctionReturn(0); 46 } 47 if (format != PETSC_VIEWER_ASCII_VTK_DEPRECATED && format != PETSC_VIEWER_ASCII_VTK_CELL_DEPRECATED && format != PETSC_VIEWER_ASCII_GLVIS) { 48 DMDALocalInfo info; 49 CHKERRQ(DMDAGetLocalInfo(da,&info)); 50 CHKERRQ(PetscViewerASCIIPushSynchronized(viewer)); 51 CHKERRQ(PetscViewerASCIISynchronizedPrintf(viewer,"Processor [%d] M %D N %D m %D n %D w %D s %D\n",rank,dd->M,dd->N,dd->m,dd->n,dd->w,dd->s)); 52 CHKERRQ(PetscViewerASCIISynchronizedPrintf(viewer,"X range of indices: %D %D, Y range of indices: %D %D\n",info.xs,info.xs+info.xm,info.ys,info.ys+info.ym)); 53 CHKERRQ(PetscViewerFlush(viewer)); 54 CHKERRQ(PetscViewerASCIIPopSynchronized(viewer)); 55 } else if (format == PETSC_VIEWER_ASCII_GLVIS) { 56 CHKERRQ(DMView_DA_GLVis(da,viewer)); 57 } else { 58 CHKERRQ(DMView_DA_VTK(da,viewer)); 59 } 60 } else if (isdraw) { 61 PetscDraw draw; 62 double ymin = -1*dd->s-1,ymax = dd->N+dd->s; 63 double xmin = -1*dd->s-1,xmax = dd->M+dd->s; 64 double x,y; 65 PetscInt base; 66 const PetscInt *idx; 67 char node[10]; 68 PetscBool isnull; 69 PetscErrorCode ierr; 70 71 CHKERRQ(PetscViewerDrawGetDraw(viewer,0,&draw)); 72 CHKERRQ(PetscDrawIsNull(draw,&isnull)); 73 if (isnull) PetscFunctionReturn(0); 74 75 CHKERRQ(PetscDrawCheckResizedWindow(draw)); 76 CHKERRQ(PetscDrawClear(draw)); 77 CHKERRQ(PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax)); 78 79 ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr); 80 /* first processor draw all node lines */ 81 if (rank == 0) { 82 ymin = 0.0; ymax = dd->N - 1; 83 for (xmin=0; xmin<dd->M; xmin++) { 84 CHKERRQ(PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK)); 85 } 86 xmin = 0.0; xmax = dd->M - 1; 87 for (ymin=0; ymin<dd->N; ymin++) { 88 CHKERRQ(PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK)); 89 } 90 } 91 ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr); 92 CHKERRQ(PetscDrawFlush(draw)); 93 CHKERRQ(PetscDrawPause(draw)); 94 95 ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr); 96 /* draw my box */ 97 xmin = dd->xs/dd->w; xmax =(dd->xe-1)/dd->w; ymin = dd->ys; ymax = dd->ye - 1; 98 CHKERRQ(PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED)); 99 CHKERRQ(PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED)); 100 CHKERRQ(PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED)); 101 CHKERRQ(PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED)); 102 /* put in numbers */ 103 base = (dd->base)/dd->w; 104 for (y=ymin; y<=ymax; y++) { 105 for (x=xmin; x<=xmax; x++) { 106 CHKERRQ(PetscSNPrintf(node,sizeof(node),"%d",(int)base++)); 107 CHKERRQ(PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node)); 108 } 109 } 110 ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr); 111 CHKERRQ(PetscDrawFlush(draw)); 112 CHKERRQ(PetscDrawPause(draw)); 113 114 ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr); 115 /* overlay ghost numbers, useful for error checking */ 116 CHKERRQ(ISLocalToGlobalMappingGetBlockIndices(da->ltogmap,&idx)); 117 base = 0; xmin = dd->Xs; xmax = dd->Xe; ymin = dd->Ys; ymax = dd->Ye; 118 for (y=ymin; y<ymax; y++) { 119 for (x=xmin; x<xmax; x++) { 120 if ((base % dd->w) == 0) { 121 CHKERRQ(PetscSNPrintf(node,sizeof(node),"%d",(int)(idx[base/dd->w]))); 122 CHKERRQ(PetscDrawString(draw,x/dd->w,y,PETSC_DRAW_BLUE,node)); 123 } 124 base++; 125 } 126 } 127 CHKERRQ(ISLocalToGlobalMappingRestoreBlockIndices(da->ltogmap,&idx)); 128 ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr); 129 CHKERRQ(PetscDrawFlush(draw)); 130 CHKERRQ(PetscDrawPause(draw)); 131 CHKERRQ(PetscDrawSave(draw)); 132 } else if (isglvis) { 133 CHKERRQ(DMView_DA_GLVis(da,viewer)); 134 } else if (isbinary) { 135 CHKERRQ(DMView_DA_Binary(da,viewer)); 136 #if defined(PETSC_HAVE_MATLAB_ENGINE) 137 } else if (ismatlab) { 138 CHKERRQ(DMView_DA_Matlab(da,viewer)); 139 #endif 140 } 141 PetscFunctionReturn(0); 142 } 143 144 #if defined(new) 145 /* 146 DMDAGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix where local 147 function lives on a DMDA 148 149 y ~= (F(u + ha) - F(u))/h, 150 where F = nonlinear function, as set by SNESSetFunction() 151 u = current iterate 152 h = difference interval 153 */ 154 PetscErrorCode DMDAGetDiagonal_MFFD(DM da,Vec U,Vec a) 155 { 156 PetscScalar h,*aa,*ww,v; 157 PetscReal epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 158 PetscInt gI,nI; 159 MatStencil stencil; 160 DMDALocalInfo info; 161 162 PetscFunctionBegin; 163 CHKERRQ((*ctx->func)(0,U,a,ctx->funcctx)); 164 CHKERRQ((*ctx->funcisetbase)(U,ctx->funcctx)); 165 166 CHKERRQ(VecGetArray(U,&ww)); 167 CHKERRQ(VecGetArray(a,&aa)); 168 169 nI = 0; 170 h = ww[gI]; 171 if (h == 0.0) h = 1.0; 172 if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin; 173 else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin; 174 h *= epsilon; 175 176 ww[gI] += h; 177 CHKERRQ((*ctx->funci)(i,w,&v,ctx->funcctx)); 178 aa[nI] = (v - aa[nI])/h; 179 ww[gI] -= h; 180 nI++; 181 182 CHKERRQ(VecRestoreArray(U,&ww)); 183 CHKERRQ(VecRestoreArray(a,&aa)); 184 PetscFunctionReturn(0); 185 } 186 #endif 187 188 PetscErrorCode DMSetUp_DA_2D(DM da) 189 { 190 DM_DA *dd = (DM_DA*)da->data; 191 const PetscInt M = dd->M; 192 const PetscInt N = dd->N; 193 PetscInt m = dd->m; 194 PetscInt n = dd->n; 195 const PetscInt dof = dd->w; 196 const PetscInt s = dd->s; 197 DMBoundaryType bx = dd->bx; 198 DMBoundaryType by = dd->by; 199 DMDAStencilType stencil_type = dd->stencil_type; 200 PetscInt *lx = dd->lx; 201 PetscInt *ly = dd->ly; 202 MPI_Comm comm; 203 PetscMPIInt rank,size; 204 PetscInt xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,IXs,IXe,IYs,IYe; 205 PetscInt up,down,left,right,i,n0,n1,n2,n3,n5,n6,n7,n8,*idx,nn; 206 PetscInt xbase,*bases,*ldims,j,x_t,y_t,s_t,base,count; 207 PetscInt s_x,s_y; /* s proportionalized to w */ 208 PetscInt sn0 = 0,sn2 = 0,sn6 = 0,sn8 = 0; 209 Vec local,global; 210 VecScatter gtol; 211 IS to,from; 212 213 PetscFunctionBegin; 214 PetscCheck(stencil_type != DMDA_STENCIL_BOX || (bx != DM_BOUNDARY_MIRROR && by != DM_BOUNDARY_MIRROR),PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Mirror boundary and box stencil"); 215 CHKERRQ(PetscObjectGetComm((PetscObject)da,&comm)); 216 #if !defined(PETSC_USE_64BIT_INDICES) 217 PetscCheck(((PetscInt64) M)*((PetscInt64) N)*((PetscInt64) dof) <= (PetscInt64) PETSC_MPI_INT_MAX,comm,PETSC_ERR_INT_OVERFLOW,"Mesh of %D by %D by %D (dof) is too large for 32 bit indices",M,N,dof); 218 #endif 219 220 CHKERRMPI(MPI_Comm_size(comm,&size)); 221 CHKERRMPI(MPI_Comm_rank(comm,&rank)); 222 223 dd->p = 1; 224 if (m != PETSC_DECIDE) { 225 PetscCheck(m >= 1,comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m); 226 else PetscCheck(m <= size,comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size); 227 } 228 if (n != PETSC_DECIDE) { 229 PetscCheck(n >= 1,comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n); 230 else PetscCheck(n <= size,comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size); 231 } 232 233 if (m == PETSC_DECIDE || n == PETSC_DECIDE) { 234 if (n != PETSC_DECIDE) { 235 m = size/n; 236 } else if (m != PETSC_DECIDE) { 237 n = size/m; 238 } else { 239 /* try for squarish distribution */ 240 m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)N))); 241 if (!m) m = 1; 242 while (m > 0) { 243 n = size/m; 244 if (m*n == size) break; 245 m--; 246 } 247 if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;} 248 } 249 PetscCheck(m*n == size,comm,PETSC_ERR_PLIB,"Unable to create partition, check the size of the communicator and input m and n "); 250 } else PetscCheck(m*n == size,comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition"); 251 252 PetscCheck(M >= m,comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m); 253 PetscCheck(N >= n,comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n); 254 255 /* 256 Determine locally owned region 257 xs is the first local node number, x is the number of local nodes 258 */ 259 if (!lx) { 260 CHKERRQ(PetscMalloc1(m, &dd->lx)); 261 lx = dd->lx; 262 for (i=0; i<m; i++) { 263 lx[i] = M/m + ((M % m) > i); 264 } 265 } 266 x = lx[rank % m]; 267 xs = 0; 268 for (i=0; i<(rank % m); i++) { 269 xs += lx[i]; 270 } 271 if (PetscDefined(USE_DEBUG)) { 272 left = xs; 273 for (i=(rank % m); i<m; i++) { 274 left += lx[i]; 275 } 276 PetscCheck(left == M,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of lx across processors not equal to M: %D %D",left,M); 277 } 278 279 /* 280 Determine locally owned region 281 ys is the first local node number, y is the number of local nodes 282 */ 283 if (!ly) { 284 CHKERRQ(PetscMalloc1(n, &dd->ly)); 285 ly = dd->ly; 286 for (i=0; i<n; i++) { 287 ly[i] = N/n + ((N % n) > i); 288 } 289 } 290 y = ly[rank/m]; 291 ys = 0; 292 for (i=0; i<(rank/m); i++) { 293 ys += ly[i]; 294 } 295 if (PetscDefined(USE_DEBUG)) { 296 left = ys; 297 for (i=(rank/m); i<n; i++) { 298 left += ly[i]; 299 } 300 PetscCheck(left == N,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of ly across processors not equal to N: %D %D",left,N); 301 } 302 303 /* 304 check if the scatter requires more than one process neighbor or wraps around 305 the domain more than once 306 */ 307 PetscCheck((x >= s) || ((m <= 1) && (bx != DM_BOUNDARY_PERIODIC)),PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local x-width of domain x %D is smaller than stencil width s %D",x,s); 308 PetscCheck((y >= s) || ((n <= 1) && (by != DM_BOUNDARY_PERIODIC)),PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local y-width of domain y %D is smaller than stencil width s %D",y,s); 309 xe = xs + x; 310 ye = ys + y; 311 312 /* determine ghost region (Xs) and region scattered into (IXs) */ 313 if (xs-s > 0) { 314 Xs = xs - s; IXs = xs - s; 315 } else { 316 if (bx) { 317 Xs = xs - s; 318 } else { 319 Xs = 0; 320 } 321 IXs = 0; 322 } 323 if (xe+s <= M) { 324 Xe = xe + s; IXe = xe + s; 325 } else { 326 if (bx) { 327 Xs = xs - s; Xe = xe + s; 328 } else { 329 Xe = M; 330 } 331 IXe = M; 332 } 333 334 if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) { 335 IXs = xs - s; 336 IXe = xe + s; 337 Xs = xs - s; 338 Xe = xe + s; 339 } 340 341 if (ys-s > 0) { 342 Ys = ys - s; IYs = ys - s; 343 } else { 344 if (by) { 345 Ys = ys - s; 346 } else { 347 Ys = 0; 348 } 349 IYs = 0; 350 } 351 if (ye+s <= N) { 352 Ye = ye + s; IYe = ye + s; 353 } else { 354 if (by) { 355 Ye = ye + s; 356 } else { 357 Ye = N; 358 } 359 IYe = N; 360 } 361 362 if (by == DM_BOUNDARY_PERIODIC || by == DM_BOUNDARY_MIRROR) { 363 IYs = ys - s; 364 IYe = ye + s; 365 Ys = ys - s; 366 Ye = ye + s; 367 } 368 369 /* stencil length in each direction */ 370 s_x = s; 371 s_y = s; 372 373 /* determine starting point of each processor */ 374 nn = x*y; 375 CHKERRQ(PetscMalloc2(size+1,&bases,size,&ldims)); 376 CHKERRMPI(MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm)); 377 bases[0] = 0; 378 for (i=1; i<=size; i++) { 379 bases[i] = ldims[i-1]; 380 } 381 for (i=1; i<=size; i++) { 382 bases[i] += bases[i-1]; 383 } 384 base = bases[rank]*dof; 385 386 /* allocate the base parallel and sequential vectors */ 387 dd->Nlocal = x*y*dof; 388 CHKERRQ(VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,NULL,&global)); 389 dd->nlocal = (Xe-Xs)*(Ye-Ys)*dof; 390 CHKERRQ(VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,NULL,&local)); 391 392 /* generate global to local vector scatter and local to global mapping*/ 393 394 /* global to local must include ghost points within the domain, 395 but not ghost points outside the domain that aren't periodic */ 396 CHKERRQ(PetscMalloc1((IXe-IXs)*(IYe-IYs),&idx)); 397 if (stencil_type == DMDA_STENCIL_BOX) { 398 left = IXs - Xs; right = left + (IXe-IXs); 399 down = IYs - Ys; up = down + (IYe-IYs); 400 count = 0; 401 for (i=down; i<up; i++) { 402 for (j=left; j<right; j++) { 403 idx[count++] = j + i*(Xe-Xs); 404 } 405 } 406 CHKERRQ(ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to)); 407 408 } else { 409 /* must drop into cross shape region */ 410 /* ---------| 411 | top | 412 |--- ---| up 413 | middle | 414 | | 415 ---- ---- down 416 | bottom | 417 ----------- 418 Xs xs xe Xe */ 419 left = xs - Xs; right = left + x; 420 down = ys - Ys; up = down + y; 421 count = 0; 422 /* bottom */ 423 for (i=(IYs-Ys); i<down; i++) { 424 for (j=left; j<right; j++) { 425 idx[count++] = j + i*(Xe-Xs); 426 } 427 } 428 /* middle */ 429 for (i=down; i<up; i++) { 430 for (j=(IXs-Xs); j<(IXe-Xs); j++) { 431 idx[count++] = j + i*(Xe-Xs); 432 } 433 } 434 /* top */ 435 for (i=up; i<up+IYe-ye; i++) { 436 for (j=left; j<right; j++) { 437 idx[count++] = j + i*(Xe-Xs); 438 } 439 } 440 CHKERRQ(ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to)); 441 } 442 443 /* determine who lies on each side of us stored in n6 n7 n8 444 n3 n5 445 n0 n1 n2 446 */ 447 448 /* Assume the Non-Periodic Case */ 449 n1 = rank - m; 450 if (rank % m) { 451 n0 = n1 - 1; 452 } else { 453 n0 = -1; 454 } 455 if ((rank+1) % m) { 456 n2 = n1 + 1; 457 n5 = rank + 1; 458 n8 = rank + m + 1; if (n8 >= m*n) n8 = -1; 459 } else { 460 n2 = -1; n5 = -1; n8 = -1; 461 } 462 if (rank % m) { 463 n3 = rank - 1; 464 n6 = n3 + m; if (n6 >= m*n) n6 = -1; 465 } else { 466 n3 = -1; n6 = -1; 467 } 468 n7 = rank + m; if (n7 >= m*n) n7 = -1; 469 470 if (bx == DM_BOUNDARY_PERIODIC && by == DM_BOUNDARY_PERIODIC) { 471 /* Modify for Periodic Cases */ 472 /* Handle all four corners */ 473 if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1; 474 if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0; 475 if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m; 476 if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1; 477 478 /* Handle Top and Bottom Sides */ 479 if (n1 < 0) n1 = rank + m * (n-1); 480 if (n7 < 0) n7 = rank - m * (n-1); 481 if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1; 482 if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1; 483 if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1; 484 if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1; 485 486 /* Handle Left and Right Sides */ 487 if (n3 < 0) n3 = rank + (m-1); 488 if (n5 < 0) n5 = rank - (m-1); 489 if ((n1 >= 0) && (n0 < 0)) n0 = rank-1; 490 if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1; 491 if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1; 492 if ((n7 >= 0) && (n8 < 0)) n8 = rank+1; 493 } else if (by == DM_BOUNDARY_PERIODIC) { /* Handle Top and Bottom Sides */ 494 if (n1 < 0) n1 = rank + m * (n-1); 495 if (n7 < 0) n7 = rank - m * (n-1); 496 if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1; 497 if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1; 498 if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1; 499 if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1; 500 } else if (bx == DM_BOUNDARY_PERIODIC) { /* Handle Left and Right Sides */ 501 if (n3 < 0) n3 = rank + (m-1); 502 if (n5 < 0) n5 = rank - (m-1); 503 if ((n1 >= 0) && (n0 < 0)) n0 = rank-1; 504 if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1; 505 if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1; 506 if ((n7 >= 0) && (n8 < 0)) n8 = rank+1; 507 } 508 509 CHKERRQ(PetscMalloc1(9,&dd->neighbors)); 510 511 dd->neighbors[0] = n0; 512 dd->neighbors[1] = n1; 513 dd->neighbors[2] = n2; 514 dd->neighbors[3] = n3; 515 dd->neighbors[4] = rank; 516 dd->neighbors[5] = n5; 517 dd->neighbors[6] = n6; 518 dd->neighbors[7] = n7; 519 dd->neighbors[8] = n8; 520 521 if (stencil_type == DMDA_STENCIL_STAR) { 522 /* save corner processor numbers */ 523 sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8; 524 n0 = n2 = n6 = n8 = -1; 525 } 526 527 CHKERRQ(PetscMalloc1((Xe-Xs)*(Ye-Ys),&idx)); 528 529 nn = 0; 530 xbase = bases[rank]; 531 for (i=1; i<=s_y; i++) { 532 if (n0 >= 0) { /* left below */ 533 x_t = lx[n0 % m]; 534 y_t = ly[(n0/m)]; 535 s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x; 536 for (j=0; j<s_x; j++) idx[nn++] = s_t++; 537 } 538 539 if (n1 >= 0) { /* directly below */ 540 x_t = x; 541 y_t = ly[(n1/m)]; 542 s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t; 543 for (j=0; j<x_t; j++) idx[nn++] = s_t++; 544 } else if (by == DM_BOUNDARY_MIRROR) { 545 for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1) + j; 546 } 547 548 if (n2 >= 0) { /* right below */ 549 x_t = lx[n2 % m]; 550 y_t = ly[(n2/m)]; 551 s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t; 552 for (j=0; j<s_x; j++) idx[nn++] = s_t++; 553 } 554 } 555 556 for (i=0; i<y; i++) { 557 if (n3 >= 0) { /* directly left */ 558 x_t = lx[n3 % m]; 559 /* y_t = y; */ 560 s_t = bases[n3] + (i+1)*x_t - s_x; 561 for (j=0; j<s_x; j++) idx[nn++] = s_t++; 562 } else if (bx == DM_BOUNDARY_MIRROR) { 563 for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j; 564 } 565 566 for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */ 567 568 if (n5 >= 0) { /* directly right */ 569 x_t = lx[n5 % m]; 570 /* y_t = y; */ 571 s_t = bases[n5] + (i)*x_t; 572 for (j=0; j<s_x; j++) idx[nn++] = s_t++; 573 } else if (bx == DM_BOUNDARY_MIRROR) { 574 for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j; 575 } 576 } 577 578 for (i=1; i<=s_y; i++) { 579 if (n6 >= 0) { /* left above */ 580 x_t = lx[n6 % m]; 581 /* y_t = ly[(n6/m)]; */ 582 s_t = bases[n6] + (i)*x_t - s_x; 583 for (j=0; j<s_x; j++) idx[nn++] = s_t++; 584 } 585 586 if (n7 >= 0) { /* directly above */ 587 x_t = x; 588 /* y_t = ly[(n7/m)]; */ 589 s_t = bases[n7] + (i-1)*x_t; 590 for (j=0; j<x_t; j++) idx[nn++] = s_t++; 591 } else if (by == DM_BOUNDARY_MIRROR) { 592 for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1) + j; 593 } 594 595 if (n8 >= 0) { /* right above */ 596 x_t = lx[n8 % m]; 597 /* y_t = ly[(n8/m)]; */ 598 s_t = bases[n8] + (i-1)*x_t; 599 for (j=0; j<s_x; j++) idx[nn++] = s_t++; 600 } 601 } 602 603 CHKERRQ(ISCreateBlock(comm,dof,nn,idx,PETSC_USE_POINTER,&from)); 604 CHKERRQ(VecScatterCreate(global,from,local,to,>ol)); 605 CHKERRQ(PetscLogObjectParent((PetscObject)da,(PetscObject)gtol)); 606 CHKERRQ(ISDestroy(&to)); 607 CHKERRQ(ISDestroy(&from)); 608 609 if (stencil_type == DMDA_STENCIL_STAR) { 610 n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8; 611 } 612 613 if (((stencil_type == DMDA_STENCIL_STAR) || (bx && bx != DM_BOUNDARY_PERIODIC) || (by && by != DM_BOUNDARY_PERIODIC))) { 614 /* 615 Recompute the local to global mappings, this time keeping the 616 information about the cross corner processor numbers and any ghosted 617 but not periodic indices. 618 */ 619 nn = 0; 620 xbase = bases[rank]; 621 for (i=1; i<=s_y; i++) { 622 if (n0 >= 0) { /* left below */ 623 x_t = lx[n0 % m]; 624 y_t = ly[(n0/m)]; 625 s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x; 626 for (j=0; j<s_x; j++) idx[nn++] = s_t++; 627 } else if (xs-Xs > 0 && ys-Ys > 0) { 628 for (j=0; j<s_x; j++) idx[nn++] = -1; 629 } 630 if (n1 >= 0) { /* directly below */ 631 x_t = x; 632 y_t = ly[(n1/m)]; 633 s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t; 634 for (j=0; j<x_t; j++) idx[nn++] = s_t++; 635 } else if (ys-Ys > 0) { 636 if (by == DM_BOUNDARY_MIRROR) { 637 for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1) + j; 638 } else { 639 for (j=0; j<x; j++) idx[nn++] = -1; 640 } 641 } 642 if (n2 >= 0) { /* right below */ 643 x_t = lx[n2 % m]; 644 y_t = ly[(n2/m)]; 645 s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t; 646 for (j=0; j<s_x; j++) idx[nn++] = s_t++; 647 } else if (Xe-xe> 0 && ys-Ys > 0) { 648 for (j=0; j<s_x; j++) idx[nn++] = -1; 649 } 650 } 651 652 for (i=0; i<y; i++) { 653 if (n3 >= 0) { /* directly left */ 654 x_t = lx[n3 % m]; 655 /* y_t = y; */ 656 s_t = bases[n3] + (i+1)*x_t - s_x; 657 for (j=0; j<s_x; j++) idx[nn++] = s_t++; 658 } else if (xs-Xs > 0) { 659 if (bx == DM_BOUNDARY_MIRROR) { 660 for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j; 661 } else { 662 for (j=0; j<s_x; j++) idx[nn++] = -1; 663 } 664 } 665 666 for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */ 667 668 if (n5 >= 0) { /* directly right */ 669 x_t = lx[n5 % m]; 670 /* y_t = y; */ 671 s_t = bases[n5] + (i)*x_t; 672 for (j=0; j<s_x; j++) idx[nn++] = s_t++; 673 } else if (Xe-xe > 0) { 674 if (bx == DM_BOUNDARY_MIRROR) { 675 for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j; 676 } else { 677 for (j=0; j<s_x; j++) idx[nn++] = -1; 678 } 679 } 680 } 681 682 for (i=1; i<=s_y; i++) { 683 if (n6 >= 0) { /* left above */ 684 x_t = lx[n6 % m]; 685 /* y_t = ly[(n6/m)]; */ 686 s_t = bases[n6] + (i)*x_t - s_x; 687 for (j=0; j<s_x; j++) idx[nn++] = s_t++; 688 } else if (xs-Xs > 0 && Ye-ye > 0) { 689 for (j=0; j<s_x; j++) idx[nn++] = -1; 690 } 691 if (n7 >= 0) { /* directly above */ 692 x_t = x; 693 /* y_t = ly[(n7/m)]; */ 694 s_t = bases[n7] + (i-1)*x_t; 695 for (j=0; j<x_t; j++) idx[nn++] = s_t++; 696 } else if (Ye-ye > 0) { 697 if (by == DM_BOUNDARY_MIRROR) { 698 for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1) + j; 699 } else { 700 for (j=0; j<x; j++) idx[nn++] = -1; 701 } 702 } 703 if (n8 >= 0) { /* right above */ 704 x_t = lx[n8 % m]; 705 /* y_t = ly[(n8/m)]; */ 706 s_t = bases[n8] + (i-1)*x_t; 707 for (j=0; j<s_x; j++) idx[nn++] = s_t++; 708 } else if (Xe-xe > 0 && Ye-ye > 0) { 709 for (j=0; j<s_x; j++) idx[nn++] = -1; 710 } 711 } 712 } 713 /* 714 Set the local to global ordering in the global vector, this allows use 715 of VecSetValuesLocal(). 716 */ 717 CHKERRQ(ISLocalToGlobalMappingCreate(comm,dof,nn,idx,PETSC_OWN_POINTER,&da->ltogmap)); 718 CHKERRQ(PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap)); 719 720 CHKERRQ(PetscFree2(bases,ldims)); 721 dd->m = m; dd->n = n; 722 /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */ 723 dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = 0; dd->ze = 1; 724 dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = 0; dd->Ze = 1; 725 726 CHKERRQ(VecDestroy(&local)); 727 CHKERRQ(VecDestroy(&global)); 728 729 dd->gtol = gtol; 730 dd->base = base; 731 da->ops->view = DMView_DA_2d; 732 dd->ltol = NULL; 733 dd->ao = NULL; 734 PetscFunctionReturn(0); 735 } 736 737 /*@C 738 DMDACreate2d - Creates an object that will manage the communication of two-dimensional 739 regular array data that is distributed across some processors. 740 741 Collective 742 743 Input Parameters: 744 + comm - MPI communicator 745 . bx,by - type of ghost nodes the array have. 746 Use one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC. 747 . stencil_type - stencil type. Use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR. 748 . M,N - global dimension in each direction of the array 749 . m,n - corresponding number of processors in each dimension 750 (or PETSC_DECIDE to have calculated) 751 . dof - number of degrees of freedom per node 752 . s - stencil width 753 - lx, ly - arrays containing the number of nodes in each cell along 754 the x and y coordinates, or NULL. If non-null, these 755 must be of length as m and n, and the corresponding 756 m and n cannot be PETSC_DECIDE. The sum of the lx[] entries 757 must be M, and the sum of the ly[] entries must be N. 758 759 Output Parameter: 760 . da - the resulting distributed array object 761 762 Options Database Key: 763 + -dm_view - Calls DMView() at the conclusion of DMDACreate2d() 764 . -da_grid_x <nx> - number of grid points in x direction 765 . -da_grid_y <ny> - number of grid points in y direction 766 . -da_processors_x <nx> - number of processors in x direction 767 . -da_processors_y <ny> - number of processors in y direction 768 . -da_refine_x <rx> - refinement ratio in x direction 769 . -da_refine_y <ry> - refinement ratio in y direction 770 - -da_refine <n> - refine the DMDA n times before creating 771 772 Level: beginner 773 774 Notes: 775 The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the 776 standard 5-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes 777 the standard 9-pt stencil. 778 779 The array data itself is NOT stored in the DMDA, it is stored in Vec objects; 780 The appropriate vector objects can be obtained with calls to DMCreateGlobalVector() 781 and DMCreateLocalVector() and calls to VecDuplicate() if more are needed. 782 783 You must call DMSetUp() after this call before using this DM. 784 785 If you wish to use the options database to change values in the DMDA call DMSetFromOptions() after this call 786 but before DMSetUp(). 787 788 .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate3d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(), 789 DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDASetRefinementFactor(), 790 DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges(), 791 DMStagCreate2d() 792 793 @*/ 794 795 PetscErrorCode DMDACreate2d(MPI_Comm comm,DMBoundaryType bx,DMBoundaryType by,DMDAStencilType stencil_type, 796 PetscInt M,PetscInt N,PetscInt m,PetscInt n,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],DM *da) 797 { 798 PetscFunctionBegin; 799 CHKERRQ(DMDACreate(comm, da)); 800 CHKERRQ(DMSetDimension(*da, 2)); 801 CHKERRQ(DMDASetSizes(*da, M, N, 1)); 802 CHKERRQ(DMDASetNumProcs(*da, m, n, PETSC_DECIDE)); 803 CHKERRQ(DMDASetBoundaryType(*da, bx, by, DM_BOUNDARY_NONE)); 804 CHKERRQ(DMDASetDof(*da, dof)); 805 CHKERRQ(DMDASetStencilType(*da, stencil_type)); 806 CHKERRQ(DMDASetStencilWidth(*da, s)); 807 CHKERRQ(DMDASetOwnershipRanges(*da, lx, ly, NULL)); 808 PetscFunctionReturn(0); 809 } 810