17d0a6c19SBarry Smith 247c6ae99SBarry Smith /* 347c6ae99SBarry Smith Code for manipulating distributed regular 3d arrays in parallel. 447c6ae99SBarry Smith File created by Peter Mell 7/14/95 547c6ae99SBarry Smith */ 647c6ae99SBarry Smith 7af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 847c6ae99SBarry Smith 99804daf3SBarry Smith #include <petscdraw.h> 1047c6ae99SBarry Smith #undef __FUNCT__ 119a42bb27SBarry Smith #define __FUNCT__ "DMView_DA_3d" 129a42bb27SBarry Smith PetscErrorCode DMView_DA_3d(DM da,PetscViewer viewer) 1347c6ae99SBarry Smith { 1447c6ae99SBarry Smith PetscErrorCode ierr; 1547c6ae99SBarry Smith PetscMPIInt rank; 169a42bb27SBarry Smith PetscBool iascii,isdraw,isbinary; 1747c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 189a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 199a42bb27SBarry Smith PetscBool ismatlab; 209a42bb27SBarry Smith #endif 2147c6ae99SBarry Smith 2247c6ae99SBarry Smith PetscFunctionBegin; 23ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank);CHKERRQ(ierr); 2447c6ae99SBarry Smith 25251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 26251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 27251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 289a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 29251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab);CHKERRQ(ierr); 309a42bb27SBarry Smith #endif 3147c6ae99SBarry Smith if (iascii) { 3247c6ae99SBarry Smith PetscViewerFormat format; 3347c6ae99SBarry Smith 341575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 3547c6ae99SBarry Smith ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 3647c6ae99SBarry Smith if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL) { 37aa219208SBarry Smith DMDALocalInfo info; 38aa219208SBarry Smith ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr); 3947c6ae99SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Processor [%d] M %D N %D P %D m %D n %D p %D w %D s %D\n",rank,dd->M,dd->N,dd->P,dd->m,dd->n,dd->p,dd->w,dd->s);CHKERRQ(ierr); 4047c6ae99SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"X range of indices: %D %D, Y range of indices: %D %D, Z range of indices: %D %D\n", 4147c6ae99SBarry Smith info.xs,info.xs+info.xm,info.ys,info.ys+info.ym,info.zs,info.zs+info.zm);CHKERRQ(ierr); 4247c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX) 436636e97aSMatthew G Knepley if (da->coordinates) { 4447c6ae99SBarry Smith PetscInt last; 4547c6ae99SBarry Smith const PetscReal *coors; 466636e97aSMatthew G Knepley ierr = VecGetArrayRead(da->coordinates,&coors);CHKERRQ(ierr); 476636e97aSMatthew G Knepley ierr = VecGetLocalSize(da->coordinates,&last);CHKERRQ(ierr); 4847c6ae99SBarry Smith last = last - 3; 4957622a8eSBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Lower left corner %g %g %g : Upper right %g %g %g\n",(double)coors[0],(double)coors[1],(double)coors[2],(double)coors[last],(double)coors[last+1],(double)coors[last+2]);CHKERRQ(ierr); 506636e97aSMatthew G Knepley ierr = VecRestoreArrayRead(da->coordinates,&coors);CHKERRQ(ierr); 5147c6ae99SBarry Smith } 5247c6ae99SBarry Smith #endif 5347c6ae99SBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 541575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 55616157d6SJed Brown } else { 56616157d6SJed Brown ierr = DMView_DA_VTK(da,viewer);CHKERRQ(ierr); 5747c6ae99SBarry Smith } 5847c6ae99SBarry Smith } else if (isdraw) { 5947c6ae99SBarry Smith PetscDraw draw; 6047c6ae99SBarry Smith PetscReal ymin = -1.0,ymax = (PetscReal)dd->N; 6147c6ae99SBarry Smith PetscReal xmin = -1.0,xmax = (PetscReal)((dd->M+2)*dd->P),x,y,ycoord,xcoord; 628ea3bf28SBarry Smith PetscInt k,plane,base; 638ea3bf28SBarry Smith const PetscInt *idx; 6447c6ae99SBarry Smith char node[10]; 6547c6ae99SBarry Smith PetscBool isnull; 6647c6ae99SBarry Smith 6747c6ae99SBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 685b399a63SLisandro Dalcin ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 695b399a63SLisandro Dalcin if (isnull) PetscFunctionReturn(0); 7047c6ae99SBarry Smith 715b399a63SLisandro Dalcin ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr); 725b399a63SLisandro Dalcin ierr = PetscDrawClear(draw);CHKERRQ(ierr); 735b399a63SLisandro Dalcin ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr); 745b399a63SLisandro Dalcin 755b399a63SLisandro Dalcin ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr); 7647c6ae99SBarry Smith /* first processor draw all node lines */ 7747c6ae99SBarry Smith if (!rank) { 7847c6ae99SBarry Smith for (k=0; k<dd->P; k++) { 7947c6ae99SBarry Smith ymin = 0.0; ymax = (PetscReal)(dd->N - 1); 8047c6ae99SBarry Smith for (xmin=(PetscReal)(k*(dd->M+1)); xmin<(PetscReal)(dd->M+(k*(dd->M+1))); xmin++) { 8147c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr); 8247c6ae99SBarry Smith } 8347c6ae99SBarry Smith xmin = (PetscReal)(k*(dd->M+1)); xmax = xmin + (PetscReal)(dd->M - 1); 8447c6ae99SBarry Smith for (ymin=0; ymin<(PetscReal)dd->N; ymin++) { 8547c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr); 8647c6ae99SBarry Smith } 8747c6ae99SBarry Smith } 8847c6ae99SBarry Smith } 895b399a63SLisandro Dalcin ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr); 905b399a63SLisandro Dalcin ierr = PetscDrawFlush(draw);CHKERRQ(ierr); 9147c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 9247c6ae99SBarry Smith 935b399a63SLisandro Dalcin ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr); 945b399a63SLisandro Dalcin /*Go through and draw for each plane*/ 955b399a63SLisandro Dalcin for (k=0; k<dd->P; k++) { 9647c6ae99SBarry Smith if ((k >= dd->zs) && (k < dd->ze)) { 9747c6ae99SBarry Smith /* draw my box */ 9847c6ae99SBarry Smith ymin = dd->ys; 9947c6ae99SBarry Smith ymax = dd->ye - 1; 10047c6ae99SBarry Smith xmin = dd->xs/dd->w + (dd->M+1)*k; 10147c6ae99SBarry Smith xmax =(dd->xe-1)/dd->w + (dd->M+1)*k; 10247c6ae99SBarry Smith 10347c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr); 10447c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 10547c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 10647c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 10747c6ae99SBarry Smith 10847c6ae99SBarry Smith xmin = dd->xs/dd->w; 10947c6ae99SBarry Smith xmax =(dd->xe-1)/dd->w; 11047c6ae99SBarry Smith 111*832b7cebSLisandro Dalcin /* identify which processor owns the box */ 112*832b7cebSLisandro Dalcin ierr = PetscSNPrintf(node,sizeof(node),"%d",(int)rank);CHKERRQ(ierr); 113*832b7cebSLisandro Dalcin ierr = PetscDrawString(draw,xmin+(dd->M+1)*k+.2,ymin+.3,PETSC_DRAW_RED,node);CHKERRQ(ierr); 11447c6ae99SBarry Smith /* put in numbers*/ 11547c6ae99SBarry Smith base = (dd->base+(dd->xe-dd->xs)*(dd->ye-dd->ys)*(k-dd->zs))/dd->w; 11647c6ae99SBarry Smith for (y=ymin; y<=ymax; y++) { 11747c6ae99SBarry Smith for (x=xmin+(dd->M+1)*k; x<=xmax+(dd->M+1)*k; x++) { 118*832b7cebSLisandro Dalcin ierr = PetscSNPrintf(node,sizeof(node),"%d",(int)base++);CHKERRQ(ierr); 11947c6ae99SBarry Smith ierr = PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);CHKERRQ(ierr); 12047c6ae99SBarry Smith } 12147c6ae99SBarry Smith } 12247c6ae99SBarry Smith 12347c6ae99SBarry Smith } 12447c6ae99SBarry Smith } 1255b399a63SLisandro Dalcin ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr); 1265b399a63SLisandro Dalcin ierr = PetscDrawFlush(draw);CHKERRQ(ierr); 12747c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 12847c6ae99SBarry Smith 1295b399a63SLisandro Dalcin ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr); 13047c6ae99SBarry Smith for (k=0-dd->s; k<dd->P+dd->s; k++) { 13147c6ae99SBarry Smith /* Go through and draw for each plane */ 13247c6ae99SBarry Smith if ((k >= dd->Zs) && (k < dd->Ze)) { 13347c6ae99SBarry Smith /* overlay ghost numbers, useful for error checking */ 134565245c5SBarry Smith base = (dd->Xe-dd->Xs)*(dd->Ye-dd->Ys)*(k-dd->Zs)/dd->w; 13545b6f7e9SBarry Smith ierr = ISLocalToGlobalMappingGetBlockIndices(da->ltogmap,&idx);CHKERRQ(ierr); 13647c6ae99SBarry Smith plane=k; 137565245c5SBarry Smith /* Keep z wrap around points on the drawing */ 1388865f1eaSKarl Rupp if (k<0) plane=dd->P+k; 1398865f1eaSKarl Rupp if (k>=dd->P) plane=k-dd->P; 14047c6ae99SBarry Smith ymin = dd->Ys; ymax = dd->Ye; 14147c6ae99SBarry Smith xmin = (dd->M+1)*plane*dd->w; 14247c6ae99SBarry Smith xmax = (dd->M+1)*plane*dd->w+dd->M*dd->w; 14347c6ae99SBarry Smith for (y=ymin; y<ymax; y++) { 14447c6ae99SBarry Smith for (x=xmin+dd->Xs; x<xmin+dd->Xe; x+=dd->w) { 1458ea3bf28SBarry Smith sprintf(node,"%d",(int)(idx[base])); 14647c6ae99SBarry Smith ycoord = y; 14747c6ae99SBarry Smith /*Keep y wrap around points on drawing */ 1488865f1eaSKarl Rupp if (y<0) ycoord = dd->N+y; 1498865f1eaSKarl Rupp if (y>=dd->N) ycoord = y-dd->N; 15047c6ae99SBarry Smith xcoord = x; /* Keep x wrap points on drawing */ 1518865f1eaSKarl Rupp if (x<xmin) xcoord = xmax - (xmin-x); 1528865f1eaSKarl Rupp if (x>=xmax) xcoord = xmin + (x-xmax); 15347c6ae99SBarry Smith ierr = PetscDrawString(draw,xcoord/dd->w,ycoord,PETSC_DRAW_BLUE,node);CHKERRQ(ierr); 154565245c5SBarry Smith base++; 15547c6ae99SBarry Smith } 15647c6ae99SBarry Smith } 15745b6f7e9SBarry Smith ierr = ISLocalToGlobalMappingRestoreBlockIndices(da->ltogmap,&idx);CHKERRQ(ierr); 15847c6ae99SBarry Smith } 15947c6ae99SBarry Smith } 1605b399a63SLisandro Dalcin ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr); 1615b399a63SLisandro Dalcin ierr = PetscDrawFlush(draw);CHKERRQ(ierr); 16247c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 163*832b7cebSLisandro Dalcin ierr = PetscDrawSave(draw);CHKERRQ(ierr); 1649a42bb27SBarry Smith } else if (isbinary) { 1659a42bb27SBarry Smith ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr); 1669a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 1679a42bb27SBarry Smith } else if (ismatlab) { 1689a42bb27SBarry Smith ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr); 1699a42bb27SBarry Smith #endif 17011aeaf0aSBarry Smith } 17147c6ae99SBarry Smith PetscFunctionReturn(0); 17247c6ae99SBarry Smith } 17347c6ae99SBarry Smith 17447c6ae99SBarry Smith #undef __FUNCT__ 1759a42bb27SBarry Smith #define __FUNCT__ "DMSetUp_DA_3D" 1767087cfbeSBarry Smith PetscErrorCode DMSetUp_DA_3D(DM da) 17747c6ae99SBarry Smith { 17847c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 17947c6ae99SBarry Smith const PetscInt M = dd->M; 18047c6ae99SBarry Smith const PetscInt N = dd->N; 18147c6ae99SBarry Smith const PetscInt P = dd->P; 18247c6ae99SBarry Smith PetscInt m = dd->m; 18347c6ae99SBarry Smith PetscInt n = dd->n; 18447c6ae99SBarry Smith PetscInt p = dd->p; 18547c6ae99SBarry Smith const PetscInt dof = dd->w; 18647c6ae99SBarry Smith const PetscInt s = dd->s; 187bff4a2f0SMatthew G. Knepley DMBoundaryType bx = dd->bx; 188bff4a2f0SMatthew G. Knepley DMBoundaryType by = dd->by; 189bff4a2f0SMatthew G. Knepley DMBoundaryType bz = dd->bz; 19019fd82e9SBarry Smith DMDAStencilType stencil_type = dd->stencil_type; 19147c6ae99SBarry Smith PetscInt *lx = dd->lx; 19247c6ae99SBarry Smith PetscInt *ly = dd->ly; 19347c6ae99SBarry Smith PetscInt *lz = dd->lz; 19447c6ae99SBarry Smith MPI_Comm comm; 19547c6ae99SBarry Smith PetscMPIInt rank,size; 196ce00eea3SSatish Balay PetscInt xs = 0,xe,ys = 0,ye,zs = 0,ze,x = 0,y = 0,z = 0; 197bd1fc5aeSBarry Smith PetscInt Xs,Xe,Ys,Ye,Zs,Ze,IXs,IXe,IYs,IYe,IZs,IZe,pm; 1988ea3bf28SBarry Smith PetscInt left,right,up,down,bottom,top,i,j,k,*idx,nn; 19947c6ae99SBarry Smith PetscInt n0,n1,n2,n3,n4,n5,n6,n7,n8,n9,n10,n11,n12,n14; 20047c6ae99SBarry Smith PetscInt n15,n16,n17,n18,n19,n20,n21,n22,n23,n24,n25,n26; 201db87c5ecSEthan Coon PetscInt *bases,*ldims,base,x_t,y_t,z_t,s_t,count,s_x,s_y,s_z; 20247c6ae99SBarry Smith PetscInt sn0 = 0,sn1 = 0,sn2 = 0,sn3 = 0,sn5 = 0,sn6 = 0,sn7 = 0; 20347c6ae99SBarry Smith PetscInt sn8 = 0,sn9 = 0,sn11 = 0,sn15 = 0,sn24 = 0,sn25 = 0,sn26 = 0; 20447c6ae99SBarry Smith PetscInt sn17 = 0,sn18 = 0,sn19 = 0,sn20 = 0,sn21 = 0,sn23 = 0; 20547c6ae99SBarry Smith Vec local,global; 206bd1fc5aeSBarry Smith VecScatter gtol; 20745b6f7e9SBarry Smith IS to,from; 2086f951b95Secoon PetscBool twod; 20947c6ae99SBarry Smith PetscErrorCode ierr; 21047c6ae99SBarry Smith 2116f951b95Secoon 21247c6ae99SBarry Smith PetscFunctionBegin; 213bff4a2f0SMatthew G. Knepley if (stencil_type == DMDA_STENCIL_BOX && (bx == DM_BOUNDARY_MIRROR || by == DM_BOUNDARY_MIRROR || bz == DM_BOUNDARY_MIRROR)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Mirror boundary and box stencil"); 21447c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject) da, &comm);CHKERRQ(ierr); 2153855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES) 2163855c12bSBarry Smith if (((Petsc64bitInt) M)*((Petsc64bitInt) N)*((Petsc64bitInt) P)*((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); 2173855c12bSBarry Smith #endif 2183855c12bSBarry Smith 21947c6ae99SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 22047c6ae99SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 22147c6ae99SBarry Smith 22247c6ae99SBarry Smith if (m != PETSC_DECIDE) { 22347c6ae99SBarry Smith if (m < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m); 22447c6ae99SBarry Smith else if (m > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size); 22547c6ae99SBarry Smith } 22647c6ae99SBarry Smith if (n != PETSC_DECIDE) { 22747c6ae99SBarry Smith if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n); 22847c6ae99SBarry Smith else if (n > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size); 22947c6ae99SBarry Smith } 23047c6ae99SBarry Smith if (p != PETSC_DECIDE) { 23147c6ae99SBarry Smith if (p < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Z direction: %D",p); 23247c6ae99SBarry Smith else if (p > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Z direction: %D %d",p,size); 23347c6ae99SBarry Smith } 2340716a85fSBarry Smith if ((m > 0) && (n > 0) && (p > 0) && (m*n*p != size)) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"m %D * n %D * p %D != size %d",m,n,p,size); 23547c6ae99SBarry Smith 23647c6ae99SBarry Smith /* Partition the array among the processors */ 23747c6ae99SBarry Smith if (m == PETSC_DECIDE && n != PETSC_DECIDE && p != PETSC_DECIDE) { 23847c6ae99SBarry Smith m = size/(n*p); 23947c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) { 24047c6ae99SBarry Smith n = size/(m*p); 24147c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) { 24247c6ae99SBarry Smith p = size/(m*n); 24347c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) { 24447c6ae99SBarry Smith /* try for squarish distribution */ 245369cc0aeSBarry Smith m = (int)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)N*p))); 24647c6ae99SBarry Smith if (!m) m = 1; 24747c6ae99SBarry Smith while (m > 0) { 24847c6ae99SBarry Smith n = size/(m*p); 24947c6ae99SBarry Smith if (m*n*p == size) break; 25047c6ae99SBarry Smith m--; 25147c6ae99SBarry Smith } 25247c6ae99SBarry Smith if (!m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad p value: p = %D",p); 25347c6ae99SBarry Smith if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;} 25447c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) { 25547c6ae99SBarry Smith /* try for squarish distribution */ 256369cc0aeSBarry Smith m = (int)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n))); 25747c6ae99SBarry Smith if (!m) m = 1; 25847c6ae99SBarry Smith while (m > 0) { 25947c6ae99SBarry Smith p = size/(m*n); 26047c6ae99SBarry Smith if (m*n*p == size) break; 26147c6ae99SBarry Smith m--; 26247c6ae99SBarry Smith } 26347c6ae99SBarry Smith if (!m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad n value: n = %D",n); 26447c6ae99SBarry Smith if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;} 26547c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) { 26647c6ae99SBarry Smith /* try for squarish distribution */ 267369cc0aeSBarry Smith n = (int)(0.5 + PetscSqrtReal(((PetscReal)N)*((PetscReal)size)/((PetscReal)P*m))); 26847c6ae99SBarry Smith if (!n) n = 1; 26947c6ae99SBarry Smith while (n > 0) { 27047c6ae99SBarry Smith p = size/(m*n); 27147c6ae99SBarry Smith if (m*n*p == size) break; 27247c6ae99SBarry Smith n--; 27347c6ae99SBarry Smith } 27447c6ae99SBarry Smith if (!n) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad m value: m = %D",n); 27547c6ae99SBarry Smith if (N > P && n < p) {PetscInt _n = n; n = p; p = _n;} 27647c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) { 27747c6ae99SBarry Smith /* try for squarish distribution */ 278369cc0aeSBarry Smith n = (PetscInt)(0.5 + PetscPowReal(((PetscReal)N*N)*((PetscReal)size)/((PetscReal)P*M),(PetscReal)(1./3.))); 27947c6ae99SBarry Smith if (!n) n = 1; 28047c6ae99SBarry Smith while (n > 0) { 28147c6ae99SBarry Smith pm = size/n; 28247c6ae99SBarry Smith if (n*pm == size) break; 28347c6ae99SBarry Smith n--; 28447c6ae99SBarry Smith } 28547c6ae99SBarry Smith if (!n) n = 1; 286369cc0aeSBarry Smith m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n))); 28747c6ae99SBarry Smith if (!m) m = 1; 28847c6ae99SBarry Smith while (m > 0) { 28947c6ae99SBarry Smith p = size/(m*n); 29047c6ae99SBarry Smith if (m*n*p == size) break; 29147c6ae99SBarry Smith m--; 29247c6ae99SBarry Smith } 29347c6ae99SBarry Smith if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;} 294ce94432eSBarry Smith } else if (m*n*p != size) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition"); 29547c6ae99SBarry Smith 296ce94432eSBarry Smith if (m*n*p != size) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_PLIB,"Could not find good partition"); 297ce94432eSBarry Smith if (M < m) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m); 298ce94432eSBarry Smith if (N < n) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n); 299ce94432eSBarry Smith if (P < p) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Partition in z direction is too fine! %D %D",P,p); 30047c6ae99SBarry Smith 30147c6ae99SBarry Smith /* 30247c6ae99SBarry Smith Determine locally owned region 30347c6ae99SBarry Smith [x, y, or z]s is the first local node number, [x, y, z] is the number of local nodes 30447c6ae99SBarry Smith */ 30547c6ae99SBarry Smith 30647c6ae99SBarry Smith if (!lx) { 307785e854fSJed Brown ierr = PetscMalloc1(m, &dd->lx);CHKERRQ(ierr); 30847c6ae99SBarry Smith lx = dd->lx; 3098865f1eaSKarl Rupp for (i=0; i<m; i++) lx[i] = M/m + ((M % m) > (i % m)); 31047c6ae99SBarry Smith } 31147c6ae99SBarry Smith x = lx[rank % m]; 31247c6ae99SBarry Smith xs = 0; 3138865f1eaSKarl Rupp for (i=0; i<(rank%m); i++) xs += lx[i]; 314bff4a2f0SMatthew G. Knepley 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); 31547c6ae99SBarry Smith 31647c6ae99SBarry Smith if (!ly) { 317785e854fSJed Brown ierr = PetscMalloc1(n, &dd->ly);CHKERRQ(ierr); 31847c6ae99SBarry Smith ly = dd->ly; 3198865f1eaSKarl Rupp for (i=0; i<n; i++) ly[i] = N/n + ((N % n) > (i % n)); 32047c6ae99SBarry Smith } 32147c6ae99SBarry Smith y = ly[(rank % (m*n))/m]; 322bff4a2f0SMatthew G. Knepley 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); 32330729d88SBarry Smith 32447c6ae99SBarry Smith ys = 0; 3258865f1eaSKarl Rupp for (i=0; i<(rank % (m*n))/m; i++) ys += ly[i]; 32647c6ae99SBarry Smith 32747c6ae99SBarry Smith if (!lz) { 328785e854fSJed Brown ierr = PetscMalloc1(p, &dd->lz);CHKERRQ(ierr); 32947c6ae99SBarry Smith lz = dd->lz; 3308865f1eaSKarl Rupp for (i=0; i<p; i++) lz[i] = P/p + ((P % p) > (i % p)); 33147c6ae99SBarry Smith } 33247c6ae99SBarry Smith z = lz[rank/(m*n)]; 333bcea557cSEthan Coon 334fdc81ce8SEthan Coon /* note this is different than x- and y-, as we will handle as an important special 335bff4a2f0SMatthew G. Knepley case when p=P=1 and DM_BOUNDARY_PERIODIC and s > z. This is to deal with 2D problems 336fdc81ce8SEthan Coon in a 3D code. Additional code for this case is noted with "2d case" comments */ 3376f951b95Secoon twod = PETSC_FALSE; 3388865f1eaSKarl Rupp if (P == 1) twod = PETSC_TRUE; 339bff4a2f0SMatthew G. Knepley else if ((z < s) && ((p > 1) || (bz == DM_BOUNDARY_PERIODIC))) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local z-width of domain z %D is smaller than stencil width s %D",z,s); 34047c6ae99SBarry Smith zs = 0; 3418865f1eaSKarl Rupp for (i=0; i<(rank/(m*n)); i++) zs += lz[i]; 34247c6ae99SBarry Smith ye = ys + y; 34347c6ae99SBarry Smith xe = xs + x; 34447c6ae99SBarry Smith ze = zs + z; 34547c6ae99SBarry Smith 34688661749SPeter Brune /* determine ghost region (Xs) and region scattered into (IXs) */ 347d9c9ebe5SPeter Brune if (xs-s > 0) { 348d9c9ebe5SPeter Brune Xs = xs - s; IXs = xs - s; 34988661749SPeter Brune } else { 3508865f1eaSKarl Rupp if (bx) Xs = xs - s; 3518865f1eaSKarl Rupp else Xs = 0; 35288661749SPeter Brune IXs = 0; 35388661749SPeter Brune } 354d9c9ebe5SPeter Brune if (xe+s <= M) { 355d9c9ebe5SPeter Brune Xe = xe + s; IXe = xe + s; 35688661749SPeter Brune } else { 35788661749SPeter Brune if (bx) { 358d9c9ebe5SPeter Brune Xs = xs - s; Xe = xe + s; 3598865f1eaSKarl Rupp } else Xe = M; 36088661749SPeter Brune IXe = M; 36188661749SPeter Brune } 36247c6ae99SBarry Smith 363bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) { 364d9c9ebe5SPeter Brune IXs = xs - s; 365d9c9ebe5SPeter Brune IXe = xe + s; 366d9c9ebe5SPeter Brune Xs = xs - s; 367d9c9ebe5SPeter Brune Xe = xe + s; 36888661749SPeter Brune } 36988661749SPeter Brune 370d9c9ebe5SPeter Brune if (ys-s > 0) { 371d9c9ebe5SPeter Brune Ys = ys - s; IYs = ys - s; 37288661749SPeter Brune } else { 3738865f1eaSKarl Rupp if (by) Ys = ys - s; 3748865f1eaSKarl Rupp else Ys = 0; 37588661749SPeter Brune IYs = 0; 37688661749SPeter Brune } 377d9c9ebe5SPeter Brune if (ye+s <= N) { 378d9c9ebe5SPeter Brune Ye = ye + s; IYe = ye + s; 37988661749SPeter Brune } else { 3808865f1eaSKarl Rupp if (by) Ye = ye + s; 3818865f1eaSKarl Rupp else Ye = N; 38288661749SPeter Brune IYe = N; 38388661749SPeter Brune } 38488661749SPeter Brune 385bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC || by == DM_BOUNDARY_MIRROR) { 386d9c9ebe5SPeter Brune IYs = ys - s; 387d9c9ebe5SPeter Brune IYe = ye + s; 388d9c9ebe5SPeter Brune Ys = ys - s; 389d9c9ebe5SPeter Brune Ye = ye + s; 39088661749SPeter Brune } 39188661749SPeter Brune 392d9c9ebe5SPeter Brune if (zs-s > 0) { 393d9c9ebe5SPeter Brune Zs = zs - s; IZs = zs - s; 39488661749SPeter Brune } else { 3958865f1eaSKarl Rupp if (bz) Zs = zs - s; 3968865f1eaSKarl Rupp else Zs = 0; 39788661749SPeter Brune IZs = 0; 39888661749SPeter Brune } 399d9c9ebe5SPeter Brune if (ze+s <= P) { 400d9c9ebe5SPeter Brune Ze = ze + s; IZe = ze + s; 40188661749SPeter Brune } else { 4028865f1eaSKarl Rupp if (bz) Ze = ze + s; 4038865f1eaSKarl Rupp else Ze = P; 40488661749SPeter Brune IZe = P; 40588661749SPeter Brune } 40688661749SPeter Brune 407bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_PERIODIC || bz == DM_BOUNDARY_MIRROR) { 408d9c9ebe5SPeter Brune IZs = zs - s; 409d9c9ebe5SPeter Brune IZe = ze + s; 410d9c9ebe5SPeter Brune Zs = zs - s; 411d9c9ebe5SPeter Brune Ze = ze + s; 41288661749SPeter Brune } 41347c6ae99SBarry Smith 41447c6ae99SBarry Smith /* Resize all X parameters to reflect w */ 415d9c9ebe5SPeter Brune s_x = s; 416d9c9ebe5SPeter Brune s_y = s; 417d9c9ebe5SPeter Brune s_z = s; 41847c6ae99SBarry Smith 41947c6ae99SBarry Smith /* determine starting point of each processor */ 42047c6ae99SBarry Smith nn = x*y*z; 421dcca6d9dSJed Brown ierr = PetscMalloc2(size+1,&bases,size,&ldims);CHKERRQ(ierr); 42247c6ae99SBarry Smith ierr = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr); 42347c6ae99SBarry Smith bases[0] = 0; 4248865f1eaSKarl Rupp for (i=1; i<=size; i++) bases[i] = ldims[i-1]; 4258865f1eaSKarl Rupp for (i=1; i<=size; i++) bases[i] += bases[i-1]; 426ce00eea3SSatish Balay base = bases[rank]*dof; 42747c6ae99SBarry Smith 42847c6ae99SBarry Smith /* allocate the base parallel and sequential vectors */ 429ce00eea3SSatish Balay dd->Nlocal = x*y*z*dof; 430b1fb7eb7SBarry Smith ierr = VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,NULL,&global);CHKERRQ(ierr); 431ce00eea3SSatish Balay dd->nlocal = (Xe-Xs)*(Ye-Ys)*(Ze-Zs)*dof; 432b1fb7eb7SBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,NULL,&local);CHKERRQ(ierr); 43347c6ae99SBarry Smith 43447c6ae99SBarry Smith /* generate appropriate vector scatters */ 43547c6ae99SBarry Smith /* local to global inserts non-ghost point region into global */ 43602fe608eSBarry Smith ierr = PetscMalloc1((IXe-IXs)*(IYe-IYs)*(IZe-IZs),&idx);CHKERRQ(ierr); 437ce00eea3SSatish Balay left = xs - Xs; right = left + x; 43847c6ae99SBarry Smith bottom = ys - Ys; top = bottom + y; 43947c6ae99SBarry Smith down = zs - Zs; up = down + z; 44047c6ae99SBarry Smith count = 0; 44147c6ae99SBarry Smith for (i=down; i<up; i++) { 44247c6ae99SBarry Smith for (j=bottom; j<top; j++) { 443ce00eea3SSatish Balay for (k=left; k<right; k++) { 444ce00eea3SSatish Balay idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 44547c6ae99SBarry Smith } 44647c6ae99SBarry Smith } 44747c6ae99SBarry Smith } 44847c6ae99SBarry Smith 449ce00eea3SSatish Balay /* global to local must include ghost points within the domain, 450ce00eea3SSatish Balay but not ghost points outside the domain that aren't periodic */ 451d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_BOX) { 452ce00eea3SSatish Balay left = IXs - Xs; right = left + (IXe-IXs); 453ce00eea3SSatish Balay bottom = IYs - Ys; top = bottom + (IYe-IYs); 454ce00eea3SSatish Balay down = IZs - Zs; up = down + (IZe-IZs); 455ce00eea3SSatish Balay count = 0; 456ce00eea3SSatish Balay for (i=down; i<up; i++) { 457ce00eea3SSatish Balay for (j=bottom; j<top; j++) { 458ce00eea3SSatish Balay for (k=left; k<right; k++) { 459ce00eea3SSatish Balay idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 460ce00eea3SSatish Balay } 461ce00eea3SSatish Balay } 462ce00eea3SSatish Balay } 463ce00eea3SSatish Balay ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr); 46447c6ae99SBarry Smith } else { 46547c6ae99SBarry Smith /* This is way ugly! We need to list the funny cross type region */ 466ce00eea3SSatish Balay left = xs - Xs; right = left + x; 46747c6ae99SBarry Smith bottom = ys - Ys; top = bottom + y; 46847c6ae99SBarry Smith down = zs - Zs; up = down + z; 46947c6ae99SBarry Smith count = 0; 470ce00eea3SSatish Balay /* the bottom chunck */ 471ce00eea3SSatish Balay for (i=(IZs-Zs); i<down; i++) { 47247c6ae99SBarry Smith for (j=bottom; j<top; j++) { 473ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 47447c6ae99SBarry Smith } 47547c6ae99SBarry Smith } 47647c6ae99SBarry Smith /* the middle piece */ 47747c6ae99SBarry Smith for (i=down; i<up; i++) { 47847c6ae99SBarry Smith /* front */ 479ce00eea3SSatish Balay for (j=(IYs-Ys); j<bottom; j++) { 480ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 48147c6ae99SBarry Smith } 48247c6ae99SBarry Smith /* middle */ 48347c6ae99SBarry Smith for (j=bottom; j<top; j++) { 484ce00eea3SSatish Balay for (k=IXs-Xs; k<IXe-Xs; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 48547c6ae99SBarry Smith } 48647c6ae99SBarry Smith /* back */ 487ce00eea3SSatish Balay for (j=top; j<top+IYe-ye; j++) { 488ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 48947c6ae99SBarry Smith } 49047c6ae99SBarry Smith } 49147c6ae99SBarry Smith /* the top piece */ 492ce00eea3SSatish Balay for (i=up; i<up+IZe-ze; i++) { 49347c6ae99SBarry Smith for (j=bottom; j<top; j++) { 494ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 49547c6ae99SBarry Smith } 49647c6ae99SBarry Smith } 49747c6ae99SBarry Smith ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr); 49847c6ae99SBarry Smith } 49947c6ae99SBarry Smith 50047c6ae99SBarry Smith /* determine who lies on each side of use stored in n24 n25 n26 50147c6ae99SBarry Smith n21 n22 n23 50247c6ae99SBarry Smith n18 n19 n20 50347c6ae99SBarry Smith 50447c6ae99SBarry Smith n15 n16 n17 50547c6ae99SBarry Smith n12 n14 50647c6ae99SBarry Smith n9 n10 n11 50747c6ae99SBarry Smith 50847c6ae99SBarry Smith n6 n7 n8 50947c6ae99SBarry Smith n3 n4 n5 51047c6ae99SBarry Smith n0 n1 n2 51147c6ae99SBarry Smith */ 51247c6ae99SBarry Smith 51347c6ae99SBarry Smith /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */ 51447c6ae99SBarry Smith /* Assume Nodes are Internal to the Cube */ 51547c6ae99SBarry Smith n0 = rank - m*n - m - 1; 51647c6ae99SBarry Smith n1 = rank - m*n - m; 51747c6ae99SBarry Smith n2 = rank - m*n - m + 1; 51847c6ae99SBarry Smith n3 = rank - m*n -1; 51947c6ae99SBarry Smith n4 = rank - m*n; 52047c6ae99SBarry Smith n5 = rank - m*n + 1; 52147c6ae99SBarry Smith n6 = rank - m*n + m - 1; 52247c6ae99SBarry Smith n7 = rank - m*n + m; 52347c6ae99SBarry Smith n8 = rank - m*n + m + 1; 52447c6ae99SBarry Smith 52547c6ae99SBarry Smith n9 = rank - m - 1; 52647c6ae99SBarry Smith n10 = rank - m; 52747c6ae99SBarry Smith n11 = rank - m + 1; 52847c6ae99SBarry Smith n12 = rank - 1; 52947c6ae99SBarry Smith n14 = rank + 1; 53047c6ae99SBarry Smith n15 = rank + m - 1; 53147c6ae99SBarry Smith n16 = rank + m; 53247c6ae99SBarry Smith n17 = rank + m + 1; 53347c6ae99SBarry Smith 53447c6ae99SBarry Smith n18 = rank + m*n - m - 1; 53547c6ae99SBarry Smith n19 = rank + m*n - m; 53647c6ae99SBarry Smith n20 = rank + m*n - m + 1; 53747c6ae99SBarry Smith n21 = rank + m*n - 1; 53847c6ae99SBarry Smith n22 = rank + m*n; 53947c6ae99SBarry Smith n23 = rank + m*n + 1; 54047c6ae99SBarry Smith n24 = rank + m*n + m - 1; 54147c6ae99SBarry Smith n25 = rank + m*n + m; 54247c6ae99SBarry Smith n26 = rank + m*n + m + 1; 54347c6ae99SBarry Smith 54447c6ae99SBarry Smith /* Assume Pieces are on Faces of Cube */ 54547c6ae99SBarry Smith 54647c6ae99SBarry Smith if (xs == 0) { /* First assume not corner or edge */ 54747c6ae99SBarry Smith n0 = rank -1 - (m*n); 54847c6ae99SBarry Smith n3 = rank + m -1 - (m*n); 54947c6ae99SBarry Smith n6 = rank + 2*m -1 - (m*n); 55047c6ae99SBarry Smith n9 = rank -1; 55147c6ae99SBarry Smith n12 = rank + m -1; 55247c6ae99SBarry Smith n15 = rank + 2*m -1; 55347c6ae99SBarry Smith n18 = rank -1 + (m*n); 55447c6ae99SBarry Smith n21 = rank + m -1 + (m*n); 55547c6ae99SBarry Smith n24 = rank + 2*m -1 + (m*n); 55647c6ae99SBarry Smith } 55747c6ae99SBarry Smith 558ce00eea3SSatish Balay if (xe == M) { /* First assume not corner or edge */ 55947c6ae99SBarry Smith n2 = rank -2*m +1 - (m*n); 56047c6ae99SBarry Smith n5 = rank - m +1 - (m*n); 56147c6ae99SBarry Smith n8 = rank +1 - (m*n); 56247c6ae99SBarry Smith n11 = rank -2*m +1; 56347c6ae99SBarry Smith n14 = rank - m +1; 56447c6ae99SBarry Smith n17 = rank +1; 56547c6ae99SBarry Smith n20 = rank -2*m +1 + (m*n); 56647c6ae99SBarry Smith n23 = rank - m +1 + (m*n); 56747c6ae99SBarry Smith n26 = rank +1 + (m*n); 56847c6ae99SBarry Smith } 56947c6ae99SBarry Smith 57047c6ae99SBarry Smith if (ys==0) { /* First assume not corner or edge */ 57147c6ae99SBarry Smith n0 = rank + m * (n-1) -1 - (m*n); 57247c6ae99SBarry Smith n1 = rank + m * (n-1) - (m*n); 57347c6ae99SBarry Smith n2 = rank + m * (n-1) +1 - (m*n); 57447c6ae99SBarry Smith n9 = rank + m * (n-1) -1; 57547c6ae99SBarry Smith n10 = rank + m * (n-1); 57647c6ae99SBarry Smith n11 = rank + m * (n-1) +1; 57747c6ae99SBarry Smith n18 = rank + m * (n-1) -1 + (m*n); 57847c6ae99SBarry Smith n19 = rank + m * (n-1) + (m*n); 57947c6ae99SBarry Smith n20 = rank + m * (n-1) +1 + (m*n); 58047c6ae99SBarry Smith } 58147c6ae99SBarry Smith 58247c6ae99SBarry Smith if (ye == N) { /* First assume not corner or edge */ 58347c6ae99SBarry Smith n6 = rank - m * (n-1) -1 - (m*n); 58447c6ae99SBarry Smith n7 = rank - m * (n-1) - (m*n); 58547c6ae99SBarry Smith n8 = rank - m * (n-1) +1 - (m*n); 58647c6ae99SBarry Smith n15 = rank - m * (n-1) -1; 58747c6ae99SBarry Smith n16 = rank - m * (n-1); 58847c6ae99SBarry Smith n17 = rank - m * (n-1) +1; 58947c6ae99SBarry Smith n24 = rank - m * (n-1) -1 + (m*n); 59047c6ae99SBarry Smith n25 = rank - m * (n-1) + (m*n); 59147c6ae99SBarry Smith n26 = rank - m * (n-1) +1 + (m*n); 59247c6ae99SBarry Smith } 59347c6ae99SBarry Smith 59447c6ae99SBarry Smith if (zs == 0) { /* First assume not corner or edge */ 59547c6ae99SBarry Smith n0 = size - (m*n) + rank - m - 1; 59647c6ae99SBarry Smith n1 = size - (m*n) + rank - m; 59747c6ae99SBarry Smith n2 = size - (m*n) + rank - m + 1; 59847c6ae99SBarry Smith n3 = size - (m*n) + rank - 1; 59947c6ae99SBarry Smith n4 = size - (m*n) + rank; 60047c6ae99SBarry Smith n5 = size - (m*n) + rank + 1; 60147c6ae99SBarry Smith n6 = size - (m*n) + rank + m - 1; 60247c6ae99SBarry Smith n7 = size - (m*n) + rank + m; 60347c6ae99SBarry Smith n8 = size - (m*n) + rank + m + 1; 60447c6ae99SBarry Smith } 60547c6ae99SBarry Smith 60647c6ae99SBarry Smith if (ze == P) { /* First assume not corner or edge */ 60747c6ae99SBarry Smith n18 = (m*n) - (size-rank) - m - 1; 60847c6ae99SBarry Smith n19 = (m*n) - (size-rank) - m; 60947c6ae99SBarry Smith n20 = (m*n) - (size-rank) - m + 1; 61047c6ae99SBarry Smith n21 = (m*n) - (size-rank) - 1; 61147c6ae99SBarry Smith n22 = (m*n) - (size-rank); 61247c6ae99SBarry Smith n23 = (m*n) - (size-rank) + 1; 61347c6ae99SBarry Smith n24 = (m*n) - (size-rank) + m - 1; 61447c6ae99SBarry Smith n25 = (m*n) - (size-rank) + m; 61547c6ae99SBarry Smith n26 = (m*n) - (size-rank) + m + 1; 61647c6ae99SBarry Smith } 61747c6ae99SBarry Smith 61847c6ae99SBarry Smith if ((xs==0) && (zs==0)) { /* Assume an edge, not corner */ 61947c6ae99SBarry Smith n0 = size - m*n + rank + m-1 - m; 62047c6ae99SBarry Smith n3 = size - m*n + rank + m-1; 62147c6ae99SBarry Smith n6 = size - m*n + rank + m-1 + m; 62247c6ae99SBarry Smith } 62347c6ae99SBarry Smith 62447c6ae99SBarry Smith if ((xs==0) && (ze==P)) { /* Assume an edge, not corner */ 62547c6ae99SBarry Smith n18 = m*n - (size - rank) + m-1 - m; 62647c6ae99SBarry Smith n21 = m*n - (size - rank) + m-1; 62747c6ae99SBarry Smith n24 = m*n - (size - rank) + m-1 + m; 62847c6ae99SBarry Smith } 62947c6ae99SBarry Smith 63047c6ae99SBarry Smith if ((xs==0) && (ys==0)) { /* Assume an edge, not corner */ 63147c6ae99SBarry Smith n0 = rank + m*n -1 - m*n; 63247c6ae99SBarry Smith n9 = rank + m*n -1; 63347c6ae99SBarry Smith n18 = rank + m*n -1 + m*n; 63447c6ae99SBarry Smith } 63547c6ae99SBarry Smith 63647c6ae99SBarry Smith if ((xs==0) && (ye==N)) { /* Assume an edge, not corner */ 63747c6ae99SBarry Smith n6 = rank - m*(n-1) + m-1 - m*n; 63847c6ae99SBarry Smith n15 = rank - m*(n-1) + m-1; 63947c6ae99SBarry Smith n24 = rank - m*(n-1) + m-1 + m*n; 64047c6ae99SBarry Smith } 64147c6ae99SBarry Smith 642ce00eea3SSatish Balay if ((xe==M) && (zs==0)) { /* Assume an edge, not corner */ 64347c6ae99SBarry Smith n2 = size - (m*n-rank) - (m-1) - m; 64447c6ae99SBarry Smith n5 = size - (m*n-rank) - (m-1); 64547c6ae99SBarry Smith n8 = size - (m*n-rank) - (m-1) + m; 64647c6ae99SBarry Smith } 64747c6ae99SBarry Smith 648ce00eea3SSatish Balay if ((xe==M) && (ze==P)) { /* Assume an edge, not corner */ 64947c6ae99SBarry Smith n20 = m*n - (size - rank) - (m-1) - m; 65047c6ae99SBarry Smith n23 = m*n - (size - rank) - (m-1); 65147c6ae99SBarry Smith n26 = m*n - (size - rank) - (m-1) + m; 65247c6ae99SBarry Smith } 65347c6ae99SBarry Smith 654ce00eea3SSatish Balay if ((xe==M) && (ys==0)) { /* Assume an edge, not corner */ 65547c6ae99SBarry Smith n2 = rank + m*(n-1) - (m-1) - m*n; 65647c6ae99SBarry Smith n11 = rank + m*(n-1) - (m-1); 65747c6ae99SBarry Smith n20 = rank + m*(n-1) - (m-1) + m*n; 65847c6ae99SBarry Smith } 65947c6ae99SBarry Smith 660ce00eea3SSatish Balay if ((xe==M) && (ye==N)) { /* Assume an edge, not corner */ 66147c6ae99SBarry Smith n8 = rank - m*n +1 - m*n; 66247c6ae99SBarry Smith n17 = rank - m*n +1; 66347c6ae99SBarry Smith n26 = rank - m*n +1 + m*n; 66447c6ae99SBarry Smith } 66547c6ae99SBarry Smith 66647c6ae99SBarry Smith if ((ys==0) && (zs==0)) { /* Assume an edge, not corner */ 66747c6ae99SBarry Smith n0 = size - m + rank -1; 66847c6ae99SBarry Smith n1 = size - m + rank; 66947c6ae99SBarry Smith n2 = size - m + rank +1; 67047c6ae99SBarry Smith } 67147c6ae99SBarry Smith 67247c6ae99SBarry Smith if ((ys==0) && (ze==P)) { /* Assume an edge, not corner */ 67347c6ae99SBarry Smith n18 = m*n - (size - rank) + m*(n-1) -1; 67447c6ae99SBarry Smith n19 = m*n - (size - rank) + m*(n-1); 67547c6ae99SBarry Smith n20 = m*n - (size - rank) + m*(n-1) +1; 67647c6ae99SBarry Smith } 67747c6ae99SBarry Smith 67847c6ae99SBarry Smith if ((ye==N) && (zs==0)) { /* Assume an edge, not corner */ 67947c6ae99SBarry Smith n6 = size - (m*n-rank) - m * (n-1) -1; 68047c6ae99SBarry Smith n7 = size - (m*n-rank) - m * (n-1); 68147c6ae99SBarry Smith n8 = size - (m*n-rank) - m * (n-1) +1; 68247c6ae99SBarry Smith } 68347c6ae99SBarry Smith 68447c6ae99SBarry Smith if ((ye==N) && (ze==P)) { /* Assume an edge, not corner */ 68547c6ae99SBarry Smith n24 = rank - (size-m) -1; 68647c6ae99SBarry Smith n25 = rank - (size-m); 68747c6ae99SBarry Smith n26 = rank - (size-m) +1; 68847c6ae99SBarry Smith } 68947c6ae99SBarry Smith 69047c6ae99SBarry Smith /* Check for Corners */ 6918865f1eaSKarl Rupp if ((xs==0) && (ys==0) && (zs==0)) n0 = size -1; 6928865f1eaSKarl Rupp if ((xs==0) && (ys==0) && (ze==P)) n18 = m*n-1; 6938865f1eaSKarl Rupp if ((xs==0) && (ye==N) && (zs==0)) n6 = (size-1)-m*(n-1); 6948865f1eaSKarl Rupp if ((xs==0) && (ye==N) && (ze==P)) n24 = m-1; 6958865f1eaSKarl Rupp if ((xe==M) && (ys==0) && (zs==0)) n2 = size-m; 6968865f1eaSKarl Rupp if ((xe==M) && (ys==0) && (ze==P)) n20 = m*n-m; 6978865f1eaSKarl Rupp if ((xe==M) && (ye==N) && (zs==0)) n8 = size-m*n; 6988865f1eaSKarl Rupp if ((xe==M) && (ye==N) && (ze==P)) n26 = 0; 69947c6ae99SBarry Smith 70047c6ae99SBarry Smith /* Check for when not X,Y, and Z Periodic */ 70147c6ae99SBarry Smith 70247c6ae99SBarry Smith /* If not X periodic */ 703bff4a2f0SMatthew G. Knepley if (bx != DM_BOUNDARY_PERIODIC) { 7048865f1eaSKarl Rupp if (xs==0) n0 = n3 = n6 = n9 = n12 = n15 = n18 = n21 = n24 = -2; 7058865f1eaSKarl Rupp if (xe==M) n2 = n5 = n8 = n11 = n14 = n17 = n20 = n23 = n26 = -2; 70647c6ae99SBarry Smith } 70747c6ae99SBarry Smith 70847c6ae99SBarry Smith /* If not Y periodic */ 709bff4a2f0SMatthew G. Knepley if (by != DM_BOUNDARY_PERIODIC) { 7108865f1eaSKarl Rupp if (ys==0) n0 = n1 = n2 = n9 = n10 = n11 = n18 = n19 = n20 = -2; 7118865f1eaSKarl Rupp if (ye==N) n6 = n7 = n8 = n15 = n16 = n17 = n24 = n25 = n26 = -2; 71247c6ae99SBarry Smith } 71347c6ae99SBarry Smith 71447c6ae99SBarry Smith /* If not Z periodic */ 715bff4a2f0SMatthew G. Knepley if (bz != DM_BOUNDARY_PERIODIC) { 7168865f1eaSKarl Rupp if (zs==0) n0 = n1 = n2 = n3 = n4 = n5 = n6 = n7 = n8 = -2; 7178865f1eaSKarl Rupp if (ze==P) n18 = n19 = n20 = n21 = n22 = n23 = n24 = n25 = n26 = -2; 71847c6ae99SBarry Smith } 71947c6ae99SBarry Smith 720785e854fSJed Brown ierr = PetscMalloc1(27,&dd->neighbors);CHKERRQ(ierr); 7218865f1eaSKarl Rupp 72247c6ae99SBarry Smith dd->neighbors[0] = n0; 72347c6ae99SBarry Smith dd->neighbors[1] = n1; 72447c6ae99SBarry Smith dd->neighbors[2] = n2; 72547c6ae99SBarry Smith dd->neighbors[3] = n3; 72647c6ae99SBarry Smith dd->neighbors[4] = n4; 72747c6ae99SBarry Smith dd->neighbors[5] = n5; 72847c6ae99SBarry Smith dd->neighbors[6] = n6; 72947c6ae99SBarry Smith dd->neighbors[7] = n7; 73047c6ae99SBarry Smith dd->neighbors[8] = n8; 73147c6ae99SBarry Smith dd->neighbors[9] = n9; 73247c6ae99SBarry Smith dd->neighbors[10] = n10; 73347c6ae99SBarry Smith dd->neighbors[11] = n11; 73447c6ae99SBarry Smith dd->neighbors[12] = n12; 73547c6ae99SBarry Smith dd->neighbors[13] = rank; 73647c6ae99SBarry Smith dd->neighbors[14] = n14; 73747c6ae99SBarry Smith dd->neighbors[15] = n15; 73847c6ae99SBarry Smith dd->neighbors[16] = n16; 73947c6ae99SBarry Smith dd->neighbors[17] = n17; 74047c6ae99SBarry Smith dd->neighbors[18] = n18; 74147c6ae99SBarry Smith dd->neighbors[19] = n19; 74247c6ae99SBarry Smith dd->neighbors[20] = n20; 74347c6ae99SBarry Smith dd->neighbors[21] = n21; 74447c6ae99SBarry Smith dd->neighbors[22] = n22; 74547c6ae99SBarry Smith dd->neighbors[23] = n23; 74647c6ae99SBarry Smith dd->neighbors[24] = n24; 74747c6ae99SBarry Smith dd->neighbors[25] = n25; 74847c6ae99SBarry Smith dd->neighbors[26] = n26; 74947c6ae99SBarry Smith 75047c6ae99SBarry Smith /* If star stencil then delete the corner neighbors */ 751d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 75247c6ae99SBarry Smith /* save information about corner neighbors */ 75347c6ae99SBarry Smith sn0 = n0; sn1 = n1; sn2 = n2; sn3 = n3; sn5 = n5; sn6 = n6; sn7 = n7; 75447c6ae99SBarry Smith sn8 = n8; sn9 = n9; sn11 = n11; sn15 = n15; sn17 = n17; sn18 = n18; 75547c6ae99SBarry Smith sn19 = n19; sn20 = n20; sn21 = n21; sn23 = n23; sn24 = n24; sn25 = n25; 75647c6ae99SBarry Smith sn26 = n26; 7578865f1eaSKarl Rupp n0 = n1 = n2 = n3 = n5 = n6 = n7 = n8 = n9 = n11 = n15 = n17 = n18 = n19 = n20 = n21 = n23 = n24 = n25 = n26 = -1; 75847c6ae99SBarry Smith } 75947c6ae99SBarry Smith 760785e854fSJed Brown ierr = PetscMalloc1((Xe-Xs)*(Ye-Ys)*(Ze-Zs),&idx);CHKERRQ(ierr); 76147c6ae99SBarry Smith 76247c6ae99SBarry Smith nn = 0; 76347c6ae99SBarry Smith /* Bottom Level */ 76447c6ae99SBarry Smith for (k=0; k<s_z; k++) { 76547c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 76647c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 767ce00eea3SSatish Balay x_t = lx[n0 % m]; 76847c6ae99SBarry Smith y_t = ly[(n0 % (m*n))/m]; 76947c6ae99SBarry Smith z_t = lz[n0 / (m*n)]; 77047c6ae99SBarry Smith s_t = bases[n0] + x_t*y_t*z_t - (s_y-i)*x_t - s_x - (s_z-k-1)*x_t*y_t; 7718865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n0] + x_t*y_t*z_t - (s_y-i)*x_t - s_x; /* 2D case */ 7728865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 77347c6ae99SBarry Smith } 77447c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 77547c6ae99SBarry Smith x_t = x; 77647c6ae99SBarry Smith y_t = ly[(n1 % (m*n))/m]; 77747c6ae99SBarry Smith z_t = lz[n1 / (m*n)]; 77847c6ae99SBarry Smith s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t; 7798865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t; /* 2D case */ 7808865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 78147c6ae99SBarry Smith } 78247c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 783ce00eea3SSatish Balay x_t = lx[n2 % m]; 78447c6ae99SBarry Smith y_t = ly[(n2 % (m*n))/m]; 78547c6ae99SBarry Smith z_t = lz[n2 / (m*n)]; 78647c6ae99SBarry Smith s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t; 7878865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t; /* 2D case */ 7888865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 78947c6ae99SBarry Smith } 79047c6ae99SBarry Smith } 79147c6ae99SBarry Smith 79247c6ae99SBarry Smith for (i=0; i<y; i++) { 79347c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 794ce00eea3SSatish Balay x_t = lx[n3 % m]; 79547c6ae99SBarry Smith y_t = y; 79647c6ae99SBarry Smith z_t = lz[n3 / (m*n)]; 79747c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 7988865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - x_t*y_t; /* 2D case */ 7998865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 80047c6ae99SBarry Smith } 80147c6ae99SBarry Smith 80247c6ae99SBarry Smith if (n4 >= 0) { /* middle */ 80347c6ae99SBarry Smith x_t = x; 80447c6ae99SBarry Smith y_t = y; 80547c6ae99SBarry Smith z_t = lz[n4 / (m*n)]; 80647c6ae99SBarry Smith s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8078865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n4] + i*x_t + x_t*y_t*z_t - x_t*y_t; /* 2D case */ 8088865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 809bff4a2f0SMatthew G. Knepley } else if (bz == DM_BOUNDARY_MIRROR) { 8108865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 81147c6ae99SBarry Smith } 81247c6ae99SBarry Smith 81347c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 814ce00eea3SSatish Balay x_t = lx[n5 % m]; 81547c6ae99SBarry Smith y_t = y; 81647c6ae99SBarry Smith z_t = lz[n5 / (m*n)]; 81747c6ae99SBarry Smith s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8188865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n5] + i*x_t + x_t*y_t*z_t - x_t*y_t; /* 2D case */ 8198865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 82047c6ae99SBarry Smith } 82147c6ae99SBarry Smith } 82247c6ae99SBarry Smith 82347c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 82447c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 825ce00eea3SSatish Balay x_t = lx[n6 % m]; 82647c6ae99SBarry Smith y_t = ly[(n6 % (m*n))/m]; 82747c6ae99SBarry Smith z_t = lz[n6 / (m*n)]; 82847c6ae99SBarry Smith s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8298865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - x_t*y_t; /* 2D case */ 8308865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 83147c6ae99SBarry Smith } 83247c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 83347c6ae99SBarry Smith x_t = x; 83447c6ae99SBarry Smith y_t = ly[(n7 % (m*n))/m]; 83547c6ae99SBarry Smith z_t = lz[n7 / (m*n)]; 83647c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8378865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - x_t*y_t; /* 2D case */ 8388865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 83947c6ae99SBarry Smith } 84047c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 841ce00eea3SSatish Balay x_t = lx[n8 % m]; 84247c6ae99SBarry Smith y_t = ly[(n8 % (m*n))/m]; 84347c6ae99SBarry Smith z_t = lz[n8 / (m*n)]; 84447c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8458865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - x_t*y_t; /* 2D case */ 8468865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 84747c6ae99SBarry Smith } 84847c6ae99SBarry Smith } 84947c6ae99SBarry Smith } 85047c6ae99SBarry Smith 85147c6ae99SBarry Smith /* Middle Level */ 85247c6ae99SBarry Smith for (k=0; k<z; k++) { 85347c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 85447c6ae99SBarry Smith if (n9 >= 0) { /* left below */ 855ce00eea3SSatish Balay x_t = lx[n9 % m]; 85647c6ae99SBarry Smith y_t = ly[(n9 % (m*n))/m]; 85747c6ae99SBarry Smith /* z_t = z; */ 85847c6ae99SBarry Smith s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 8598865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 86047c6ae99SBarry Smith } 86147c6ae99SBarry Smith if (n10 >= 0) { /* directly below */ 86247c6ae99SBarry Smith x_t = x; 86347c6ae99SBarry Smith y_t = ly[(n10 % (m*n))/m]; 86447c6ae99SBarry Smith /* z_t = z; */ 86547c6ae99SBarry Smith s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 8668865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 867bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_MIRROR) { 8688865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 86947c6ae99SBarry Smith } 87047c6ae99SBarry Smith if (n11 >= 0) { /* right below */ 871ce00eea3SSatish Balay x_t = lx[n11 % m]; 87247c6ae99SBarry Smith y_t = ly[(n11 % (m*n))/m]; 87347c6ae99SBarry Smith /* z_t = z; */ 87447c6ae99SBarry Smith s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 8758865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 87647c6ae99SBarry Smith } 87747c6ae99SBarry Smith } 87847c6ae99SBarry Smith 87947c6ae99SBarry Smith for (i=0; i<y; i++) { 88047c6ae99SBarry Smith if (n12 >= 0) { /* directly left */ 881ce00eea3SSatish Balay x_t = lx[n12 % m]; 88247c6ae99SBarry Smith y_t = y; 88347c6ae99SBarry Smith /* z_t = z; */ 88447c6ae99SBarry Smith s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t; 8858865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 886bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { 8878865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = 0; 88847c6ae99SBarry Smith } 88947c6ae99SBarry Smith 89047c6ae99SBarry Smith /* Interior */ 89147c6ae99SBarry Smith s_t = bases[rank] + i*x + k*x*y; 8928865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = s_t++; 89347c6ae99SBarry Smith 89447c6ae99SBarry Smith if (n14 >= 0) { /* directly right */ 895ce00eea3SSatish Balay x_t = lx[n14 % m]; 89647c6ae99SBarry Smith y_t = y; 89747c6ae99SBarry Smith /* z_t = z; */ 89847c6ae99SBarry Smith s_t = bases[n14] + i*x_t + k*x_t*y_t; 8998865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 900bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { 9018865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = 0; 90247c6ae99SBarry Smith } 90347c6ae99SBarry Smith } 90447c6ae99SBarry Smith 90547c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 90647c6ae99SBarry Smith if (n15 >= 0) { /* left above */ 907ce00eea3SSatish Balay x_t = lx[n15 % m]; 90847c6ae99SBarry Smith y_t = ly[(n15 % (m*n))/m]; 90947c6ae99SBarry Smith /* z_t = z; */ 91047c6ae99SBarry Smith s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t; 9118865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 91247c6ae99SBarry Smith } 91347c6ae99SBarry Smith if (n16 >= 0) { /* directly above */ 91447c6ae99SBarry Smith x_t = x; 91547c6ae99SBarry Smith y_t = ly[(n16 % (m*n))/m]; 91647c6ae99SBarry Smith /* z_t = z; */ 91747c6ae99SBarry Smith s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t; 9188865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 919bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_MIRROR) { 9208865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 92147c6ae99SBarry Smith } 92247c6ae99SBarry Smith if (n17 >= 0) { /* right above */ 923ce00eea3SSatish Balay x_t = lx[n17 % m]; 92447c6ae99SBarry Smith y_t = ly[(n17 % (m*n))/m]; 92547c6ae99SBarry Smith /* z_t = z; */ 92647c6ae99SBarry Smith s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t; 9278865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 92847c6ae99SBarry Smith } 92947c6ae99SBarry Smith } 93047c6ae99SBarry Smith } 93147c6ae99SBarry Smith 93247c6ae99SBarry Smith /* Upper Level */ 93347c6ae99SBarry Smith for (k=0; k<s_z; k++) { 93447c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 93547c6ae99SBarry Smith if (n18 >= 0) { /* left below */ 936ce00eea3SSatish Balay x_t = lx[n18 % m]; 93747c6ae99SBarry Smith y_t = ly[(n18 % (m*n))/m]; 93847c6ae99SBarry Smith /* z_t = lz[n18 / (m*n)]; */ 93947c6ae99SBarry Smith s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 9408865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n18] - (s_y-i)*x_t -s_x + x_t*y_t; /* 2d case */ 9418865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 94247c6ae99SBarry Smith } 94347c6ae99SBarry Smith if (n19 >= 0) { /* directly below */ 94447c6ae99SBarry Smith x_t = x; 94547c6ae99SBarry Smith y_t = ly[(n19 % (m*n))/m]; 94647c6ae99SBarry Smith /* z_t = lz[n19 / (m*n)]; */ 94747c6ae99SBarry Smith s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 9488865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n19] - (s_y+1-i)*x_t + x_t*y_t; /* 2d case */ 9498865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 95047c6ae99SBarry Smith } 95147c6ae99SBarry Smith if (n20 >= 0) { /* right below */ 952ce00eea3SSatish Balay x_t = lx[n20 % m]; 95347c6ae99SBarry Smith y_t = ly[(n20 % (m*n))/m]; 95447c6ae99SBarry Smith /* z_t = lz[n20 / (m*n)]; */ 95547c6ae99SBarry Smith s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 9568865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n20] - (s_y+1-i)*x_t + x_t*y_t; /* 2d case */ 9578865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 95847c6ae99SBarry Smith } 95947c6ae99SBarry Smith } 96047c6ae99SBarry Smith 96147c6ae99SBarry Smith for (i=0; i<y; i++) { 96247c6ae99SBarry Smith if (n21 >= 0) { /* directly left */ 963ce00eea3SSatish Balay x_t = lx[n21 % m]; 96447c6ae99SBarry Smith y_t = y; 96547c6ae99SBarry Smith /* z_t = lz[n21 / (m*n)]; */ 96647c6ae99SBarry Smith s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t; 9678865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n21] + (i+1)*x_t - s_x; /* 2d case */ 9688865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 96947c6ae99SBarry Smith } 97047c6ae99SBarry Smith 97147c6ae99SBarry Smith if (n22 >= 0) { /* middle */ 97247c6ae99SBarry Smith x_t = x; 97347c6ae99SBarry Smith y_t = y; 97447c6ae99SBarry Smith /* z_t = lz[n22 / (m*n)]; */ 97547c6ae99SBarry Smith s_t = bases[n22] + i*x_t + k*x_t*y_t; 9768865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n22] + i*x_t; /* 2d case */ 9778865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 978bff4a2f0SMatthew G. Knepley } else if (bz == DM_BOUNDARY_MIRROR) { 9798865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 98047c6ae99SBarry Smith } 98147c6ae99SBarry Smith 98247c6ae99SBarry Smith if (n23 >= 0) { /* directly right */ 983ce00eea3SSatish Balay x_t = lx[n23 % m]; 98447c6ae99SBarry Smith y_t = y; 98547c6ae99SBarry Smith /* z_t = lz[n23 / (m*n)]; */ 98647c6ae99SBarry Smith s_t = bases[n23] + i*x_t + k*x_t*y_t; 9878865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n23] + i*x_t; /* 2d case */ 9888865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 98947c6ae99SBarry Smith } 99047c6ae99SBarry Smith } 99147c6ae99SBarry Smith 99247c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 99347c6ae99SBarry Smith if (n24 >= 0) { /* left above */ 994ce00eea3SSatish Balay x_t = lx[n24 % m]; 99547c6ae99SBarry Smith y_t = ly[(n24 % (m*n))/m]; 99647c6ae99SBarry Smith /* z_t = lz[n24 / (m*n)]; */ 99747c6ae99SBarry Smith s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t; 9988865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n24] + i*x_t - s_x; /* 2d case */ 9998865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 100047c6ae99SBarry Smith } 100147c6ae99SBarry Smith if (n25 >= 0) { /* directly above */ 100247c6ae99SBarry Smith x_t = x; 100347c6ae99SBarry Smith y_t = ly[(n25 % (m*n))/m]; 100447c6ae99SBarry Smith /* z_t = lz[n25 / (m*n)]; */ 100547c6ae99SBarry Smith s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t; 10068865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n25] + (i-1)*x_t; /* 2d case */ 10078865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 100847c6ae99SBarry Smith } 100947c6ae99SBarry Smith if (n26 >= 0) { /* right above */ 1010ce00eea3SSatish Balay x_t = lx[n26 % m]; 101147c6ae99SBarry Smith y_t = ly[(n26 % (m*n))/m]; 101247c6ae99SBarry Smith /* z_t = lz[n26 / (m*n)]; */ 101347c6ae99SBarry Smith s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t; 10148865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n26] + (i-1)*x_t; /* 2d case */ 10158865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 101647c6ae99SBarry Smith } 101747c6ae99SBarry Smith } 101847c6ae99SBarry Smith } 1019ce00eea3SSatish Balay 1020b1fb7eb7SBarry Smith ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_USE_POINTER,&from);CHKERRQ(ierr); 102147c6ae99SBarry Smith ierr = VecScatterCreate(global,from,local,to,>ol);CHKERRQ(ierr); 10223bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)gtol);CHKERRQ(ierr); 1023fcfd50ebSBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 1024fcfd50ebSBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 102547c6ae99SBarry Smith 1026d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 102747c6ae99SBarry Smith n0 = sn0; n1 = sn1; n2 = sn2; n3 = sn3; n5 = sn5; n6 = sn6; n7 = sn7; 102847c6ae99SBarry Smith n8 = sn8; n9 = sn9; n11 = sn11; n15 = sn15; n17 = sn17; n18 = sn18; 102947c6ae99SBarry Smith n19 = sn19; n20 = sn20; n21 = sn21; n23 = sn23; n24 = sn24; n25 = sn25; 103047c6ae99SBarry Smith n26 = sn26; 1031ce00eea3SSatish Balay } 103247c6ae99SBarry Smith 103388661749SPeter Brune if (((stencil_type == DMDA_STENCIL_STAR) || 1034bff4a2f0SMatthew G. Knepley (bx != DM_BOUNDARY_PERIODIC && bx) || 1035bff4a2f0SMatthew G. Knepley (by != DM_BOUNDARY_PERIODIC && by) || 1036bff4a2f0SMatthew G. Knepley (bz != DM_BOUNDARY_PERIODIC && bz))) { 1037ce00eea3SSatish Balay /* 1038ce00eea3SSatish Balay Recompute the local to global mappings, this time keeping the 1039ce00eea3SSatish Balay information about the cross corner processor numbers. 1040ce00eea3SSatish Balay */ 104147c6ae99SBarry Smith nn = 0; 104247c6ae99SBarry Smith /* Bottom Level */ 104347c6ae99SBarry Smith for (k=0; k<s_z; k++) { 104447c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 104547c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 1046ce00eea3SSatish Balay x_t = lx[n0 % m]; 104747c6ae99SBarry Smith y_t = ly[(n0 % (m*n))/m]; 104847c6ae99SBarry Smith z_t = lz[n0 / (m*n)]; 104947c6ae99SBarry Smith s_t = bases[n0] + x_t*y_t*z_t - (s_y-i)*x_t - s_x - (s_z-k-1)*x_t*y_t; 10508865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1051ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0 && Zs-zs < 0) { 10528865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 105347c6ae99SBarry Smith } 105447c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 105547c6ae99SBarry Smith x_t = x; 105647c6ae99SBarry Smith y_t = ly[(n1 % (m*n))/m]; 105747c6ae99SBarry Smith z_t = lz[n1 / (m*n)]; 105847c6ae99SBarry Smith s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t; 10598865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1060ce00eea3SSatish Balay } else if (Ys-ys < 0 && Zs-zs < 0) { 10618865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 106247c6ae99SBarry Smith } 106347c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 1064ce00eea3SSatish Balay x_t = lx[n2 % m]; 106547c6ae99SBarry Smith y_t = ly[(n2 % (m*n))/m]; 106647c6ae99SBarry Smith z_t = lz[n2 / (m*n)]; 106747c6ae99SBarry Smith s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t; 10688865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1069ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0 && Zs-zs < 0) { 10708865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 107147c6ae99SBarry Smith } 107247c6ae99SBarry Smith } 107347c6ae99SBarry Smith 107447c6ae99SBarry Smith for (i=0; i<y; i++) { 107547c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 1076ce00eea3SSatish Balay x_t = lx[n3 % m]; 107747c6ae99SBarry Smith y_t = y; 107847c6ae99SBarry Smith z_t = lz[n3 / (m*n)]; 107947c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 10808865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1081ce00eea3SSatish Balay } else if (Xs-xs < 0 && Zs-zs < 0) { 10828865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 108347c6ae99SBarry Smith } 108447c6ae99SBarry Smith 108547c6ae99SBarry Smith if (n4 >= 0) { /* middle */ 108647c6ae99SBarry Smith x_t = x; 108747c6ae99SBarry Smith y_t = y; 108847c6ae99SBarry Smith z_t = lz[n4 / (m*n)]; 108947c6ae99SBarry Smith s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 10908865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1091ce00eea3SSatish Balay } else if (Zs-zs < 0) { 1092bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_MIRROR) { 10938865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 1094c2859e5eSBarry Smith } else { 10958865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 109647c6ae99SBarry Smith } 1097c2859e5eSBarry Smith } 109847c6ae99SBarry Smith 109947c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 1100ce00eea3SSatish Balay x_t = lx[n5 % m]; 110147c6ae99SBarry Smith y_t = y; 110247c6ae99SBarry Smith z_t = lz[n5 / (m*n)]; 110347c6ae99SBarry Smith s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11048865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1105ce00eea3SSatish Balay } else if (xe-Xe < 0 && Zs-zs < 0) { 11068865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 110747c6ae99SBarry Smith } 110847c6ae99SBarry Smith } 110947c6ae99SBarry Smith 111047c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 111147c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 1112ce00eea3SSatish Balay x_t = lx[n6 % m]; 111347c6ae99SBarry Smith y_t = ly[(n6 % (m*n))/m]; 111447c6ae99SBarry Smith z_t = lz[n6 / (m*n)]; 111547c6ae99SBarry Smith s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11168865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1117ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0 && Zs-zs < 0) { 11188865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 111947c6ae99SBarry Smith } 112047c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 112147c6ae99SBarry Smith x_t = x; 112247c6ae99SBarry Smith y_t = ly[(n7 % (m*n))/m]; 112347c6ae99SBarry Smith z_t = lz[n7 / (m*n)]; 112447c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11258865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1126ce00eea3SSatish Balay } else if (ye-Ye < 0 && Zs-zs < 0) { 11278865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 112847c6ae99SBarry Smith } 112947c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 1130ce00eea3SSatish Balay x_t = lx[n8 % m]; 113147c6ae99SBarry Smith y_t = ly[(n8 % (m*n))/m]; 113247c6ae99SBarry Smith z_t = lz[n8 / (m*n)]; 113347c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11348865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1135ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0 && Zs-zs < 0) { 11368865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 113747c6ae99SBarry Smith } 113847c6ae99SBarry Smith } 113947c6ae99SBarry Smith } 114047c6ae99SBarry Smith 114147c6ae99SBarry Smith /* Middle Level */ 114247c6ae99SBarry Smith for (k=0; k<z; k++) { 114347c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 114447c6ae99SBarry Smith if (n9 >= 0) { /* left below */ 1145ce00eea3SSatish Balay x_t = lx[n9 % m]; 114647c6ae99SBarry Smith y_t = ly[(n9 % (m*n))/m]; 114747c6ae99SBarry Smith /* z_t = z; */ 114847c6ae99SBarry Smith s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 11498865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1150ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0) { 11518865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 115247c6ae99SBarry Smith } 115347c6ae99SBarry Smith if (n10 >= 0) { /* directly below */ 115447c6ae99SBarry Smith x_t = x; 115547c6ae99SBarry Smith y_t = ly[(n10 % (m*n))/m]; 115647c6ae99SBarry Smith /* z_t = z; */ 115747c6ae99SBarry Smith s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 11588865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1159ce00eea3SSatish Balay } else if (Ys-ys < 0) { 1160bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_MIRROR) { 11618865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 1162c2859e5eSBarry Smith } else { 11638865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 1164c2859e5eSBarry Smith } 116547c6ae99SBarry Smith } 116647c6ae99SBarry Smith if (n11 >= 0) { /* right below */ 1167ce00eea3SSatish Balay x_t = lx[n11 % m]; 116847c6ae99SBarry Smith y_t = ly[(n11 % (m*n))/m]; 116947c6ae99SBarry Smith /* z_t = z; */ 117047c6ae99SBarry Smith s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 11718865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1172ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0) { 11738865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 117447c6ae99SBarry Smith } 117547c6ae99SBarry Smith } 117647c6ae99SBarry Smith 117747c6ae99SBarry Smith for (i=0; i<y; i++) { 117847c6ae99SBarry Smith if (n12 >= 0) { /* directly left */ 1179ce00eea3SSatish Balay x_t = lx[n12 % m]; 118047c6ae99SBarry Smith y_t = y; 118147c6ae99SBarry Smith /* z_t = z; */ 118247c6ae99SBarry Smith s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t; 11838865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1184ce00eea3SSatish Balay } else if (Xs-xs < 0) { 1185bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_MIRROR) { 11868865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = 0; 1187c2859e5eSBarry Smith } else { 11888865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 118947c6ae99SBarry Smith } 1190c2859e5eSBarry Smith } 119147c6ae99SBarry Smith 119247c6ae99SBarry Smith /* Interior */ 119347c6ae99SBarry Smith s_t = bases[rank] + i*x + k*x*y; 11948865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = s_t++; 119547c6ae99SBarry Smith 119647c6ae99SBarry Smith if (n14 >= 0) { /* directly right */ 1197ce00eea3SSatish Balay x_t = lx[n14 % m]; 119847c6ae99SBarry Smith y_t = y; 119947c6ae99SBarry Smith /* z_t = z; */ 120047c6ae99SBarry Smith s_t = bases[n14] + i*x_t + k*x_t*y_t; 12018865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1202ce00eea3SSatish Balay } else if (xe-Xe < 0) { 1203bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_MIRROR) { 12048865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = 0; 1205c2859e5eSBarry Smith } else { 12068865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 120747c6ae99SBarry Smith } 120847c6ae99SBarry Smith } 1209c2859e5eSBarry Smith } 121047c6ae99SBarry Smith 121147c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 121247c6ae99SBarry Smith if (n15 >= 0) { /* left above */ 1213ce00eea3SSatish Balay x_t = lx[n15 % m]; 121447c6ae99SBarry Smith y_t = ly[(n15 % (m*n))/m]; 121547c6ae99SBarry Smith /* z_t = z; */ 121647c6ae99SBarry Smith s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t; 12178865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1218ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0) { 12198865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 122047c6ae99SBarry Smith } 122147c6ae99SBarry Smith if (n16 >= 0) { /* directly above */ 122247c6ae99SBarry Smith x_t = x; 122347c6ae99SBarry Smith y_t = ly[(n16 % (m*n))/m]; 122447c6ae99SBarry Smith /* z_t = z; */ 122547c6ae99SBarry Smith s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t; 12268865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1227ce00eea3SSatish Balay } else if (ye-Ye < 0) { 1228bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_MIRROR) { 12298865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 1230c2859e5eSBarry Smith } else { 12318865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 123247c6ae99SBarry Smith } 1233c2859e5eSBarry Smith } 123447c6ae99SBarry Smith if (n17 >= 0) { /* right above */ 1235ce00eea3SSatish Balay x_t = lx[n17 % m]; 123647c6ae99SBarry Smith y_t = ly[(n17 % (m*n))/m]; 123747c6ae99SBarry Smith /* z_t = z; */ 123847c6ae99SBarry Smith s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t; 12398865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1240ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0) { 12418865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 124247c6ae99SBarry Smith } 124347c6ae99SBarry Smith } 124447c6ae99SBarry Smith } 124547c6ae99SBarry Smith 124647c6ae99SBarry Smith /* Upper Level */ 124747c6ae99SBarry Smith for (k=0; k<s_z; k++) { 124847c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 124947c6ae99SBarry Smith if (n18 >= 0) { /* left below */ 1250ce00eea3SSatish Balay x_t = lx[n18 % m]; 125147c6ae99SBarry Smith y_t = ly[(n18 % (m*n))/m]; 125247c6ae99SBarry Smith /* z_t = lz[n18 / (m*n)]; */ 125347c6ae99SBarry Smith s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 12548865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1255ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0 && ze-Ze < 0) { 12568865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 125747c6ae99SBarry Smith } 125847c6ae99SBarry Smith if (n19 >= 0) { /* directly below */ 125947c6ae99SBarry Smith x_t = x; 126047c6ae99SBarry Smith y_t = ly[(n19 % (m*n))/m]; 126147c6ae99SBarry Smith /* z_t = lz[n19 / (m*n)]; */ 126247c6ae99SBarry Smith s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 12638865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1264ce00eea3SSatish Balay } else if (Ys-ys < 0 && ze-Ze < 0) { 12658865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 126647c6ae99SBarry Smith } 126747c6ae99SBarry Smith if (n20 >= 0) { /* right below */ 1268ce00eea3SSatish Balay x_t = lx[n20 % m]; 126947c6ae99SBarry Smith y_t = ly[(n20 % (m*n))/m]; 127047c6ae99SBarry Smith /* z_t = lz[n20 / (m*n)]; */ 127147c6ae99SBarry Smith s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 12728865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1273ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0 && ze-Ze < 0) { 12748865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 127547c6ae99SBarry Smith } 127647c6ae99SBarry Smith } 127747c6ae99SBarry Smith 127847c6ae99SBarry Smith for (i=0; i<y; i++) { 127947c6ae99SBarry Smith if (n21 >= 0) { /* directly left */ 1280ce00eea3SSatish Balay x_t = lx[n21 % m]; 128147c6ae99SBarry Smith y_t = y; 128247c6ae99SBarry Smith /* z_t = lz[n21 / (m*n)]; */ 128347c6ae99SBarry Smith s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t; 12848865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1285ce00eea3SSatish Balay } else if (Xs-xs < 0 && ze-Ze < 0) { 12868865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 128747c6ae99SBarry Smith } 128847c6ae99SBarry Smith 128947c6ae99SBarry Smith if (n22 >= 0) { /* middle */ 129047c6ae99SBarry Smith x_t = x; 129147c6ae99SBarry Smith y_t = y; 129247c6ae99SBarry Smith /* z_t = lz[n22 / (m*n)]; */ 129347c6ae99SBarry Smith s_t = bases[n22] + i*x_t + k*x_t*y_t; 12948865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1295ce00eea3SSatish Balay } else if (ze-Ze < 0) { 1296bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_MIRROR) { 12978865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 1298c2859e5eSBarry Smith } else { 12998865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 130047c6ae99SBarry Smith } 1301c2859e5eSBarry Smith } 130247c6ae99SBarry Smith 130347c6ae99SBarry Smith if (n23 >= 0) { /* directly right */ 1304ce00eea3SSatish Balay x_t = lx[n23 % m]; 130547c6ae99SBarry Smith y_t = y; 130647c6ae99SBarry Smith /* z_t = lz[n23 / (m*n)]; */ 130747c6ae99SBarry Smith s_t = bases[n23] + i*x_t + k*x_t*y_t; 13088865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1309ce00eea3SSatish Balay } else if (xe-Xe < 0 && ze-Ze < 0) { 13108865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 131147c6ae99SBarry Smith } 131247c6ae99SBarry Smith } 131347c6ae99SBarry Smith 131447c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 131547c6ae99SBarry Smith if (n24 >= 0) { /* left above */ 1316ce00eea3SSatish Balay x_t = lx[n24 % m]; 131747c6ae99SBarry Smith y_t = ly[(n24 % (m*n))/m]; 131847c6ae99SBarry Smith /* z_t = lz[n24 / (m*n)]; */ 131947c6ae99SBarry Smith s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t; 13208865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1321ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0 && ze-Ze < 0) { 13228865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 132347c6ae99SBarry Smith } 132447c6ae99SBarry Smith if (n25 >= 0) { /* directly above */ 132547c6ae99SBarry Smith x_t = x; 132647c6ae99SBarry Smith y_t = ly[(n25 % (m*n))/m]; 132747c6ae99SBarry Smith /* z_t = lz[n25 / (m*n)]; */ 132847c6ae99SBarry Smith s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t; 13298865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1330ce00eea3SSatish Balay } else if (ye-Ye < 0 && ze-Ze < 0) { 13318865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 133247c6ae99SBarry Smith } 133347c6ae99SBarry Smith if (n26 >= 0) { /* right above */ 1334ce00eea3SSatish Balay x_t = lx[n26 % m]; 133547c6ae99SBarry Smith y_t = ly[(n26 % (m*n))/m]; 133647c6ae99SBarry Smith /* z_t = lz[n26 / (m*n)]; */ 133747c6ae99SBarry Smith s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t; 13388865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1339ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0 && ze-Ze < 0) { 13408865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 134147c6ae99SBarry Smith } 134247c6ae99SBarry Smith } 134347c6ae99SBarry Smith } 134447c6ae99SBarry Smith } 134547c6ae99SBarry Smith /* 134647c6ae99SBarry Smith Set the local to global ordering in the global vector, this allows use 134747c6ae99SBarry Smith of VecSetValuesLocal(). 134847c6ae99SBarry Smith */ 134945b6f7e9SBarry Smith ierr = ISLocalToGlobalMappingCreate(comm,dof,nn,idx,PETSC_OWN_POINTER,&da->ltogmap);CHKERRQ(ierr); 13503bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap);CHKERRQ(ierr); 135147c6ae99SBarry Smith 1352ce00eea3SSatish Balay ierr = PetscFree2(bases,ldims);CHKERRQ(ierr); 1353ce00eea3SSatish Balay dd->m = m; dd->n = n; dd->p = p; 1354ce00eea3SSatish Balay /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */ 1355ce00eea3SSatish Balay dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = zs; dd->ze = ze; 1356ce00eea3SSatish Balay dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = Zs; dd->Ze = Ze; 1357ce00eea3SSatish Balay 1358fcfd50ebSBarry Smith ierr = VecDestroy(&local);CHKERRQ(ierr); 1359fcfd50ebSBarry Smith ierr = VecDestroy(&global);CHKERRQ(ierr); 1360ce00eea3SSatish Balay 1361ce00eea3SSatish Balay dd->gtol = gtol; 1362ce00eea3SSatish Balay dd->base = base; 1363ce00eea3SSatish Balay da->ops->view = DMView_DA_3d; 13640298fd71SBarry Smith dd->ltol = NULL; 13650298fd71SBarry Smith dd->ao = NULL; 136647c6ae99SBarry Smith PetscFunctionReturn(0); 136747c6ae99SBarry Smith } 1368564755cdSBarry Smith 136947c6ae99SBarry Smith 137047c6ae99SBarry Smith #undef __FUNCT__ 1371aa219208SBarry Smith #define __FUNCT__ "DMDACreate3d" 137247c6ae99SBarry Smith /*@C 1373aa219208SBarry Smith DMDACreate3d - Creates an object that will manage the communication of three-dimensional 137447c6ae99SBarry Smith regular array data that is distributed across some processors. 137547c6ae99SBarry Smith 137647c6ae99SBarry Smith Collective on MPI_Comm 137747c6ae99SBarry Smith 137847c6ae99SBarry Smith Input Parameters: 137947c6ae99SBarry Smith + comm - MPI communicator 13801321219cSEthan Coon . bx,by,bz - type of ghost nodes the array have. 1381bff4a2f0SMatthew G. Knepley Use one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC. 1382aa219208SBarry Smith . stencil_type - Type of stencil (DMDA_STENCIL_STAR or DMDA_STENCIL_BOX) 138347c6ae99SBarry Smith . M,N,P - global dimension in each direction of the array (use -M, -N, and or -P to indicate that it may be set to a different value 138447c6ae99SBarry Smith from the command line with -da_grid_x <M> -da_grid_y <N> -da_grid_z <P>) 138547c6ae99SBarry Smith . m,n,p - corresponding number of processors in each dimension 138647c6ae99SBarry Smith (or PETSC_DECIDE to have calculated) 138747c6ae99SBarry Smith . dof - number of degrees of freedom per node 138810d7c030SMatthew G Knepley . s - stencil width 138910d7c030SMatthew G Knepley - lx, ly, lz - arrays containing the number of nodes in each cell along 13900298fd71SBarry Smith the x, y, and z coordinates, or NULL. If non-null, these 139147c6ae99SBarry Smith must be of length as m,n,p and the corresponding 139247c6ae99SBarry Smith m,n, or p cannot be PETSC_DECIDE. Sum of the lx[] entries must be M, sum of 139347c6ae99SBarry Smith the ly[] must N, sum of the lz[] must be P 139447c6ae99SBarry Smith 139547c6ae99SBarry Smith Output Parameter: 139647c6ae99SBarry Smith . da - the resulting distributed array object 139747c6ae99SBarry Smith 139847c6ae99SBarry Smith Options Database Key: 1399706a11cbSBarry Smith + -dm_view - Calls DMView() at the conclusion of DMDACreate3d() 140047c6ae99SBarry Smith . -da_grid_x <nx> - number of grid points in x direction, if M < 0 140147c6ae99SBarry Smith . -da_grid_y <ny> - number of grid points in y direction, if N < 0 140247c6ae99SBarry Smith . -da_grid_z <nz> - number of grid points in z direction, if P < 0 140347c6ae99SBarry Smith . -da_processors_x <MX> - number of processors in x direction 140447c6ae99SBarry Smith . -da_processors_y <MY> - number of processors in y direction 140547c6ae99SBarry Smith . -da_processors_z <MZ> - number of processors in z direction 1406e0f5d30fSBarry Smith . -da_refine_x <rx> - refinement ratio in x direction 1407e0f5d30fSBarry Smith . -da_refine_y <ry> - refinement ratio in y direction 1408e0f5d30fSBarry Smith . -da_refine_z <rz>- refinement ratio in z directio 1409e0f5d30fSBarry Smith - -da_refine <n> - refine the DMDA n times before creating it, , if M, N, or P < 0 141047c6ae99SBarry Smith 141147c6ae99SBarry Smith Level: beginner 141247c6ae99SBarry Smith 141347c6ae99SBarry Smith Notes: 1414aa219208SBarry Smith The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the 1415aa219208SBarry Smith standard 7-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes 141647c6ae99SBarry Smith the standard 27-pt stencil. 141747c6ae99SBarry Smith 1418aa219208SBarry Smith The array data itself is NOT stored in the DMDA, it is stored in Vec objects; 1419564755cdSBarry Smith The appropriate vector objects can be obtained with calls to DMCreateGlobalVector() 1420564755cdSBarry Smith and DMCreateLocalVector() and calls to VecDuplicate() if more are needed. 142147c6ae99SBarry Smith 142247c6ae99SBarry Smith .keywords: distributed array, create, three-dimensional 142347c6ae99SBarry Smith 1424aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate2d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(), 142599f0b487SRichard Tran Mills DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDASetRefinementFactor(), 1426d461ba97SBarry Smith DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges() 142747c6ae99SBarry Smith 142847c6ae99SBarry Smith @*/ 1429bff4a2f0SMatthew G. Knepley PetscErrorCode DMDACreate3d(MPI_Comm comm,DMBoundaryType bx,DMBoundaryType by,DMBoundaryType bz,DMDAStencilType stencil_type,PetscInt M, 14309a42bb27SBarry Smith PetscInt N,PetscInt P,PetscInt m,PetscInt n,PetscInt p,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],const PetscInt lz[],DM *da) 143147c6ae99SBarry Smith { 143247c6ae99SBarry Smith PetscErrorCode ierr; 143347c6ae99SBarry Smith 143447c6ae99SBarry Smith PetscFunctionBegin; 1435aa219208SBarry Smith ierr = DMDACreate(comm, da);CHKERRQ(ierr); 1436c73cfb54SMatthew G. Knepley ierr = DMSetDimension(*da, 3);CHKERRQ(ierr); 1437aa219208SBarry Smith ierr = DMDASetSizes(*da, M, N, P);CHKERRQ(ierr); 1438aa219208SBarry Smith ierr = DMDASetNumProcs(*da, m, n, p);CHKERRQ(ierr); 14391321219cSEthan Coon ierr = DMDASetBoundaryType(*da, bx, by, bz);CHKERRQ(ierr); 1440aa219208SBarry Smith ierr = DMDASetDof(*da, dof);CHKERRQ(ierr); 1441aa219208SBarry Smith ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr); 1442aa219208SBarry Smith ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr); 1443aa219208SBarry Smith ierr = DMDASetOwnershipRanges(*da, lx, ly, lz);CHKERRQ(ierr); 144447c6ae99SBarry Smith /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */ 14459a42bb27SBarry Smith ierr = DMSetFromOptions(*da);CHKERRQ(ierr); 14469a42bb27SBarry Smith ierr = DMSetUp(*da);CHKERRQ(ierr); 144747c6ae99SBarry Smith PetscFunctionReturn(0); 144847c6ae99SBarry Smith } 1449