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> 10e0877f53SBarry Smith static PetscErrorCode DMView_DA_3d(DM da,PetscViewer viewer) 1147c6ae99SBarry Smith { 1247c6ae99SBarry Smith PetscErrorCode ierr; 1347c6ae99SBarry Smith PetscMPIInt rank; 14c9493c35SStefano Zampini PetscBool iascii,isdraw,isglvis,isbinary; 1547c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 169a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 179a42bb27SBarry Smith PetscBool ismatlab; 189a42bb27SBarry Smith #endif 1947c6ae99SBarry Smith 2047c6ae99SBarry Smith PetscFunctionBegin; 21ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank);CHKERRQ(ierr); 2247c6ae99SBarry Smith 23251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 24251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 25c9493c35SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis);CHKERRQ(ierr); 26251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 279a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 28251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab);CHKERRQ(ierr); 299a42bb27SBarry Smith #endif 3047c6ae99SBarry Smith if (iascii) { 3147c6ae99SBarry Smith PetscViewerFormat format; 3247c6ae99SBarry Smith 331575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 3447c6ae99SBarry Smith ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 358135c375SStefano Zampini if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL && format != PETSC_VIEWER_ASCII_GLVIS) { 36aa219208SBarry Smith DMDALocalInfo info; 37aa219208SBarry Smith ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr); 3847c6ae99SBarry 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); 3947c6ae99SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"X range of indices: %D %D, Y range of indices: %D %D, Z range of indices: %D %D\n", 4047c6ae99SBarry Smith info.xs,info.xs+info.xm,info.ys,info.ys+info.ym,info.zs,info.zs+info.zm);CHKERRQ(ierr); 4147c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX) 426636e97aSMatthew G Knepley if (da->coordinates) { 4347c6ae99SBarry Smith PetscInt last; 4447c6ae99SBarry Smith const PetscReal *coors; 456636e97aSMatthew G Knepley ierr = VecGetArrayRead(da->coordinates,&coors);CHKERRQ(ierr); 466636e97aSMatthew G Knepley ierr = VecGetLocalSize(da->coordinates,&last);CHKERRQ(ierr); 4747c6ae99SBarry Smith last = last - 3; 4857622a8eSBarry 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); 496636e97aSMatthew G Knepley ierr = VecRestoreArrayRead(da->coordinates,&coors);CHKERRQ(ierr); 5047c6ae99SBarry Smith } 5147c6ae99SBarry Smith #endif 5247c6ae99SBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 531575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 548135c375SStefano Zampini } else if (format == PETSC_VIEWER_ASCII_GLVIS) { 558135c375SStefano Zampini ierr = DMView_DA_GLVis(da,viewer);CHKERRQ(ierr); 56616157d6SJed Brown } else { 57616157d6SJed Brown ierr = DMView_DA_VTK(da,viewer);CHKERRQ(ierr); 5847c6ae99SBarry Smith } 5947c6ae99SBarry Smith } else if (isdraw) { 6047c6ae99SBarry Smith PetscDraw draw; 6147c6ae99SBarry Smith PetscReal ymin = -1.0,ymax = (PetscReal)dd->N; 6247c6ae99SBarry Smith PetscReal xmin = -1.0,xmax = (PetscReal)((dd->M+2)*dd->P),x,y,ycoord,xcoord; 638ea3bf28SBarry Smith PetscInt k,plane,base; 648ea3bf28SBarry Smith const PetscInt *idx; 6547c6ae99SBarry Smith char node[10]; 6647c6ae99SBarry Smith PetscBool isnull; 6747c6ae99SBarry Smith 6847c6ae99SBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 695b399a63SLisandro Dalcin ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 705b399a63SLisandro Dalcin if (isnull) PetscFunctionReturn(0); 7147c6ae99SBarry Smith 725b399a63SLisandro Dalcin ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr); 735b399a63SLisandro Dalcin ierr = PetscDrawClear(draw);CHKERRQ(ierr); 745b399a63SLisandro Dalcin ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr); 755b399a63SLisandro Dalcin 765b399a63SLisandro Dalcin ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr); 7747c6ae99SBarry Smith /* first processor draw all node lines */ 7847c6ae99SBarry Smith if (!rank) { 7947c6ae99SBarry Smith for (k=0; k<dd->P; k++) { 8047c6ae99SBarry Smith ymin = 0.0; ymax = (PetscReal)(dd->N - 1); 8147c6ae99SBarry Smith for (xmin=(PetscReal)(k*(dd->M+1)); xmin<(PetscReal)(dd->M+(k*(dd->M+1))); xmin++) { 8247c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr); 8347c6ae99SBarry Smith } 8447c6ae99SBarry Smith xmin = (PetscReal)(k*(dd->M+1)); xmax = xmin + (PetscReal)(dd->M - 1); 8547c6ae99SBarry Smith for (ymin=0; ymin<(PetscReal)dd->N; ymin++) { 8647c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr); 8747c6ae99SBarry Smith } 8847c6ae99SBarry Smith } 8947c6ae99SBarry Smith } 905b399a63SLisandro Dalcin ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr); 915b399a63SLisandro Dalcin ierr = PetscDrawFlush(draw);CHKERRQ(ierr); 9247c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 9347c6ae99SBarry Smith 945b399a63SLisandro Dalcin ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr); 955b399a63SLisandro Dalcin /*Go through and draw for each plane*/ 965b399a63SLisandro Dalcin for (k=0; k<dd->P; k++) { 9747c6ae99SBarry Smith if ((k >= dd->zs) && (k < dd->ze)) { 9847c6ae99SBarry Smith /* draw my box */ 9947c6ae99SBarry Smith ymin = dd->ys; 10047c6ae99SBarry Smith ymax = dd->ye - 1; 10147c6ae99SBarry Smith xmin = dd->xs/dd->w + (dd->M+1)*k; 10247c6ae99SBarry Smith xmax =(dd->xe-1)/dd->w + (dd->M+1)*k; 10347c6ae99SBarry Smith 10447c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr); 10547c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 10647c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 10747c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 10847c6ae99SBarry Smith 10947c6ae99SBarry Smith xmin = dd->xs/dd->w; 11047c6ae99SBarry Smith xmax =(dd->xe-1)/dd->w; 11147c6ae99SBarry Smith 112832b7cebSLisandro Dalcin /* identify which processor owns the box */ 113832b7cebSLisandro Dalcin ierr = PetscSNPrintf(node,sizeof(node),"%d",(int)rank);CHKERRQ(ierr); 114832b7cebSLisandro Dalcin ierr = PetscDrawString(draw,xmin+(dd->M+1)*k+.2,ymin+.3,PETSC_DRAW_RED,node);CHKERRQ(ierr); 11547c6ae99SBarry Smith /* put in numbers*/ 11647c6ae99SBarry Smith base = (dd->base+(dd->xe-dd->xs)*(dd->ye-dd->ys)*(k-dd->zs))/dd->w; 11747c6ae99SBarry Smith for (y=ymin; y<=ymax; y++) { 11847c6ae99SBarry Smith for (x=xmin+(dd->M+1)*k; x<=xmax+(dd->M+1)*k; x++) { 119832b7cebSLisandro Dalcin ierr = PetscSNPrintf(node,sizeof(node),"%d",(int)base++);CHKERRQ(ierr); 12047c6ae99SBarry Smith ierr = PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);CHKERRQ(ierr); 12147c6ae99SBarry Smith } 12247c6ae99SBarry Smith } 12347c6ae99SBarry Smith 12447c6ae99SBarry Smith } 12547c6ae99SBarry Smith } 1265b399a63SLisandro Dalcin ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr); 1275b399a63SLisandro Dalcin ierr = PetscDrawFlush(draw);CHKERRQ(ierr); 12847c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 12947c6ae99SBarry Smith 1305b399a63SLisandro Dalcin ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr); 13147c6ae99SBarry Smith for (k=0-dd->s; k<dd->P+dd->s; k++) { 13247c6ae99SBarry Smith /* Go through and draw for each plane */ 13347c6ae99SBarry Smith if ((k >= dd->Zs) && (k < dd->Ze)) { 13447c6ae99SBarry Smith /* overlay ghost numbers, useful for error checking */ 135565245c5SBarry Smith base = (dd->Xe-dd->Xs)*(dd->Ye-dd->Ys)*(k-dd->Zs)/dd->w; 13645b6f7e9SBarry Smith ierr = ISLocalToGlobalMappingGetBlockIndices(da->ltogmap,&idx);CHKERRQ(ierr); 13747c6ae99SBarry Smith plane=k; 138565245c5SBarry Smith /* Keep z wrap around points on the drawing */ 1398865f1eaSKarl Rupp if (k<0) plane=dd->P+k; 1408865f1eaSKarl Rupp if (k>=dd->P) plane=k-dd->P; 14147c6ae99SBarry Smith ymin = dd->Ys; ymax = dd->Ye; 14247c6ae99SBarry Smith xmin = (dd->M+1)*plane*dd->w; 14347c6ae99SBarry Smith xmax = (dd->M+1)*plane*dd->w+dd->M*dd->w; 14447c6ae99SBarry Smith for (y=ymin; y<ymax; y++) { 14547c6ae99SBarry Smith for (x=xmin+dd->Xs; x<xmin+dd->Xe; x+=dd->w) { 1468ea3bf28SBarry Smith sprintf(node,"%d",(int)(idx[base])); 14747c6ae99SBarry Smith ycoord = y; 14847c6ae99SBarry Smith /*Keep y wrap around points on drawing */ 1498865f1eaSKarl Rupp if (y<0) ycoord = dd->N+y; 1508865f1eaSKarl Rupp if (y>=dd->N) ycoord = y-dd->N; 15147c6ae99SBarry Smith xcoord = x; /* Keep x wrap points on drawing */ 1528865f1eaSKarl Rupp if (x<xmin) xcoord = xmax - (xmin-x); 1538865f1eaSKarl Rupp if (x>=xmax) xcoord = xmin + (x-xmax); 15447c6ae99SBarry Smith ierr = PetscDrawString(draw,xcoord/dd->w,ycoord,PETSC_DRAW_BLUE,node);CHKERRQ(ierr); 155565245c5SBarry Smith base++; 15647c6ae99SBarry Smith } 15747c6ae99SBarry Smith } 15845b6f7e9SBarry Smith ierr = ISLocalToGlobalMappingRestoreBlockIndices(da->ltogmap,&idx);CHKERRQ(ierr); 15947c6ae99SBarry Smith } 16047c6ae99SBarry Smith } 1615b399a63SLisandro Dalcin ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr); 1625b399a63SLisandro Dalcin ierr = PetscDrawFlush(draw);CHKERRQ(ierr); 16347c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 164832b7cebSLisandro Dalcin ierr = PetscDrawSave(draw);CHKERRQ(ierr); 165c9493c35SStefano Zampini } else if (isglvis) { 166c9493c35SStefano Zampini ierr = DMView_DA_GLVis(da,viewer);CHKERRQ(ierr); 1679a42bb27SBarry Smith } else if (isbinary) { 1689a42bb27SBarry Smith ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr); 1699a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 1709a42bb27SBarry Smith } else if (ismatlab) { 1719a42bb27SBarry Smith ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr); 1729a42bb27SBarry Smith #endif 17311aeaf0aSBarry Smith } 17447c6ae99SBarry Smith PetscFunctionReturn(0); 17547c6ae99SBarry Smith } 17647c6ae99SBarry Smith 1777087cfbeSBarry Smith PetscErrorCode DMSetUp_DA_3D(DM da) 17847c6ae99SBarry Smith { 17947c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 18047c6ae99SBarry Smith const PetscInt M = dd->M; 18147c6ae99SBarry Smith const PetscInt N = dd->N; 18247c6ae99SBarry Smith const PetscInt P = dd->P; 18347c6ae99SBarry Smith PetscInt m = dd->m; 18447c6ae99SBarry Smith PetscInt n = dd->n; 18547c6ae99SBarry Smith PetscInt p = dd->p; 18647c6ae99SBarry Smith const PetscInt dof = dd->w; 18747c6ae99SBarry Smith const PetscInt s = dd->s; 188bff4a2f0SMatthew G. Knepley DMBoundaryType bx = dd->bx; 189bff4a2f0SMatthew G. Knepley DMBoundaryType by = dd->by; 190bff4a2f0SMatthew G. Knepley DMBoundaryType bz = dd->bz; 19119fd82e9SBarry Smith DMDAStencilType stencil_type = dd->stencil_type; 19247c6ae99SBarry Smith PetscInt *lx = dd->lx; 19347c6ae99SBarry Smith PetscInt *ly = dd->ly; 19447c6ae99SBarry Smith PetscInt *lz = dd->lz; 19547c6ae99SBarry Smith MPI_Comm comm; 19647c6ae99SBarry Smith PetscMPIInt rank,size; 197ce00eea3SSatish Balay PetscInt xs = 0,xe,ys = 0,ye,zs = 0,ze,x = 0,y = 0,z = 0; 198bd1fc5aeSBarry Smith PetscInt Xs,Xe,Ys,Ye,Zs,Ze,IXs,IXe,IYs,IYe,IZs,IZe,pm; 1998ea3bf28SBarry Smith PetscInt left,right,up,down,bottom,top,i,j,k,*idx,nn; 20047c6ae99SBarry Smith PetscInt n0,n1,n2,n3,n4,n5,n6,n7,n8,n9,n10,n11,n12,n14; 20147c6ae99SBarry Smith PetscInt n15,n16,n17,n18,n19,n20,n21,n22,n23,n24,n25,n26; 202db87c5ecSEthan Coon PetscInt *bases,*ldims,base,x_t,y_t,z_t,s_t,count,s_x,s_y,s_z; 20347c6ae99SBarry Smith PetscInt sn0 = 0,sn1 = 0,sn2 = 0,sn3 = 0,sn5 = 0,sn6 = 0,sn7 = 0; 20447c6ae99SBarry Smith PetscInt sn8 = 0,sn9 = 0,sn11 = 0,sn15 = 0,sn24 = 0,sn25 = 0,sn26 = 0; 20547c6ae99SBarry Smith PetscInt sn17 = 0,sn18 = 0,sn19 = 0,sn20 = 0,sn21 = 0,sn23 = 0; 20647c6ae99SBarry Smith Vec local,global; 207bd1fc5aeSBarry Smith VecScatter gtol; 20845b6f7e9SBarry Smith IS to,from; 2096f951b95Secoon PetscBool twod; 21047c6ae99SBarry Smith PetscErrorCode ierr; 21147c6ae99SBarry Smith 2126f951b95Secoon 21347c6ae99SBarry Smith PetscFunctionBegin; 214bff4a2f0SMatthew 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"); 21540088605SBarry Smith if ((bx == DM_BOUNDARY_MIRROR) || (by == DM_BOUNDARY_MIRROR) || (bz == DM_BOUNDARY_MIRROR)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Mirror boundary not supported yet in 3d"); 21647c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject) da, &comm);CHKERRQ(ierr); 2173855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES) 218bafee8b4SSatish Balay if (((PetscInt64) M)*((PetscInt64) N)*((PetscInt64) P)*((PetscInt64) dof) > (PetscInt64) 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); 2193855c12bSBarry Smith #endif 2203855c12bSBarry Smith 22147c6ae99SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 22247c6ae99SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 22347c6ae99SBarry Smith 22447c6ae99SBarry Smith if (m != PETSC_DECIDE) { 22547c6ae99SBarry Smith if (m < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m); 22647c6ae99SBarry Smith else if (m > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size); 22747c6ae99SBarry Smith } 22847c6ae99SBarry Smith if (n != PETSC_DECIDE) { 22947c6ae99SBarry Smith if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n); 23047c6ae99SBarry Smith else if (n > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size); 23147c6ae99SBarry Smith } 23247c6ae99SBarry Smith if (p != PETSC_DECIDE) { 23347c6ae99SBarry Smith if (p < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Z direction: %D",p); 23447c6ae99SBarry Smith else if (p > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Z direction: %D %d",p,size); 23547c6ae99SBarry Smith } 2360716a85fSBarry 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); 23747c6ae99SBarry Smith 23847c6ae99SBarry Smith /* Partition the array among the processors */ 23947c6ae99SBarry Smith if (m == PETSC_DECIDE && n != PETSC_DECIDE && p != PETSC_DECIDE) { 24047c6ae99SBarry Smith m = size/(n*p); 24147c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) { 24247c6ae99SBarry Smith n = size/(m*p); 24347c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) { 24447c6ae99SBarry Smith p = size/(m*n); 24547c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) { 24647c6ae99SBarry Smith /* try for squarish distribution */ 247369cc0aeSBarry Smith m = (int)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)N*p))); 24847c6ae99SBarry Smith if (!m) m = 1; 24947c6ae99SBarry Smith while (m > 0) { 25047c6ae99SBarry Smith n = size/(m*p); 25147c6ae99SBarry Smith if (m*n*p == size) break; 25247c6ae99SBarry Smith m--; 25347c6ae99SBarry Smith } 25447c6ae99SBarry Smith if (!m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad p value: p = %D",p); 25547c6ae99SBarry Smith if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;} 25647c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) { 25747c6ae99SBarry Smith /* try for squarish distribution */ 258369cc0aeSBarry Smith m = (int)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n))); 25947c6ae99SBarry Smith if (!m) m = 1; 26047c6ae99SBarry Smith while (m > 0) { 26147c6ae99SBarry Smith p = size/(m*n); 26247c6ae99SBarry Smith if (m*n*p == size) break; 26347c6ae99SBarry Smith m--; 26447c6ae99SBarry Smith } 26547c6ae99SBarry Smith if (!m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad n value: n = %D",n); 26647c6ae99SBarry Smith if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;} 26747c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) { 26847c6ae99SBarry Smith /* try for squarish distribution */ 269369cc0aeSBarry Smith n = (int)(0.5 + PetscSqrtReal(((PetscReal)N)*((PetscReal)size)/((PetscReal)P*m))); 27047c6ae99SBarry Smith if (!n) n = 1; 27147c6ae99SBarry Smith while (n > 0) { 27247c6ae99SBarry Smith p = size/(m*n); 27347c6ae99SBarry Smith if (m*n*p == size) break; 27447c6ae99SBarry Smith n--; 27547c6ae99SBarry Smith } 27647c6ae99SBarry Smith if (!n) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad m value: m = %D",n); 27747c6ae99SBarry Smith if (N > P && n < p) {PetscInt _n = n; n = p; p = _n;} 27847c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) { 27947c6ae99SBarry Smith /* try for squarish distribution */ 280369cc0aeSBarry Smith n = (PetscInt)(0.5 + PetscPowReal(((PetscReal)N*N)*((PetscReal)size)/((PetscReal)P*M),(PetscReal)(1./3.))); 28147c6ae99SBarry Smith if (!n) n = 1; 28247c6ae99SBarry Smith while (n > 0) { 28347c6ae99SBarry Smith pm = size/n; 28447c6ae99SBarry Smith if (n*pm == size) break; 28547c6ae99SBarry Smith n--; 28647c6ae99SBarry Smith } 28747c6ae99SBarry Smith if (!n) n = 1; 288369cc0aeSBarry Smith m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n))); 28947c6ae99SBarry Smith if (!m) m = 1; 29047c6ae99SBarry Smith while (m > 0) { 29147c6ae99SBarry Smith p = size/(m*n); 29247c6ae99SBarry Smith if (m*n*p == size) break; 29347c6ae99SBarry Smith m--; 29447c6ae99SBarry Smith } 29547c6ae99SBarry Smith if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;} 296ce94432eSBarry Smith } else if (m*n*p != size) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition"); 29747c6ae99SBarry Smith 298ce94432eSBarry Smith if (m*n*p != size) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_PLIB,"Could not find good partition"); 299ce94432eSBarry Smith if (M < m) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m); 300ce94432eSBarry Smith if (N < n) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n); 301ce94432eSBarry Smith if (P < p) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Partition in z direction is too fine! %D %D",P,p); 30247c6ae99SBarry Smith 30347c6ae99SBarry Smith /* 30447c6ae99SBarry Smith Determine locally owned region 30547c6ae99SBarry Smith [x, y, or z]s is the first local node number, [x, y, z] is the number of local nodes 30647c6ae99SBarry Smith */ 30747c6ae99SBarry Smith 30847c6ae99SBarry Smith if (!lx) { 309785e854fSJed Brown ierr = PetscMalloc1(m, &dd->lx);CHKERRQ(ierr); 31047c6ae99SBarry Smith lx = dd->lx; 3118865f1eaSKarl Rupp for (i=0; i<m; i++) lx[i] = M/m + ((M % m) > (i % m)); 31247c6ae99SBarry Smith } 31347c6ae99SBarry Smith x = lx[rank % m]; 31447c6ae99SBarry Smith xs = 0; 3158865f1eaSKarl Rupp for (i=0; i<(rank%m); i++) xs += lx[i]; 316bff4a2f0SMatthew 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); 31747c6ae99SBarry Smith 31847c6ae99SBarry Smith if (!ly) { 319785e854fSJed Brown ierr = PetscMalloc1(n, &dd->ly);CHKERRQ(ierr); 32047c6ae99SBarry Smith ly = dd->ly; 3218865f1eaSKarl Rupp for (i=0; i<n; i++) ly[i] = N/n + ((N % n) > (i % n)); 32247c6ae99SBarry Smith } 32347c6ae99SBarry Smith y = ly[(rank % (m*n))/m]; 324bff4a2f0SMatthew 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); 32530729d88SBarry Smith 32647c6ae99SBarry Smith ys = 0; 3278865f1eaSKarl Rupp for (i=0; i<(rank % (m*n))/m; i++) ys += ly[i]; 32847c6ae99SBarry Smith 32947c6ae99SBarry Smith if (!lz) { 330785e854fSJed Brown ierr = PetscMalloc1(p, &dd->lz);CHKERRQ(ierr); 33147c6ae99SBarry Smith lz = dd->lz; 3328865f1eaSKarl Rupp for (i=0; i<p; i++) lz[i] = P/p + ((P % p) > (i % p)); 33347c6ae99SBarry Smith } 33447c6ae99SBarry Smith z = lz[rank/(m*n)]; 335bcea557cSEthan Coon 336fdc81ce8SEthan Coon /* note this is different than x- and y-, as we will handle as an important special 337bff4a2f0SMatthew G. Knepley case when p=P=1 and DM_BOUNDARY_PERIODIC and s > z. This is to deal with 2D problems 338fdc81ce8SEthan Coon in a 3D code. Additional code for this case is noted with "2d case" comments */ 3396f951b95Secoon twod = PETSC_FALSE; 3408865f1eaSKarl Rupp if (P == 1) twod = PETSC_TRUE; 341bff4a2f0SMatthew 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); 34247c6ae99SBarry Smith zs = 0; 3438865f1eaSKarl Rupp for (i=0; i<(rank/(m*n)); i++) zs += lz[i]; 34447c6ae99SBarry Smith ye = ys + y; 34547c6ae99SBarry Smith xe = xs + x; 34647c6ae99SBarry Smith ze = zs + z; 34747c6ae99SBarry Smith 34888661749SPeter Brune /* determine ghost region (Xs) and region scattered into (IXs) */ 349d9c9ebe5SPeter Brune if (xs-s > 0) { 350d9c9ebe5SPeter Brune Xs = xs - s; IXs = xs - s; 35188661749SPeter Brune } else { 3528865f1eaSKarl Rupp if (bx) Xs = xs - s; 3538865f1eaSKarl Rupp else Xs = 0; 35488661749SPeter Brune IXs = 0; 35588661749SPeter Brune } 356d9c9ebe5SPeter Brune if (xe+s <= M) { 357d9c9ebe5SPeter Brune Xe = xe + s; IXe = xe + s; 35888661749SPeter Brune } else { 35988661749SPeter Brune if (bx) { 360d9c9ebe5SPeter Brune Xs = xs - s; Xe = xe + s; 3618865f1eaSKarl Rupp } else Xe = M; 36288661749SPeter Brune IXe = M; 36388661749SPeter Brune } 36447c6ae99SBarry Smith 365bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) { 366d9c9ebe5SPeter Brune IXs = xs - s; 367d9c9ebe5SPeter Brune IXe = xe + s; 368d9c9ebe5SPeter Brune Xs = xs - s; 369d9c9ebe5SPeter Brune Xe = xe + s; 37088661749SPeter Brune } 37188661749SPeter Brune 372d9c9ebe5SPeter Brune if (ys-s > 0) { 373d9c9ebe5SPeter Brune Ys = ys - s; IYs = ys - s; 37488661749SPeter Brune } else { 3758865f1eaSKarl Rupp if (by) Ys = ys - s; 3768865f1eaSKarl Rupp else Ys = 0; 37788661749SPeter Brune IYs = 0; 37888661749SPeter Brune } 379d9c9ebe5SPeter Brune if (ye+s <= N) { 380d9c9ebe5SPeter Brune Ye = ye + s; IYe = ye + s; 38188661749SPeter Brune } else { 3828865f1eaSKarl Rupp if (by) Ye = ye + s; 3838865f1eaSKarl Rupp else Ye = N; 38488661749SPeter Brune IYe = N; 38588661749SPeter Brune } 38688661749SPeter Brune 387bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC || by == DM_BOUNDARY_MIRROR) { 388d9c9ebe5SPeter Brune IYs = ys - s; 389d9c9ebe5SPeter Brune IYe = ye + s; 390d9c9ebe5SPeter Brune Ys = ys - s; 391d9c9ebe5SPeter Brune Ye = ye + s; 39288661749SPeter Brune } 39388661749SPeter Brune 394d9c9ebe5SPeter Brune if (zs-s > 0) { 395d9c9ebe5SPeter Brune Zs = zs - s; IZs = zs - s; 39688661749SPeter Brune } else { 3978865f1eaSKarl Rupp if (bz) Zs = zs - s; 3988865f1eaSKarl Rupp else Zs = 0; 39988661749SPeter Brune IZs = 0; 40088661749SPeter Brune } 401d9c9ebe5SPeter Brune if (ze+s <= P) { 402d9c9ebe5SPeter Brune Ze = ze + s; IZe = ze + s; 40388661749SPeter Brune } else { 4048865f1eaSKarl Rupp if (bz) Ze = ze + s; 4058865f1eaSKarl Rupp else Ze = P; 40688661749SPeter Brune IZe = P; 40788661749SPeter Brune } 40888661749SPeter Brune 409bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_PERIODIC || bz == DM_BOUNDARY_MIRROR) { 410d9c9ebe5SPeter Brune IZs = zs - s; 411d9c9ebe5SPeter Brune IZe = ze + s; 412d9c9ebe5SPeter Brune Zs = zs - s; 413d9c9ebe5SPeter Brune Ze = ze + s; 41488661749SPeter Brune } 41547c6ae99SBarry Smith 41647c6ae99SBarry Smith /* Resize all X parameters to reflect w */ 417d9c9ebe5SPeter Brune s_x = s; 418d9c9ebe5SPeter Brune s_y = s; 419d9c9ebe5SPeter Brune s_z = s; 42047c6ae99SBarry Smith 42147c6ae99SBarry Smith /* determine starting point of each processor */ 42247c6ae99SBarry Smith nn = x*y*z; 423dcca6d9dSJed Brown ierr = PetscMalloc2(size+1,&bases,size,&ldims);CHKERRQ(ierr); 42447c6ae99SBarry Smith ierr = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr); 42547c6ae99SBarry Smith bases[0] = 0; 4268865f1eaSKarl Rupp for (i=1; i<=size; i++) bases[i] = ldims[i-1]; 4278865f1eaSKarl Rupp for (i=1; i<=size; i++) bases[i] += bases[i-1]; 428ce00eea3SSatish Balay base = bases[rank]*dof; 42947c6ae99SBarry Smith 43047c6ae99SBarry Smith /* allocate the base parallel and sequential vectors */ 431ce00eea3SSatish Balay dd->Nlocal = x*y*z*dof; 432b1fb7eb7SBarry Smith ierr = VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,NULL,&global);CHKERRQ(ierr); 433ce00eea3SSatish Balay dd->nlocal = (Xe-Xs)*(Ye-Ys)*(Ze-Zs)*dof; 434b1fb7eb7SBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,NULL,&local);CHKERRQ(ierr); 43547c6ae99SBarry Smith 436*4104a7a0SPatrick Sanan /* generate global to local vector scatter and local to global mapping*/ 43747c6ae99SBarry Smith 438ce00eea3SSatish Balay /* global to local must include ghost points within the domain, 439ce00eea3SSatish Balay but not ghost points outside the domain that aren't periodic */ 440*4104a7a0SPatrick Sanan ierr = PetscMalloc1((IXe-IXs)*(IYe-IYs)*(IZe-IZs),&idx);CHKERRQ(ierr); 441d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_BOX) { 442ce00eea3SSatish Balay left = IXs - Xs; right = left + (IXe-IXs); 443ce00eea3SSatish Balay bottom = IYs - Ys; top = bottom + (IYe-IYs); 444ce00eea3SSatish Balay down = IZs - Zs; up = down + (IZe-IZs); 445ce00eea3SSatish Balay count = 0; 446ce00eea3SSatish Balay for (i=down; i<up; i++) { 447ce00eea3SSatish Balay for (j=bottom; j<top; j++) { 448ce00eea3SSatish Balay for (k=left; k<right; k++) { 449ce00eea3SSatish Balay idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 450ce00eea3SSatish Balay } 451ce00eea3SSatish Balay } 452ce00eea3SSatish Balay } 453ce00eea3SSatish Balay ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr); 45447c6ae99SBarry Smith } else { 45547c6ae99SBarry Smith /* This is way ugly! We need to list the funny cross type region */ 456ce00eea3SSatish Balay left = xs - Xs; right = left + x; 45747c6ae99SBarry Smith bottom = ys - Ys; top = bottom + y; 45847c6ae99SBarry Smith down = zs - Zs; up = down + z; 45947c6ae99SBarry Smith count = 0; 460ce00eea3SSatish Balay /* the bottom chunck */ 461ce00eea3SSatish Balay for (i=(IZs-Zs); i<down; i++) { 46247c6ae99SBarry Smith for (j=bottom; j<top; j++) { 463ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 46447c6ae99SBarry Smith } 46547c6ae99SBarry Smith } 46647c6ae99SBarry Smith /* the middle piece */ 46747c6ae99SBarry Smith for (i=down; i<up; i++) { 46847c6ae99SBarry Smith /* front */ 469ce00eea3SSatish Balay for (j=(IYs-Ys); j<bottom; j++) { 470ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 47147c6ae99SBarry Smith } 47247c6ae99SBarry Smith /* middle */ 47347c6ae99SBarry Smith for (j=bottom; j<top; j++) { 474ce00eea3SSatish Balay for (k=IXs-Xs; k<IXe-Xs; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 47547c6ae99SBarry Smith } 47647c6ae99SBarry Smith /* back */ 477ce00eea3SSatish Balay for (j=top; j<top+IYe-ye; j++) { 478ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 47947c6ae99SBarry Smith } 48047c6ae99SBarry Smith } 48147c6ae99SBarry Smith /* the top piece */ 482ce00eea3SSatish Balay for (i=up; i<up+IZe-ze; i++) { 48347c6ae99SBarry Smith for (j=bottom; j<top; j++) { 484ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 48547c6ae99SBarry Smith } 48647c6ae99SBarry Smith } 48747c6ae99SBarry Smith ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr); 48847c6ae99SBarry Smith } 48947c6ae99SBarry Smith 49047c6ae99SBarry Smith /* determine who lies on each side of use stored in n24 n25 n26 49147c6ae99SBarry Smith n21 n22 n23 49247c6ae99SBarry Smith n18 n19 n20 49347c6ae99SBarry Smith 49447c6ae99SBarry Smith n15 n16 n17 49547c6ae99SBarry Smith n12 n14 49647c6ae99SBarry Smith n9 n10 n11 49747c6ae99SBarry Smith 49847c6ae99SBarry Smith n6 n7 n8 49947c6ae99SBarry Smith n3 n4 n5 50047c6ae99SBarry Smith n0 n1 n2 50147c6ae99SBarry Smith */ 50247c6ae99SBarry Smith 50347c6ae99SBarry Smith /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */ 50447c6ae99SBarry Smith /* Assume Nodes are Internal to the Cube */ 50547c6ae99SBarry Smith n0 = rank - m*n - m - 1; 50647c6ae99SBarry Smith n1 = rank - m*n - m; 50747c6ae99SBarry Smith n2 = rank - m*n - m + 1; 50847c6ae99SBarry Smith n3 = rank - m*n -1; 50947c6ae99SBarry Smith n4 = rank - m*n; 51047c6ae99SBarry Smith n5 = rank - m*n + 1; 51147c6ae99SBarry Smith n6 = rank - m*n + m - 1; 51247c6ae99SBarry Smith n7 = rank - m*n + m; 51347c6ae99SBarry Smith n8 = rank - m*n + m + 1; 51447c6ae99SBarry Smith 51547c6ae99SBarry Smith n9 = rank - m - 1; 51647c6ae99SBarry Smith n10 = rank - m; 51747c6ae99SBarry Smith n11 = rank - m + 1; 51847c6ae99SBarry Smith n12 = rank - 1; 51947c6ae99SBarry Smith n14 = rank + 1; 52047c6ae99SBarry Smith n15 = rank + m - 1; 52147c6ae99SBarry Smith n16 = rank + m; 52247c6ae99SBarry Smith n17 = rank + m + 1; 52347c6ae99SBarry Smith 52447c6ae99SBarry Smith n18 = rank + m*n - m - 1; 52547c6ae99SBarry Smith n19 = rank + m*n - m; 52647c6ae99SBarry Smith n20 = rank + m*n - m + 1; 52747c6ae99SBarry Smith n21 = rank + m*n - 1; 52847c6ae99SBarry Smith n22 = rank + m*n; 52947c6ae99SBarry Smith n23 = rank + m*n + 1; 53047c6ae99SBarry Smith n24 = rank + m*n + m - 1; 53147c6ae99SBarry Smith n25 = rank + m*n + m; 53247c6ae99SBarry Smith n26 = rank + m*n + m + 1; 53347c6ae99SBarry Smith 53447c6ae99SBarry Smith /* Assume Pieces are on Faces of Cube */ 53547c6ae99SBarry Smith 53647c6ae99SBarry Smith if (xs == 0) { /* First assume not corner or edge */ 53747c6ae99SBarry Smith n0 = rank -1 - (m*n); 53847c6ae99SBarry Smith n3 = rank + m -1 - (m*n); 53947c6ae99SBarry Smith n6 = rank + 2*m -1 - (m*n); 54047c6ae99SBarry Smith n9 = rank -1; 54147c6ae99SBarry Smith n12 = rank + m -1; 54247c6ae99SBarry Smith n15 = rank + 2*m -1; 54347c6ae99SBarry Smith n18 = rank -1 + (m*n); 54447c6ae99SBarry Smith n21 = rank + m -1 + (m*n); 54547c6ae99SBarry Smith n24 = rank + 2*m -1 + (m*n); 54647c6ae99SBarry Smith } 54747c6ae99SBarry Smith 548ce00eea3SSatish Balay if (xe == M) { /* First assume not corner or edge */ 54947c6ae99SBarry Smith n2 = rank -2*m +1 - (m*n); 55047c6ae99SBarry Smith n5 = rank - m +1 - (m*n); 55147c6ae99SBarry Smith n8 = rank +1 - (m*n); 55247c6ae99SBarry Smith n11 = rank -2*m +1; 55347c6ae99SBarry Smith n14 = rank - m +1; 55447c6ae99SBarry Smith n17 = rank +1; 55547c6ae99SBarry Smith n20 = rank -2*m +1 + (m*n); 55647c6ae99SBarry Smith n23 = rank - m +1 + (m*n); 55747c6ae99SBarry Smith n26 = rank +1 + (m*n); 55847c6ae99SBarry Smith } 55947c6ae99SBarry Smith 56047c6ae99SBarry Smith if (ys==0) { /* First assume not corner or edge */ 56147c6ae99SBarry Smith n0 = rank + m * (n-1) -1 - (m*n); 56247c6ae99SBarry Smith n1 = rank + m * (n-1) - (m*n); 56347c6ae99SBarry Smith n2 = rank + m * (n-1) +1 - (m*n); 56447c6ae99SBarry Smith n9 = rank + m * (n-1) -1; 56547c6ae99SBarry Smith n10 = rank + m * (n-1); 56647c6ae99SBarry Smith n11 = rank + m * (n-1) +1; 56747c6ae99SBarry Smith n18 = rank + m * (n-1) -1 + (m*n); 56847c6ae99SBarry Smith n19 = rank + m * (n-1) + (m*n); 56947c6ae99SBarry Smith n20 = rank + m * (n-1) +1 + (m*n); 57047c6ae99SBarry Smith } 57147c6ae99SBarry Smith 57247c6ae99SBarry Smith if (ye == N) { /* First assume not corner or edge */ 57347c6ae99SBarry Smith n6 = rank - m * (n-1) -1 - (m*n); 57447c6ae99SBarry Smith n7 = rank - m * (n-1) - (m*n); 57547c6ae99SBarry Smith n8 = rank - m * (n-1) +1 - (m*n); 57647c6ae99SBarry Smith n15 = rank - m * (n-1) -1; 57747c6ae99SBarry Smith n16 = rank - m * (n-1); 57847c6ae99SBarry Smith n17 = rank - m * (n-1) +1; 57947c6ae99SBarry Smith n24 = rank - m * (n-1) -1 + (m*n); 58047c6ae99SBarry Smith n25 = rank - m * (n-1) + (m*n); 58147c6ae99SBarry Smith n26 = rank - m * (n-1) +1 + (m*n); 58247c6ae99SBarry Smith } 58347c6ae99SBarry Smith 58447c6ae99SBarry Smith if (zs == 0) { /* First assume not corner or edge */ 58547c6ae99SBarry Smith n0 = size - (m*n) + rank - m - 1; 58647c6ae99SBarry Smith n1 = size - (m*n) + rank - m; 58747c6ae99SBarry Smith n2 = size - (m*n) + rank - m + 1; 58847c6ae99SBarry Smith n3 = size - (m*n) + rank - 1; 58947c6ae99SBarry Smith n4 = size - (m*n) + rank; 59047c6ae99SBarry Smith n5 = size - (m*n) + rank + 1; 59147c6ae99SBarry Smith n6 = size - (m*n) + rank + m - 1; 59247c6ae99SBarry Smith n7 = size - (m*n) + rank + m; 59347c6ae99SBarry Smith n8 = size - (m*n) + rank + m + 1; 59447c6ae99SBarry Smith } 59547c6ae99SBarry Smith 59647c6ae99SBarry Smith if (ze == P) { /* First assume not corner or edge */ 59747c6ae99SBarry Smith n18 = (m*n) - (size-rank) - m - 1; 59847c6ae99SBarry Smith n19 = (m*n) - (size-rank) - m; 59947c6ae99SBarry Smith n20 = (m*n) - (size-rank) - m + 1; 60047c6ae99SBarry Smith n21 = (m*n) - (size-rank) - 1; 60147c6ae99SBarry Smith n22 = (m*n) - (size-rank); 60247c6ae99SBarry Smith n23 = (m*n) - (size-rank) + 1; 60347c6ae99SBarry Smith n24 = (m*n) - (size-rank) + m - 1; 60447c6ae99SBarry Smith n25 = (m*n) - (size-rank) + m; 60547c6ae99SBarry Smith n26 = (m*n) - (size-rank) + m + 1; 60647c6ae99SBarry Smith } 60747c6ae99SBarry Smith 60847c6ae99SBarry Smith if ((xs==0) && (zs==0)) { /* Assume an edge, not corner */ 60947c6ae99SBarry Smith n0 = size - m*n + rank + m-1 - m; 61047c6ae99SBarry Smith n3 = size - m*n + rank + m-1; 61147c6ae99SBarry Smith n6 = size - m*n + rank + m-1 + m; 61247c6ae99SBarry Smith } 61347c6ae99SBarry Smith 61447c6ae99SBarry Smith if ((xs==0) && (ze==P)) { /* Assume an edge, not corner */ 61547c6ae99SBarry Smith n18 = m*n - (size - rank) + m-1 - m; 61647c6ae99SBarry Smith n21 = m*n - (size - rank) + m-1; 61747c6ae99SBarry Smith n24 = m*n - (size - rank) + m-1 + m; 61847c6ae99SBarry Smith } 61947c6ae99SBarry Smith 62047c6ae99SBarry Smith if ((xs==0) && (ys==0)) { /* Assume an edge, not corner */ 62147c6ae99SBarry Smith n0 = rank + m*n -1 - m*n; 62247c6ae99SBarry Smith n9 = rank + m*n -1; 62347c6ae99SBarry Smith n18 = rank + m*n -1 + m*n; 62447c6ae99SBarry Smith } 62547c6ae99SBarry Smith 62647c6ae99SBarry Smith if ((xs==0) && (ye==N)) { /* Assume an edge, not corner */ 62747c6ae99SBarry Smith n6 = rank - m*(n-1) + m-1 - m*n; 62847c6ae99SBarry Smith n15 = rank - m*(n-1) + m-1; 62947c6ae99SBarry Smith n24 = rank - m*(n-1) + m-1 + m*n; 63047c6ae99SBarry Smith } 63147c6ae99SBarry Smith 632ce00eea3SSatish Balay if ((xe==M) && (zs==0)) { /* Assume an edge, not corner */ 63347c6ae99SBarry Smith n2 = size - (m*n-rank) - (m-1) - m; 63447c6ae99SBarry Smith n5 = size - (m*n-rank) - (m-1); 63547c6ae99SBarry Smith n8 = size - (m*n-rank) - (m-1) + m; 63647c6ae99SBarry Smith } 63747c6ae99SBarry Smith 638ce00eea3SSatish Balay if ((xe==M) && (ze==P)) { /* Assume an edge, not corner */ 63947c6ae99SBarry Smith n20 = m*n - (size - rank) - (m-1) - m; 64047c6ae99SBarry Smith n23 = m*n - (size - rank) - (m-1); 64147c6ae99SBarry Smith n26 = m*n - (size - rank) - (m-1) + m; 64247c6ae99SBarry Smith } 64347c6ae99SBarry Smith 644ce00eea3SSatish Balay if ((xe==M) && (ys==0)) { /* Assume an edge, not corner */ 64547c6ae99SBarry Smith n2 = rank + m*(n-1) - (m-1) - m*n; 64647c6ae99SBarry Smith n11 = rank + m*(n-1) - (m-1); 64747c6ae99SBarry Smith n20 = rank + m*(n-1) - (m-1) + m*n; 64847c6ae99SBarry Smith } 64947c6ae99SBarry Smith 650ce00eea3SSatish Balay if ((xe==M) && (ye==N)) { /* Assume an edge, not corner */ 65147c6ae99SBarry Smith n8 = rank - m*n +1 - m*n; 65247c6ae99SBarry Smith n17 = rank - m*n +1; 65347c6ae99SBarry Smith n26 = rank - m*n +1 + m*n; 65447c6ae99SBarry Smith } 65547c6ae99SBarry Smith 65647c6ae99SBarry Smith if ((ys==0) && (zs==0)) { /* Assume an edge, not corner */ 65747c6ae99SBarry Smith n0 = size - m + rank -1; 65847c6ae99SBarry Smith n1 = size - m + rank; 65947c6ae99SBarry Smith n2 = size - m + rank +1; 66047c6ae99SBarry Smith } 66147c6ae99SBarry Smith 66247c6ae99SBarry Smith if ((ys==0) && (ze==P)) { /* Assume an edge, not corner */ 66347c6ae99SBarry Smith n18 = m*n - (size - rank) + m*(n-1) -1; 66447c6ae99SBarry Smith n19 = m*n - (size - rank) + m*(n-1); 66547c6ae99SBarry Smith n20 = m*n - (size - rank) + m*(n-1) +1; 66647c6ae99SBarry Smith } 66747c6ae99SBarry Smith 66847c6ae99SBarry Smith if ((ye==N) && (zs==0)) { /* Assume an edge, not corner */ 66947c6ae99SBarry Smith n6 = size - (m*n-rank) - m * (n-1) -1; 67047c6ae99SBarry Smith n7 = size - (m*n-rank) - m * (n-1); 67147c6ae99SBarry Smith n8 = size - (m*n-rank) - m * (n-1) +1; 67247c6ae99SBarry Smith } 67347c6ae99SBarry Smith 67447c6ae99SBarry Smith if ((ye==N) && (ze==P)) { /* Assume an edge, not corner */ 67547c6ae99SBarry Smith n24 = rank - (size-m) -1; 67647c6ae99SBarry Smith n25 = rank - (size-m); 67747c6ae99SBarry Smith n26 = rank - (size-m) +1; 67847c6ae99SBarry Smith } 67947c6ae99SBarry Smith 68047c6ae99SBarry Smith /* Check for Corners */ 6818865f1eaSKarl Rupp if ((xs==0) && (ys==0) && (zs==0)) n0 = size -1; 6828865f1eaSKarl Rupp if ((xs==0) && (ys==0) && (ze==P)) n18 = m*n-1; 6838865f1eaSKarl Rupp if ((xs==0) && (ye==N) && (zs==0)) n6 = (size-1)-m*(n-1); 6848865f1eaSKarl Rupp if ((xs==0) && (ye==N) && (ze==P)) n24 = m-1; 6858865f1eaSKarl Rupp if ((xe==M) && (ys==0) && (zs==0)) n2 = size-m; 6868865f1eaSKarl Rupp if ((xe==M) && (ys==0) && (ze==P)) n20 = m*n-m; 6878865f1eaSKarl Rupp if ((xe==M) && (ye==N) && (zs==0)) n8 = size-m*n; 6888865f1eaSKarl Rupp if ((xe==M) && (ye==N) && (ze==P)) n26 = 0; 68947c6ae99SBarry Smith 69047c6ae99SBarry Smith /* Check for when not X,Y, and Z Periodic */ 69147c6ae99SBarry Smith 69247c6ae99SBarry Smith /* If not X periodic */ 693bff4a2f0SMatthew G. Knepley if (bx != DM_BOUNDARY_PERIODIC) { 6948865f1eaSKarl Rupp if (xs==0) n0 = n3 = n6 = n9 = n12 = n15 = n18 = n21 = n24 = -2; 6958865f1eaSKarl Rupp if (xe==M) n2 = n5 = n8 = n11 = n14 = n17 = n20 = n23 = n26 = -2; 69647c6ae99SBarry Smith } 69747c6ae99SBarry Smith 69847c6ae99SBarry Smith /* If not Y periodic */ 699bff4a2f0SMatthew G. Knepley if (by != DM_BOUNDARY_PERIODIC) { 7008865f1eaSKarl Rupp if (ys==0) n0 = n1 = n2 = n9 = n10 = n11 = n18 = n19 = n20 = -2; 7018865f1eaSKarl Rupp if (ye==N) n6 = n7 = n8 = n15 = n16 = n17 = n24 = n25 = n26 = -2; 70247c6ae99SBarry Smith } 70347c6ae99SBarry Smith 70447c6ae99SBarry Smith /* If not Z periodic */ 705bff4a2f0SMatthew G. Knepley if (bz != DM_BOUNDARY_PERIODIC) { 7068865f1eaSKarl Rupp if (zs==0) n0 = n1 = n2 = n3 = n4 = n5 = n6 = n7 = n8 = -2; 7078865f1eaSKarl Rupp if (ze==P) n18 = n19 = n20 = n21 = n22 = n23 = n24 = n25 = n26 = -2; 70847c6ae99SBarry Smith } 70947c6ae99SBarry Smith 710785e854fSJed Brown ierr = PetscMalloc1(27,&dd->neighbors);CHKERRQ(ierr); 7118865f1eaSKarl Rupp 71247c6ae99SBarry Smith dd->neighbors[0] = n0; 71347c6ae99SBarry Smith dd->neighbors[1] = n1; 71447c6ae99SBarry Smith dd->neighbors[2] = n2; 71547c6ae99SBarry Smith dd->neighbors[3] = n3; 71647c6ae99SBarry Smith dd->neighbors[4] = n4; 71747c6ae99SBarry Smith dd->neighbors[5] = n5; 71847c6ae99SBarry Smith dd->neighbors[6] = n6; 71947c6ae99SBarry Smith dd->neighbors[7] = n7; 72047c6ae99SBarry Smith dd->neighbors[8] = n8; 72147c6ae99SBarry Smith dd->neighbors[9] = n9; 72247c6ae99SBarry Smith dd->neighbors[10] = n10; 72347c6ae99SBarry Smith dd->neighbors[11] = n11; 72447c6ae99SBarry Smith dd->neighbors[12] = n12; 72547c6ae99SBarry Smith dd->neighbors[13] = rank; 72647c6ae99SBarry Smith dd->neighbors[14] = n14; 72747c6ae99SBarry Smith dd->neighbors[15] = n15; 72847c6ae99SBarry Smith dd->neighbors[16] = n16; 72947c6ae99SBarry Smith dd->neighbors[17] = n17; 73047c6ae99SBarry Smith dd->neighbors[18] = n18; 73147c6ae99SBarry Smith dd->neighbors[19] = n19; 73247c6ae99SBarry Smith dd->neighbors[20] = n20; 73347c6ae99SBarry Smith dd->neighbors[21] = n21; 73447c6ae99SBarry Smith dd->neighbors[22] = n22; 73547c6ae99SBarry Smith dd->neighbors[23] = n23; 73647c6ae99SBarry Smith dd->neighbors[24] = n24; 73747c6ae99SBarry Smith dd->neighbors[25] = n25; 73847c6ae99SBarry Smith dd->neighbors[26] = n26; 73947c6ae99SBarry Smith 74047c6ae99SBarry Smith /* If star stencil then delete the corner neighbors */ 741d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 74247c6ae99SBarry Smith /* save information about corner neighbors */ 74347c6ae99SBarry Smith sn0 = n0; sn1 = n1; sn2 = n2; sn3 = n3; sn5 = n5; sn6 = n6; sn7 = n7; 74447c6ae99SBarry Smith sn8 = n8; sn9 = n9; sn11 = n11; sn15 = n15; sn17 = n17; sn18 = n18; 74547c6ae99SBarry Smith sn19 = n19; sn20 = n20; sn21 = n21; sn23 = n23; sn24 = n24; sn25 = n25; 74647c6ae99SBarry Smith sn26 = n26; 7478865f1eaSKarl Rupp n0 = n1 = n2 = n3 = n5 = n6 = n7 = n8 = n9 = n11 = n15 = n17 = n18 = n19 = n20 = n21 = n23 = n24 = n25 = n26 = -1; 74847c6ae99SBarry Smith } 74947c6ae99SBarry Smith 750785e854fSJed Brown ierr = PetscMalloc1((Xe-Xs)*(Ye-Ys)*(Ze-Zs),&idx);CHKERRQ(ierr); 75147c6ae99SBarry Smith 75247c6ae99SBarry Smith nn = 0; 75347c6ae99SBarry Smith /* Bottom Level */ 75447c6ae99SBarry Smith for (k=0; k<s_z; k++) { 75547c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 75647c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 757ce00eea3SSatish Balay x_t = lx[n0 % m]; 75847c6ae99SBarry Smith y_t = ly[(n0 % (m*n))/m]; 75947c6ae99SBarry Smith z_t = lz[n0 / (m*n)]; 76047c6ae99SBarry 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; 7618865f1eaSKarl 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 */ 7628865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 76347c6ae99SBarry Smith } 76447c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 76547c6ae99SBarry Smith x_t = x; 76647c6ae99SBarry Smith y_t = ly[(n1 % (m*n))/m]; 76747c6ae99SBarry Smith z_t = lz[n1 / (m*n)]; 76847c6ae99SBarry 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; 7698865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t; /* 2D case */ 7708865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 77147c6ae99SBarry Smith } 77247c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 773ce00eea3SSatish Balay x_t = lx[n2 % m]; 77447c6ae99SBarry Smith y_t = ly[(n2 % (m*n))/m]; 77547c6ae99SBarry Smith z_t = lz[n2 / (m*n)]; 77647c6ae99SBarry 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; 7778865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t; /* 2D case */ 7788865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 77947c6ae99SBarry Smith } 78047c6ae99SBarry Smith } 78147c6ae99SBarry Smith 78247c6ae99SBarry Smith for (i=0; i<y; i++) { 78347c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 784ce00eea3SSatish Balay x_t = lx[n3 % m]; 78547c6ae99SBarry Smith y_t = y; 78647c6ae99SBarry Smith z_t = lz[n3 / (m*n)]; 78747c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 7888865f1eaSKarl 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 */ 7898865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 79047c6ae99SBarry Smith } 79147c6ae99SBarry Smith 79247c6ae99SBarry Smith if (n4 >= 0) { /* middle */ 79347c6ae99SBarry Smith x_t = x; 79447c6ae99SBarry Smith y_t = y; 79547c6ae99SBarry Smith z_t = lz[n4 / (m*n)]; 79647c6ae99SBarry Smith s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 7978865f1eaSKarl 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 */ 7988865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 799bff4a2f0SMatthew G. Knepley } else if (bz == DM_BOUNDARY_MIRROR) { 8008865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 80147c6ae99SBarry Smith } 80247c6ae99SBarry Smith 80347c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 804ce00eea3SSatish Balay x_t = lx[n5 % m]; 80547c6ae99SBarry Smith y_t = y; 80647c6ae99SBarry Smith z_t = lz[n5 / (m*n)]; 80747c6ae99SBarry Smith s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8088865f1eaSKarl 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 */ 8098865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 81047c6ae99SBarry Smith } 81147c6ae99SBarry Smith } 81247c6ae99SBarry Smith 81347c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 81447c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 815ce00eea3SSatish Balay x_t = lx[n6 % m]; 81647c6ae99SBarry Smith y_t = ly[(n6 % (m*n))/m]; 81747c6ae99SBarry Smith z_t = lz[n6 / (m*n)]; 81847c6ae99SBarry Smith s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8198865f1eaSKarl 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 */ 8208865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 82147c6ae99SBarry Smith } 82247c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 82347c6ae99SBarry Smith x_t = x; 82447c6ae99SBarry Smith y_t = ly[(n7 % (m*n))/m]; 82547c6ae99SBarry Smith z_t = lz[n7 / (m*n)]; 82647c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8278865f1eaSKarl 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 */ 8288865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 82947c6ae99SBarry Smith } 83047c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 831ce00eea3SSatish Balay x_t = lx[n8 % m]; 83247c6ae99SBarry Smith y_t = ly[(n8 % (m*n))/m]; 83347c6ae99SBarry Smith z_t = lz[n8 / (m*n)]; 83447c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8358865f1eaSKarl 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 */ 8368865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 83747c6ae99SBarry Smith } 83847c6ae99SBarry Smith } 83947c6ae99SBarry Smith } 84047c6ae99SBarry Smith 84147c6ae99SBarry Smith /* Middle Level */ 84247c6ae99SBarry Smith for (k=0; k<z; k++) { 84347c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 84447c6ae99SBarry Smith if (n9 >= 0) { /* left below */ 845ce00eea3SSatish Balay x_t = lx[n9 % m]; 84647c6ae99SBarry Smith y_t = ly[(n9 % (m*n))/m]; 84747c6ae99SBarry Smith /* z_t = z; */ 84847c6ae99SBarry Smith s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 8498865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 85047c6ae99SBarry Smith } 85147c6ae99SBarry Smith if (n10 >= 0) { /* directly below */ 85247c6ae99SBarry Smith x_t = x; 85347c6ae99SBarry Smith y_t = ly[(n10 % (m*n))/m]; 85447c6ae99SBarry Smith /* z_t = z; */ 85547c6ae99SBarry Smith s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 8568865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 857bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_MIRROR) { 8588865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 85947c6ae99SBarry Smith } 86047c6ae99SBarry Smith if (n11 >= 0) { /* right below */ 861ce00eea3SSatish Balay x_t = lx[n11 % m]; 86247c6ae99SBarry Smith y_t = ly[(n11 % (m*n))/m]; 86347c6ae99SBarry Smith /* z_t = z; */ 86447c6ae99SBarry Smith s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 8658865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 86647c6ae99SBarry Smith } 86747c6ae99SBarry Smith } 86847c6ae99SBarry Smith 86947c6ae99SBarry Smith for (i=0; i<y; i++) { 87047c6ae99SBarry Smith if (n12 >= 0) { /* directly left */ 871ce00eea3SSatish Balay x_t = lx[n12 % m]; 87247c6ae99SBarry Smith y_t = y; 87347c6ae99SBarry Smith /* z_t = z; */ 87447c6ae99SBarry Smith s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t; 8758865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 876bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { 8778865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = 0; 87847c6ae99SBarry Smith } 87947c6ae99SBarry Smith 88047c6ae99SBarry Smith /* Interior */ 88147c6ae99SBarry Smith s_t = bases[rank] + i*x + k*x*y; 8828865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = s_t++; 88347c6ae99SBarry Smith 88447c6ae99SBarry Smith if (n14 >= 0) { /* directly right */ 885ce00eea3SSatish Balay x_t = lx[n14 % m]; 88647c6ae99SBarry Smith y_t = y; 88747c6ae99SBarry Smith /* z_t = z; */ 88847c6ae99SBarry Smith s_t = bases[n14] + i*x_t + k*x_t*y_t; 8898865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 890bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { 8918865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = 0; 89247c6ae99SBarry Smith } 89347c6ae99SBarry Smith } 89447c6ae99SBarry Smith 89547c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 89647c6ae99SBarry Smith if (n15 >= 0) { /* left above */ 897ce00eea3SSatish Balay x_t = lx[n15 % m]; 89847c6ae99SBarry Smith y_t = ly[(n15 % (m*n))/m]; 89947c6ae99SBarry Smith /* z_t = z; */ 90047c6ae99SBarry Smith s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t; 9018865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 90247c6ae99SBarry Smith } 90347c6ae99SBarry Smith if (n16 >= 0) { /* directly above */ 90447c6ae99SBarry Smith x_t = x; 90547c6ae99SBarry Smith y_t = ly[(n16 % (m*n))/m]; 90647c6ae99SBarry Smith /* z_t = z; */ 90747c6ae99SBarry Smith s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t; 9088865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 909bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_MIRROR) { 9108865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 91147c6ae99SBarry Smith } 91247c6ae99SBarry Smith if (n17 >= 0) { /* right above */ 913ce00eea3SSatish Balay x_t = lx[n17 % m]; 91447c6ae99SBarry Smith y_t = ly[(n17 % (m*n))/m]; 91547c6ae99SBarry Smith /* z_t = z; */ 91647c6ae99SBarry Smith s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t; 9178865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 91847c6ae99SBarry Smith } 91947c6ae99SBarry Smith } 92047c6ae99SBarry Smith } 92147c6ae99SBarry Smith 92247c6ae99SBarry Smith /* Upper Level */ 92347c6ae99SBarry Smith for (k=0; k<s_z; k++) { 92447c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 92547c6ae99SBarry Smith if (n18 >= 0) { /* left below */ 926ce00eea3SSatish Balay x_t = lx[n18 % m]; 92747c6ae99SBarry Smith y_t = ly[(n18 % (m*n))/m]; 92847c6ae99SBarry Smith /* z_t = lz[n18 / (m*n)]; */ 92947c6ae99SBarry Smith s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 9308865f1eaSKarl 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 */ 9318865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 93247c6ae99SBarry Smith } 93347c6ae99SBarry Smith if (n19 >= 0) { /* directly below */ 93447c6ae99SBarry Smith x_t = x; 93547c6ae99SBarry Smith y_t = ly[(n19 % (m*n))/m]; 93647c6ae99SBarry Smith /* z_t = lz[n19 / (m*n)]; */ 93747c6ae99SBarry Smith s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 9388865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n19] - (s_y+1-i)*x_t + x_t*y_t; /* 2d case */ 9398865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 94047c6ae99SBarry Smith } 94147c6ae99SBarry Smith if (n20 >= 0) { /* right below */ 942ce00eea3SSatish Balay x_t = lx[n20 % m]; 94347c6ae99SBarry Smith y_t = ly[(n20 % (m*n))/m]; 94447c6ae99SBarry Smith /* z_t = lz[n20 / (m*n)]; */ 94547c6ae99SBarry Smith s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 9468865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n20] - (s_y+1-i)*x_t + x_t*y_t; /* 2d case */ 9478865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 94847c6ae99SBarry Smith } 94947c6ae99SBarry Smith } 95047c6ae99SBarry Smith 95147c6ae99SBarry Smith for (i=0; i<y; i++) { 95247c6ae99SBarry Smith if (n21 >= 0) { /* directly left */ 953ce00eea3SSatish Balay x_t = lx[n21 % m]; 95447c6ae99SBarry Smith y_t = y; 95547c6ae99SBarry Smith /* z_t = lz[n21 / (m*n)]; */ 95647c6ae99SBarry Smith s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t; 9578865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n21] + (i+1)*x_t - s_x; /* 2d case */ 9588865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 95947c6ae99SBarry Smith } 96047c6ae99SBarry Smith 96147c6ae99SBarry Smith if (n22 >= 0) { /* middle */ 96247c6ae99SBarry Smith x_t = x; 96347c6ae99SBarry Smith y_t = y; 96447c6ae99SBarry Smith /* z_t = lz[n22 / (m*n)]; */ 96547c6ae99SBarry Smith s_t = bases[n22] + i*x_t + k*x_t*y_t; 9668865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n22] + i*x_t; /* 2d case */ 9678865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 968bff4a2f0SMatthew G. Knepley } else if (bz == DM_BOUNDARY_MIRROR) { 9698865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 97047c6ae99SBarry Smith } 97147c6ae99SBarry Smith 97247c6ae99SBarry Smith if (n23 >= 0) { /* directly right */ 973ce00eea3SSatish Balay x_t = lx[n23 % m]; 97447c6ae99SBarry Smith y_t = y; 97547c6ae99SBarry Smith /* z_t = lz[n23 / (m*n)]; */ 97647c6ae99SBarry Smith s_t = bases[n23] + i*x_t + k*x_t*y_t; 9778865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n23] + i*x_t; /* 2d case */ 9788865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 97947c6ae99SBarry Smith } 98047c6ae99SBarry Smith } 98147c6ae99SBarry Smith 98247c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 98347c6ae99SBarry Smith if (n24 >= 0) { /* left above */ 984ce00eea3SSatish Balay x_t = lx[n24 % m]; 98547c6ae99SBarry Smith y_t = ly[(n24 % (m*n))/m]; 98647c6ae99SBarry Smith /* z_t = lz[n24 / (m*n)]; */ 98747c6ae99SBarry Smith s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t; 9888865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n24] + i*x_t - s_x; /* 2d case */ 9898865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 99047c6ae99SBarry Smith } 99147c6ae99SBarry Smith if (n25 >= 0) { /* directly above */ 99247c6ae99SBarry Smith x_t = x; 99347c6ae99SBarry Smith y_t = ly[(n25 % (m*n))/m]; 99447c6ae99SBarry Smith /* z_t = lz[n25 / (m*n)]; */ 99547c6ae99SBarry Smith s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t; 9968865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n25] + (i-1)*x_t; /* 2d case */ 9978865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 99847c6ae99SBarry Smith } 99947c6ae99SBarry Smith if (n26 >= 0) { /* right above */ 1000ce00eea3SSatish Balay x_t = lx[n26 % m]; 100147c6ae99SBarry Smith y_t = ly[(n26 % (m*n))/m]; 100247c6ae99SBarry Smith /* z_t = lz[n26 / (m*n)]; */ 100347c6ae99SBarry Smith s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t; 10048865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n26] + (i-1)*x_t; /* 2d case */ 10058865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 100647c6ae99SBarry Smith } 100747c6ae99SBarry Smith } 100847c6ae99SBarry Smith } 1009ce00eea3SSatish Balay 1010b1fb7eb7SBarry Smith ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_USE_POINTER,&from);CHKERRQ(ierr); 101147c6ae99SBarry Smith ierr = VecScatterCreate(global,from,local,to,>ol);CHKERRQ(ierr); 10123bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)gtol);CHKERRQ(ierr); 1013fcfd50ebSBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 1014fcfd50ebSBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 101547c6ae99SBarry Smith 1016d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 101747c6ae99SBarry Smith n0 = sn0; n1 = sn1; n2 = sn2; n3 = sn3; n5 = sn5; n6 = sn6; n7 = sn7; 101847c6ae99SBarry Smith n8 = sn8; n9 = sn9; n11 = sn11; n15 = sn15; n17 = sn17; n18 = sn18; 101947c6ae99SBarry Smith n19 = sn19; n20 = sn20; n21 = sn21; n23 = sn23; n24 = sn24; n25 = sn25; 102047c6ae99SBarry Smith n26 = sn26; 1021ce00eea3SSatish Balay } 102247c6ae99SBarry Smith 1023288e7d53SBarry Smith if (((stencil_type == DMDA_STENCIL_STAR) || (bx != DM_BOUNDARY_PERIODIC && bx) || (by != DM_BOUNDARY_PERIODIC && by) || (bz != DM_BOUNDARY_PERIODIC && bz))) { 1024ce00eea3SSatish Balay /* 1025ce00eea3SSatish Balay Recompute the local to global mappings, this time keeping the 1026ce00eea3SSatish Balay information about the cross corner processor numbers. 1027ce00eea3SSatish Balay */ 102847c6ae99SBarry Smith nn = 0; 102947c6ae99SBarry Smith /* Bottom Level */ 103047c6ae99SBarry Smith for (k=0; k<s_z; k++) { 103147c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 103247c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 1033ce00eea3SSatish Balay x_t = lx[n0 % m]; 103447c6ae99SBarry Smith y_t = ly[(n0 % (m*n))/m]; 103547c6ae99SBarry Smith z_t = lz[n0 / (m*n)]; 103647c6ae99SBarry 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; 10378865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1038ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0 && Zs-zs < 0) { 10398865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 104047c6ae99SBarry Smith } 104147c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 104247c6ae99SBarry Smith x_t = x; 104347c6ae99SBarry Smith y_t = ly[(n1 % (m*n))/m]; 104447c6ae99SBarry Smith z_t = lz[n1 / (m*n)]; 104547c6ae99SBarry 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; 10468865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1047ce00eea3SSatish Balay } else if (Ys-ys < 0 && Zs-zs < 0) { 10488865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 104947c6ae99SBarry Smith } 105047c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 1051ce00eea3SSatish Balay x_t = lx[n2 % m]; 105247c6ae99SBarry Smith y_t = ly[(n2 % (m*n))/m]; 105347c6ae99SBarry Smith z_t = lz[n2 / (m*n)]; 105447c6ae99SBarry 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; 10558865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1056ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0 && Zs-zs < 0) { 10578865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 105847c6ae99SBarry Smith } 105947c6ae99SBarry Smith } 106047c6ae99SBarry Smith 106147c6ae99SBarry Smith for (i=0; i<y; i++) { 106247c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 1063ce00eea3SSatish Balay x_t = lx[n3 % m]; 106447c6ae99SBarry Smith y_t = y; 106547c6ae99SBarry Smith z_t = lz[n3 / (m*n)]; 106647c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 10678865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1068ce00eea3SSatish Balay } else if (Xs-xs < 0 && Zs-zs < 0) { 10698865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 107047c6ae99SBarry Smith } 107147c6ae99SBarry Smith 107247c6ae99SBarry Smith if (n4 >= 0) { /* middle */ 107347c6ae99SBarry Smith x_t = x; 107447c6ae99SBarry Smith y_t = y; 107547c6ae99SBarry Smith z_t = lz[n4 / (m*n)]; 107647c6ae99SBarry Smith s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 10778865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1078ce00eea3SSatish Balay } else if (Zs-zs < 0) { 1079bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_MIRROR) { 10808865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 1081c2859e5eSBarry Smith } else { 10828865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 108347c6ae99SBarry Smith } 1084c2859e5eSBarry Smith } 108547c6ae99SBarry Smith 108647c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 1087ce00eea3SSatish Balay x_t = lx[n5 % m]; 108847c6ae99SBarry Smith y_t = y; 108947c6ae99SBarry Smith z_t = lz[n5 / (m*n)]; 109047c6ae99SBarry Smith s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 10918865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1092ce00eea3SSatish Balay } else if (xe-Xe < 0 && Zs-zs < 0) { 10938865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 109447c6ae99SBarry Smith } 109547c6ae99SBarry Smith } 109647c6ae99SBarry Smith 109747c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 109847c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 1099ce00eea3SSatish Balay x_t = lx[n6 % m]; 110047c6ae99SBarry Smith y_t = ly[(n6 % (m*n))/m]; 110147c6ae99SBarry Smith z_t = lz[n6 / (m*n)]; 110247c6ae99SBarry Smith s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11038865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1104ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0 && Zs-zs < 0) { 11058865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 110647c6ae99SBarry Smith } 110747c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 110847c6ae99SBarry Smith x_t = x; 110947c6ae99SBarry Smith y_t = ly[(n7 % (m*n))/m]; 111047c6ae99SBarry Smith z_t = lz[n7 / (m*n)]; 111147c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11128865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1113ce00eea3SSatish Balay } else if (ye-Ye < 0 && Zs-zs < 0) { 11148865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 111547c6ae99SBarry Smith } 111647c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 1117ce00eea3SSatish Balay x_t = lx[n8 % m]; 111847c6ae99SBarry Smith y_t = ly[(n8 % (m*n))/m]; 111947c6ae99SBarry Smith z_t = lz[n8 / (m*n)]; 112047c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11218865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1122ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0 && Zs-zs < 0) { 11238865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 112447c6ae99SBarry Smith } 112547c6ae99SBarry Smith } 112647c6ae99SBarry Smith } 112747c6ae99SBarry Smith 112847c6ae99SBarry Smith /* Middle Level */ 112947c6ae99SBarry Smith for (k=0; k<z; k++) { 113047c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 113147c6ae99SBarry Smith if (n9 >= 0) { /* left below */ 1132ce00eea3SSatish Balay x_t = lx[n9 % m]; 113347c6ae99SBarry Smith y_t = ly[(n9 % (m*n))/m]; 113447c6ae99SBarry Smith /* z_t = z; */ 113547c6ae99SBarry Smith s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 11368865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1137ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0) { 11388865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 113947c6ae99SBarry Smith } 114047c6ae99SBarry Smith if (n10 >= 0) { /* directly below */ 114147c6ae99SBarry Smith x_t = x; 114247c6ae99SBarry Smith y_t = ly[(n10 % (m*n))/m]; 114347c6ae99SBarry Smith /* z_t = z; */ 114447c6ae99SBarry Smith s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 11458865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1146ce00eea3SSatish Balay } else if (Ys-ys < 0) { 1147bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_MIRROR) { 11488865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 1149c2859e5eSBarry Smith } else { 11508865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 1151c2859e5eSBarry Smith } 115247c6ae99SBarry Smith } 115347c6ae99SBarry Smith if (n11 >= 0) { /* right below */ 1154ce00eea3SSatish Balay x_t = lx[n11 % m]; 115547c6ae99SBarry Smith y_t = ly[(n11 % (m*n))/m]; 115647c6ae99SBarry Smith /* z_t = z; */ 115747c6ae99SBarry Smith s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 11588865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1159ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0) { 11608865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 116147c6ae99SBarry Smith } 116247c6ae99SBarry Smith } 116347c6ae99SBarry Smith 116447c6ae99SBarry Smith for (i=0; i<y; i++) { 116547c6ae99SBarry Smith if (n12 >= 0) { /* directly left */ 1166ce00eea3SSatish Balay x_t = lx[n12 % m]; 116747c6ae99SBarry Smith y_t = y; 116847c6ae99SBarry Smith /* z_t = z; */ 116947c6ae99SBarry Smith s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t; 11708865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1171ce00eea3SSatish Balay } else if (Xs-xs < 0) { 1172bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_MIRROR) { 11738865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = 0; 1174c2859e5eSBarry Smith } else { 11758865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 117647c6ae99SBarry Smith } 1177c2859e5eSBarry Smith } 117847c6ae99SBarry Smith 117947c6ae99SBarry Smith /* Interior */ 118047c6ae99SBarry Smith s_t = bases[rank] + i*x + k*x*y; 11818865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = s_t++; 118247c6ae99SBarry Smith 118347c6ae99SBarry Smith if (n14 >= 0) { /* directly right */ 1184ce00eea3SSatish Balay x_t = lx[n14 % m]; 118547c6ae99SBarry Smith y_t = y; 118647c6ae99SBarry Smith /* z_t = z; */ 118747c6ae99SBarry Smith s_t = bases[n14] + i*x_t + k*x_t*y_t; 11888865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1189ce00eea3SSatish Balay } else if (xe-Xe < 0) { 1190bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_MIRROR) { 11918865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = 0; 1192c2859e5eSBarry Smith } else { 11938865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 119447c6ae99SBarry Smith } 119547c6ae99SBarry Smith } 1196c2859e5eSBarry Smith } 119747c6ae99SBarry Smith 119847c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 119947c6ae99SBarry Smith if (n15 >= 0) { /* left above */ 1200ce00eea3SSatish Balay x_t = lx[n15 % m]; 120147c6ae99SBarry Smith y_t = ly[(n15 % (m*n))/m]; 120247c6ae99SBarry Smith /* z_t = z; */ 120347c6ae99SBarry Smith s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t; 12048865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1205ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0) { 12068865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 120747c6ae99SBarry Smith } 120847c6ae99SBarry Smith if (n16 >= 0) { /* directly above */ 120947c6ae99SBarry Smith x_t = x; 121047c6ae99SBarry Smith y_t = ly[(n16 % (m*n))/m]; 121147c6ae99SBarry Smith /* z_t = z; */ 121247c6ae99SBarry Smith s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t; 12138865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1214ce00eea3SSatish Balay } else if (ye-Ye < 0) { 1215bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_MIRROR) { 12168865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 1217c2859e5eSBarry Smith } else { 12188865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 121947c6ae99SBarry Smith } 1220c2859e5eSBarry Smith } 122147c6ae99SBarry Smith if (n17 >= 0) { /* right above */ 1222ce00eea3SSatish Balay x_t = lx[n17 % m]; 122347c6ae99SBarry Smith y_t = ly[(n17 % (m*n))/m]; 122447c6ae99SBarry Smith /* z_t = z; */ 122547c6ae99SBarry Smith s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t; 12268865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1227ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0) { 12288865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 122947c6ae99SBarry Smith } 123047c6ae99SBarry Smith } 123147c6ae99SBarry Smith } 123247c6ae99SBarry Smith 123347c6ae99SBarry Smith /* Upper Level */ 123447c6ae99SBarry Smith for (k=0; k<s_z; k++) { 123547c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 123647c6ae99SBarry Smith if (n18 >= 0) { /* left below */ 1237ce00eea3SSatish Balay x_t = lx[n18 % m]; 123847c6ae99SBarry Smith y_t = ly[(n18 % (m*n))/m]; 123947c6ae99SBarry Smith /* z_t = lz[n18 / (m*n)]; */ 124047c6ae99SBarry Smith s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 12418865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1242ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0 && ze-Ze < 0) { 12438865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 124447c6ae99SBarry Smith } 124547c6ae99SBarry Smith if (n19 >= 0) { /* directly below */ 124647c6ae99SBarry Smith x_t = x; 124747c6ae99SBarry Smith y_t = ly[(n19 % (m*n))/m]; 124847c6ae99SBarry Smith /* z_t = lz[n19 / (m*n)]; */ 124947c6ae99SBarry Smith s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 12508865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1251ce00eea3SSatish Balay } else if (Ys-ys < 0 && ze-Ze < 0) { 12528865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 125347c6ae99SBarry Smith } 125447c6ae99SBarry Smith if (n20 >= 0) { /* right below */ 1255ce00eea3SSatish Balay x_t = lx[n20 % m]; 125647c6ae99SBarry Smith y_t = ly[(n20 % (m*n))/m]; 125747c6ae99SBarry Smith /* z_t = lz[n20 / (m*n)]; */ 125847c6ae99SBarry Smith s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 12598865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1260ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0 && ze-Ze < 0) { 12618865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 126247c6ae99SBarry Smith } 126347c6ae99SBarry Smith } 126447c6ae99SBarry Smith 126547c6ae99SBarry Smith for (i=0; i<y; i++) { 126647c6ae99SBarry Smith if (n21 >= 0) { /* directly left */ 1267ce00eea3SSatish Balay x_t = lx[n21 % m]; 126847c6ae99SBarry Smith y_t = y; 126947c6ae99SBarry Smith /* z_t = lz[n21 / (m*n)]; */ 127047c6ae99SBarry Smith s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t; 12718865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1272ce00eea3SSatish Balay } else if (Xs-xs < 0 && ze-Ze < 0) { 12738865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 127447c6ae99SBarry Smith } 127547c6ae99SBarry Smith 127647c6ae99SBarry Smith if (n22 >= 0) { /* middle */ 127747c6ae99SBarry Smith x_t = x; 127847c6ae99SBarry Smith y_t = y; 127947c6ae99SBarry Smith /* z_t = lz[n22 / (m*n)]; */ 128047c6ae99SBarry Smith s_t = bases[n22] + i*x_t + k*x_t*y_t; 12818865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1282ce00eea3SSatish Balay } else if (ze-Ze < 0) { 1283bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_MIRROR) { 12848865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 1285c2859e5eSBarry Smith } else { 12868865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 128747c6ae99SBarry Smith } 1288c2859e5eSBarry Smith } 128947c6ae99SBarry Smith 129047c6ae99SBarry Smith if (n23 >= 0) { /* directly right */ 1291ce00eea3SSatish Balay x_t = lx[n23 % m]; 129247c6ae99SBarry Smith y_t = y; 129347c6ae99SBarry Smith /* z_t = lz[n23 / (m*n)]; */ 129447c6ae99SBarry Smith s_t = bases[n23] + i*x_t + k*x_t*y_t; 12958865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1296ce00eea3SSatish Balay } else if (xe-Xe < 0 && ze-Ze < 0) { 12978865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 129847c6ae99SBarry Smith } 129947c6ae99SBarry Smith } 130047c6ae99SBarry Smith 130147c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 130247c6ae99SBarry Smith if (n24 >= 0) { /* left above */ 1303ce00eea3SSatish Balay x_t = lx[n24 % m]; 130447c6ae99SBarry Smith y_t = ly[(n24 % (m*n))/m]; 130547c6ae99SBarry Smith /* z_t = lz[n24 / (m*n)]; */ 130647c6ae99SBarry Smith s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t; 13078865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1308ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0 && ze-Ze < 0) { 13098865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 131047c6ae99SBarry Smith } 131147c6ae99SBarry Smith if (n25 >= 0) { /* directly above */ 131247c6ae99SBarry Smith x_t = x; 131347c6ae99SBarry Smith y_t = ly[(n25 % (m*n))/m]; 131447c6ae99SBarry Smith /* z_t = lz[n25 / (m*n)]; */ 131547c6ae99SBarry Smith s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t; 13168865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1317ce00eea3SSatish Balay } else if (ye-Ye < 0 && ze-Ze < 0) { 13188865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 131947c6ae99SBarry Smith } 132047c6ae99SBarry Smith if (n26 >= 0) { /* right above */ 1321ce00eea3SSatish Balay x_t = lx[n26 % m]; 132247c6ae99SBarry Smith y_t = ly[(n26 % (m*n))/m]; 132347c6ae99SBarry Smith /* z_t = lz[n26 / (m*n)]; */ 132447c6ae99SBarry Smith s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t; 13258865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1326ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0 && ze-Ze < 0) { 13278865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 132847c6ae99SBarry Smith } 132947c6ae99SBarry Smith } 133047c6ae99SBarry Smith } 133147c6ae99SBarry Smith } 133247c6ae99SBarry Smith /* 133347c6ae99SBarry Smith Set the local to global ordering in the global vector, this allows use 133447c6ae99SBarry Smith of VecSetValuesLocal(). 133547c6ae99SBarry Smith */ 133645b6f7e9SBarry Smith ierr = ISLocalToGlobalMappingCreate(comm,dof,nn,idx,PETSC_OWN_POINTER,&da->ltogmap);CHKERRQ(ierr); 13373bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap);CHKERRQ(ierr); 133847c6ae99SBarry Smith 1339ce00eea3SSatish Balay ierr = PetscFree2(bases,ldims);CHKERRQ(ierr); 1340ce00eea3SSatish Balay dd->m = m; dd->n = n; dd->p = p; 1341ce00eea3SSatish Balay /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */ 1342ce00eea3SSatish Balay dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = zs; dd->ze = ze; 1343ce00eea3SSatish Balay dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = Zs; dd->Ze = Ze; 1344ce00eea3SSatish Balay 1345fcfd50ebSBarry Smith ierr = VecDestroy(&local);CHKERRQ(ierr); 1346fcfd50ebSBarry Smith ierr = VecDestroy(&global);CHKERRQ(ierr); 1347ce00eea3SSatish Balay 1348ce00eea3SSatish Balay dd->gtol = gtol; 1349ce00eea3SSatish Balay dd->base = base; 1350ce00eea3SSatish Balay da->ops->view = DMView_DA_3d; 13510298fd71SBarry Smith dd->ltol = NULL; 13520298fd71SBarry Smith dd->ao = NULL; 135347c6ae99SBarry Smith PetscFunctionReturn(0); 135447c6ae99SBarry Smith } 1355564755cdSBarry Smith 135647c6ae99SBarry Smith 135747c6ae99SBarry Smith /*@C 1358aa219208SBarry Smith DMDACreate3d - Creates an object that will manage the communication of three-dimensional 135947c6ae99SBarry Smith regular array data that is distributed across some processors. 136047c6ae99SBarry Smith 136147c6ae99SBarry Smith Collective on MPI_Comm 136247c6ae99SBarry Smith 136347c6ae99SBarry Smith Input Parameters: 136447c6ae99SBarry Smith + comm - MPI communicator 13651321219cSEthan Coon . bx,by,bz - type of ghost nodes the array have. 1366bff4a2f0SMatthew G. Knepley Use one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC. 1367aa219208SBarry Smith . stencil_type - Type of stencil (DMDA_STENCIL_STAR or DMDA_STENCIL_BOX) 1368897f7067SBarry Smith . M,N,P - global dimension in each direction of the array 136947c6ae99SBarry Smith . m,n,p - corresponding number of processors in each dimension 137047c6ae99SBarry Smith (or PETSC_DECIDE to have calculated) 137147c6ae99SBarry Smith . dof - number of degrees of freedom per node 137210d7c030SMatthew G Knepley . s - stencil width 137310d7c030SMatthew G Knepley - lx, ly, lz - arrays containing the number of nodes in each cell along 13740298fd71SBarry Smith the x, y, and z coordinates, or NULL. If non-null, these 137547c6ae99SBarry Smith must be of length as m,n,p and the corresponding 137647c6ae99SBarry Smith m,n, or p cannot be PETSC_DECIDE. Sum of the lx[] entries must be M, sum of 137747c6ae99SBarry Smith the ly[] must N, sum of the lz[] must be P 137847c6ae99SBarry Smith 137947c6ae99SBarry Smith Output Parameter: 138047c6ae99SBarry Smith . da - the resulting distributed array object 138147c6ae99SBarry Smith 138247c6ae99SBarry Smith Options Database Key: 1383706a11cbSBarry Smith + -dm_view - Calls DMView() at the conclusion of DMDACreate3d() 138447c6ae99SBarry Smith . -da_grid_x <nx> - number of grid points in x direction, if M < 0 138547c6ae99SBarry Smith . -da_grid_y <ny> - number of grid points in y direction, if N < 0 138647c6ae99SBarry Smith . -da_grid_z <nz> - number of grid points in z direction, if P < 0 138747c6ae99SBarry Smith . -da_processors_x <MX> - number of processors in x direction 138847c6ae99SBarry Smith . -da_processors_y <MY> - number of processors in y direction 138947c6ae99SBarry Smith . -da_processors_z <MZ> - number of processors in z direction 1390e0f5d30fSBarry Smith . -da_refine_x <rx> - refinement ratio in x direction 1391e0f5d30fSBarry Smith . -da_refine_y <ry> - refinement ratio in y direction 1392e0f5d30fSBarry Smith . -da_refine_z <rz>- refinement ratio in z directio 1393e0f5d30fSBarry Smith - -da_refine <n> - refine the DMDA n times before creating it, , if M, N, or P < 0 139447c6ae99SBarry Smith 139547c6ae99SBarry Smith Level: beginner 139647c6ae99SBarry Smith 139747c6ae99SBarry Smith Notes: 1398aa219208SBarry Smith The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the 1399aa219208SBarry Smith standard 7-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes 140047c6ae99SBarry Smith the standard 27-pt stencil. 140147c6ae99SBarry Smith 1402aa219208SBarry Smith The array data itself is NOT stored in the DMDA, it is stored in Vec objects; 1403564755cdSBarry Smith The appropriate vector objects can be obtained with calls to DMCreateGlobalVector() 1404564755cdSBarry Smith and DMCreateLocalVector() and calls to VecDuplicate() if more are needed. 140547c6ae99SBarry Smith 1406897f7067SBarry Smith You must call DMSetUp() after this call before using this DM. 1407897f7067SBarry Smith 1408897f7067SBarry Smith If you wish to use the options database to change values in the DMDA call DMSetFromOptions() after this call 1409897f7067SBarry Smith but before DMSetUp(). 1410897f7067SBarry Smith 141147c6ae99SBarry Smith .keywords: distributed array, create, three-dimensional 141247c6ae99SBarry Smith 1413aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate2d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(), 141499f0b487SRichard Tran Mills DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDASetRefinementFactor(), 1415d461ba97SBarry Smith DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges() 141647c6ae99SBarry Smith 141747c6ae99SBarry Smith @*/ 1418bff4a2f0SMatthew G. Knepley PetscErrorCode DMDACreate3d(MPI_Comm comm,DMBoundaryType bx,DMBoundaryType by,DMBoundaryType bz,DMDAStencilType stencil_type,PetscInt M, 14199a42bb27SBarry 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) 142047c6ae99SBarry Smith { 142147c6ae99SBarry Smith PetscErrorCode ierr; 142247c6ae99SBarry Smith 142347c6ae99SBarry Smith PetscFunctionBegin; 1424aa219208SBarry Smith ierr = DMDACreate(comm, da);CHKERRQ(ierr); 1425c73cfb54SMatthew G. Knepley ierr = DMSetDimension(*da, 3);CHKERRQ(ierr); 1426aa219208SBarry Smith ierr = DMDASetSizes(*da, M, N, P);CHKERRQ(ierr); 1427aa219208SBarry Smith ierr = DMDASetNumProcs(*da, m, n, p);CHKERRQ(ierr); 14281321219cSEthan Coon ierr = DMDASetBoundaryType(*da, bx, by, bz);CHKERRQ(ierr); 1429aa219208SBarry Smith ierr = DMDASetDof(*da, dof);CHKERRQ(ierr); 1430aa219208SBarry Smith ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr); 1431aa219208SBarry Smith ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr); 1432aa219208SBarry Smith ierr = DMDASetOwnershipRanges(*da, lx, ly, lz);CHKERRQ(ierr); 143347c6ae99SBarry Smith PetscFunctionReturn(0); 143447c6ae99SBarry Smith } 1435