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