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 74035e84dSBarry 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 347b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);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; 4947c6ae99SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Lower left corner %G %G %G : Upper right %G %G %G\n",coors[0],coors[1],coors[2],coors[last],coors[last+1],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); 547b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);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; 6247c6ae99SBarry Smith PetscInt k,plane,base,*idx; 6347c6ae99SBarry Smith char node[10]; 6447c6ae99SBarry Smith PetscBool isnull; 6547c6ae99SBarry Smith 6647c6ae99SBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 6747c6ae99SBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 6847c6ae99SBarry Smith ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr); 6947c6ae99SBarry Smith ierr = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr); 7047c6ae99SBarry Smith 7147c6ae99SBarry Smith /* first processor draw all node lines */ 7247c6ae99SBarry Smith if (!rank) { 7347c6ae99SBarry Smith for (k=0; k<dd->P; k++) { 7447c6ae99SBarry Smith ymin = 0.0; ymax = (PetscReal)(dd->N - 1); 7547c6ae99SBarry Smith for (xmin=(PetscReal)(k*(dd->M+1)); xmin<(PetscReal)(dd->M+(k*(dd->M+1))); xmin++) { 7647c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr); 7747c6ae99SBarry Smith } 7847c6ae99SBarry Smith 7947c6ae99SBarry Smith xmin = (PetscReal)(k*(dd->M+1)); xmax = xmin + (PetscReal)(dd->M - 1); 8047c6ae99SBarry Smith for (ymin=0; ymin<(PetscReal)dd->N; ymin++) { 8147c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr); 8247c6ae99SBarry Smith } 8347c6ae99SBarry Smith } 8447c6ae99SBarry Smith } 8547c6ae99SBarry Smith ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr); 8647c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 8747c6ae99SBarry Smith 8847c6ae99SBarry Smith for (k=0; k<dd->P; k++) { /*Go through and draw for each plane*/ 8947c6ae99SBarry Smith if ((k >= dd->zs) && (k < dd->ze)) { 9047c6ae99SBarry Smith /* draw my box */ 9147c6ae99SBarry Smith ymin = dd->ys; 9247c6ae99SBarry Smith ymax = dd->ye - 1; 9347c6ae99SBarry Smith xmin = dd->xs/dd->w + (dd->M+1)*k; 9447c6ae99SBarry Smith xmax =(dd->xe-1)/dd->w + (dd->M+1)*k; 9547c6ae99SBarry Smith 9647c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr); 9747c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 9847c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 9947c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 10047c6ae99SBarry Smith 10147c6ae99SBarry Smith xmin = dd->xs/dd->w; 10247c6ae99SBarry Smith xmax =(dd->xe-1)/dd->w; 10347c6ae99SBarry Smith 10447c6ae99SBarry Smith /* put in numbers*/ 10547c6ae99SBarry Smith base = (dd->base+(dd->xe-dd->xs)*(dd->ye-dd->ys)*(k-dd->zs))/dd->w; 10647c6ae99SBarry Smith 10747c6ae99SBarry Smith /* Identify which processor owns the box */ 10847c6ae99SBarry Smith sprintf(node,"%d",rank); 10947c6ae99SBarry Smith ierr = PetscDrawString(draw,xmin+(dd->M+1)*k+.2,ymin+.3,PETSC_DRAW_RED,node);CHKERRQ(ierr); 11047c6ae99SBarry Smith 11147c6ae99SBarry Smith for (y=ymin; y<=ymax; y++) { 11247c6ae99SBarry Smith for (x=xmin+(dd->M+1)*k; x<=xmax+(dd->M+1)*k; x++) { 11347c6ae99SBarry Smith sprintf(node,"%d",(int)base++); 11447c6ae99SBarry Smith ierr = PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);CHKERRQ(ierr); 11547c6ae99SBarry Smith } 11647c6ae99SBarry Smith } 11747c6ae99SBarry Smith 11847c6ae99SBarry Smith } 11947c6ae99SBarry Smith } 12047c6ae99SBarry Smith ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr); 12147c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 12247c6ae99SBarry Smith 12347c6ae99SBarry Smith for (k=0-dd->s; k<dd->P+dd->s; k++) { 12447c6ae99SBarry Smith /* Go through and draw for each plane */ 12547c6ae99SBarry Smith if ((k >= dd->Zs) && (k < dd->Ze)) { 12647c6ae99SBarry Smith 12747c6ae99SBarry Smith /* overlay ghost numbers, useful for error checking */ 12847c6ae99SBarry Smith base = (dd->Xe-dd->Xs)*(dd->Ye-dd->Ys)*(k-dd->Zs); idx = dd->idx; 12947c6ae99SBarry Smith plane=k; 13047c6ae99SBarry Smith /* Keep z wrap around points on the dradrawg */ 1318865f1eaSKarl Rupp if (k<0) plane=dd->P+k; 1328865f1eaSKarl Rupp if (k>=dd->P) plane=k-dd->P; 13347c6ae99SBarry Smith ymin = dd->Ys; ymax = dd->Ye; 13447c6ae99SBarry Smith xmin = (dd->M+1)*plane*dd->w; 13547c6ae99SBarry Smith xmax = (dd->M+1)*plane*dd->w+dd->M*dd->w; 13647c6ae99SBarry Smith for (y=ymin; y<ymax; y++) { 13747c6ae99SBarry Smith for (x=xmin+dd->Xs; x<xmin+dd->Xe; x+=dd->w) { 13847c6ae99SBarry Smith sprintf(node,"%d",(int)(idx[base]/dd->w)); 13947c6ae99SBarry Smith ycoord = y; 14047c6ae99SBarry Smith /*Keep y wrap around points on drawing */ 1418865f1eaSKarl Rupp if (y<0) ycoord = dd->N+y; 14247c6ae99SBarry Smith 1438865f1eaSKarl Rupp if (y>=dd->N) ycoord = y-dd->N; 14447c6ae99SBarry Smith xcoord = x; /* Keep x wrap points on drawing */ 14547c6ae99SBarry Smith 1468865f1eaSKarl Rupp if (x<xmin) xcoord = xmax - (xmin-x); 1478865f1eaSKarl Rupp if (x>=xmax) xcoord = xmin + (x-xmax); 14847c6ae99SBarry Smith ierr = PetscDrawString(draw,xcoord/dd->w,ycoord,PETSC_DRAW_BLUE,node);CHKERRQ(ierr); 14947c6ae99SBarry Smith base+=dd->w; 15047c6ae99SBarry Smith } 15147c6ae99SBarry Smith } 15247c6ae99SBarry Smith } 15347c6ae99SBarry Smith } 15447c6ae99SBarry Smith ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr); 15547c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 1569a42bb27SBarry Smith } else if (isbinary) { 1579a42bb27SBarry Smith ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr); 1589a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 1599a42bb27SBarry Smith } else if (ismatlab) { 1609a42bb27SBarry Smith ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr); 1619a42bb27SBarry Smith #endif 16211aeaf0aSBarry Smith } 16347c6ae99SBarry Smith PetscFunctionReturn(0); 16447c6ae99SBarry Smith } 16547c6ae99SBarry Smith 16647c6ae99SBarry Smith #undef __FUNCT__ 1679a42bb27SBarry Smith #define __FUNCT__ "DMSetUp_DA_3D" 1687087cfbeSBarry Smith PetscErrorCode DMSetUp_DA_3D(DM da) 16947c6ae99SBarry Smith { 17047c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 17147c6ae99SBarry Smith const PetscInt M = dd->M; 17247c6ae99SBarry Smith const PetscInt N = dd->N; 17347c6ae99SBarry Smith const PetscInt P = dd->P; 17447c6ae99SBarry Smith PetscInt m = dd->m; 17547c6ae99SBarry Smith PetscInt n = dd->n; 17647c6ae99SBarry Smith PetscInt p = dd->p; 17747c6ae99SBarry Smith const PetscInt dof = dd->w; 17847c6ae99SBarry Smith const PetscInt s = dd->s; 17919fd82e9SBarry Smith DMDABoundaryType bx = dd->bx; 18019fd82e9SBarry Smith DMDABoundaryType by = dd->by; 18119fd82e9SBarry Smith DMDABoundaryType bz = dd->bz; 18219fd82e9SBarry Smith DMDAStencilType stencil_type = dd->stencil_type; 18347c6ae99SBarry Smith PetscInt *lx = dd->lx; 18447c6ae99SBarry Smith PetscInt *ly = dd->ly; 18547c6ae99SBarry Smith PetscInt *lz = dd->lz; 18647c6ae99SBarry Smith MPI_Comm comm; 18747c6ae99SBarry Smith PetscMPIInt rank,size; 188ce00eea3SSatish Balay PetscInt xs = 0,xe,ys = 0,ye,zs = 0,ze,x = 0,y = 0,z = 0; 189ce00eea3SSatish Balay PetscInt Xs,Xe,Ys,Ye,Zs,Ze,IXs,IXe,IYs,IYe,IZs,IZe,start,end,pm; 190ce00eea3SSatish Balay PetscInt left,right,up,down,bottom,top,i,j,k,*idx,*idx_cpy,nn; 191ce00eea3SSatish Balay const PetscInt *idx_full; 19247c6ae99SBarry Smith PetscInt n0,n1,n2,n3,n4,n5,n6,n7,n8,n9,n10,n11,n12,n14; 19347c6ae99SBarry Smith PetscInt n15,n16,n17,n18,n19,n20,n21,n22,n23,n24,n25,n26; 194db87c5ecSEthan Coon PetscInt *bases,*ldims,base,x_t,y_t,z_t,s_t,count,s_x,s_y,s_z; 19547c6ae99SBarry Smith PetscInt sn0 = 0,sn1 = 0,sn2 = 0,sn3 = 0,sn5 = 0,sn6 = 0,sn7 = 0; 19647c6ae99SBarry Smith PetscInt sn8 = 0,sn9 = 0,sn11 = 0,sn15 = 0,sn24 = 0,sn25 = 0,sn26 = 0; 19747c6ae99SBarry Smith PetscInt sn17 = 0,sn18 = 0,sn19 = 0,sn20 = 0,sn21 = 0,sn23 = 0; 19847c6ae99SBarry Smith Vec local,global; 19947c6ae99SBarry Smith VecScatter ltog,gtol; 200ce00eea3SSatish Balay IS to,from,ltogis; 2016f951b95Secoon PetscBool twod; 20247c6ae99SBarry Smith PetscErrorCode ierr; 20347c6ae99SBarry Smith 2046f951b95Secoon 20547c6ae99SBarry Smith PetscFunctionBegin; 206ce94432eSBarry Smith if (stencil_type == DMDA_STENCIL_BOX && (bx == DMDA_BOUNDARY_MIRROR || by == DMDA_BOUNDARY_MIRROR || bz == DMDA_BOUNDARY_MIRROR)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Mirror boundary and box stencil"); 20747c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject) da, &comm);CHKERRQ(ierr); 2083855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES) 2093855c12bSBarry 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); 2103855c12bSBarry Smith #endif 2113855c12bSBarry Smith 21247c6ae99SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 21347c6ae99SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 21447c6ae99SBarry Smith 21547c6ae99SBarry Smith if (m != PETSC_DECIDE) { 21647c6ae99SBarry Smith if (m < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m); 21747c6ae99SBarry Smith else if (m > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size); 21847c6ae99SBarry Smith } 21947c6ae99SBarry Smith if (n != PETSC_DECIDE) { 22047c6ae99SBarry Smith if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n); 22147c6ae99SBarry Smith else if (n > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size); 22247c6ae99SBarry Smith } 22347c6ae99SBarry Smith if (p != PETSC_DECIDE) { 22447c6ae99SBarry Smith if (p < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Z direction: %D",p); 22547c6ae99SBarry Smith else if (p > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Z direction: %D %d",p,size); 22647c6ae99SBarry Smith } 2270716a85fSBarry 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); 22847c6ae99SBarry Smith 22947c6ae99SBarry Smith /* Partition the array among the processors */ 23047c6ae99SBarry Smith if (m == PETSC_DECIDE && n != PETSC_DECIDE && p != PETSC_DECIDE) { 23147c6ae99SBarry Smith m = size/(n*p); 23247c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) { 23347c6ae99SBarry Smith n = size/(m*p); 23447c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) { 23547c6ae99SBarry Smith p = size/(m*n); 23647c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) { 23747c6ae99SBarry Smith /* try for squarish distribution */ 238369cc0aeSBarry Smith m = (int)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)N*p))); 23947c6ae99SBarry Smith if (!m) m = 1; 24047c6ae99SBarry Smith while (m > 0) { 24147c6ae99SBarry Smith n = size/(m*p); 24247c6ae99SBarry Smith if (m*n*p == size) break; 24347c6ae99SBarry Smith m--; 24447c6ae99SBarry Smith } 24547c6ae99SBarry Smith if (!m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad p value: p = %D",p); 24647c6ae99SBarry Smith if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;} 24747c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) { 24847c6ae99SBarry Smith /* try for squarish distribution */ 249369cc0aeSBarry Smith m = (int)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n))); 25047c6ae99SBarry Smith if (!m) m = 1; 25147c6ae99SBarry Smith while (m > 0) { 25247c6ae99SBarry Smith p = size/(m*n); 25347c6ae99SBarry Smith if (m*n*p == size) break; 25447c6ae99SBarry Smith m--; 25547c6ae99SBarry Smith } 25647c6ae99SBarry Smith if (!m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad n value: n = %D",n); 25747c6ae99SBarry Smith if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;} 25847c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) { 25947c6ae99SBarry Smith /* try for squarish distribution */ 260369cc0aeSBarry Smith n = (int)(0.5 + PetscSqrtReal(((PetscReal)N)*((PetscReal)size)/((PetscReal)P*m))); 26147c6ae99SBarry Smith if (!n) n = 1; 26247c6ae99SBarry Smith while (n > 0) { 26347c6ae99SBarry Smith p = size/(m*n); 26447c6ae99SBarry Smith if (m*n*p == size) break; 26547c6ae99SBarry Smith n--; 26647c6ae99SBarry Smith } 26747c6ae99SBarry Smith if (!n) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad m value: m = %D",n); 26847c6ae99SBarry Smith if (N > P && n < p) {PetscInt _n = n; n = p; p = _n;} 26947c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) { 27047c6ae99SBarry Smith /* try for squarish distribution */ 271369cc0aeSBarry Smith n = (PetscInt)(0.5 + PetscPowReal(((PetscReal)N*N)*((PetscReal)size)/((PetscReal)P*M),(PetscReal)(1./3.))); 27247c6ae99SBarry Smith if (!n) n = 1; 27347c6ae99SBarry Smith while (n > 0) { 27447c6ae99SBarry Smith pm = size/n; 27547c6ae99SBarry Smith if (n*pm == size) break; 27647c6ae99SBarry Smith n--; 27747c6ae99SBarry Smith } 27847c6ae99SBarry Smith if (!n) n = 1; 279369cc0aeSBarry Smith m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n))); 28047c6ae99SBarry Smith if (!m) m = 1; 28147c6ae99SBarry Smith while (m > 0) { 28247c6ae99SBarry Smith p = size/(m*n); 28347c6ae99SBarry Smith if (m*n*p == size) break; 28447c6ae99SBarry Smith m--; 28547c6ae99SBarry Smith } 28647c6ae99SBarry Smith if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;} 287ce94432eSBarry Smith } else if (m*n*p != size) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition"); 28847c6ae99SBarry Smith 289ce94432eSBarry Smith if (m*n*p != size) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_PLIB,"Could not find good partition"); 290ce94432eSBarry Smith if (M < m) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m); 291ce94432eSBarry Smith if (N < n) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n); 292ce94432eSBarry Smith if (P < p) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Partition in z direction is too fine! %D %D",P,p); 29347c6ae99SBarry Smith 29447c6ae99SBarry Smith /* 29547c6ae99SBarry Smith Determine locally owned region 29647c6ae99SBarry Smith [x, y, or z]s is the first local node number, [x, y, z] is the number of local nodes 29747c6ae99SBarry Smith */ 29847c6ae99SBarry Smith 29947c6ae99SBarry Smith if (!lx) { 30047c6ae99SBarry Smith ierr = PetscMalloc(m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr); 30147c6ae99SBarry Smith lx = dd->lx; 3028865f1eaSKarl Rupp for (i=0; i<m; i++) lx[i] = M/m + ((M % m) > (i % m)); 30347c6ae99SBarry Smith } 30447c6ae99SBarry Smith x = lx[rank % m]; 30547c6ae99SBarry Smith xs = 0; 3068865f1eaSKarl Rupp for (i=0; i<(rank%m); i++) xs += lx[i]; 307d9c9ebe5SPeter Brune if ((x < s) && ((m > 1) || (bx == DMDA_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); 30847c6ae99SBarry Smith 30947c6ae99SBarry Smith if (!ly) { 31047c6ae99SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr); 31147c6ae99SBarry Smith ly = dd->ly; 3128865f1eaSKarl Rupp for (i=0; i<n; i++) ly[i] = N/n + ((N % n) > (i % n)); 31347c6ae99SBarry Smith } 31447c6ae99SBarry Smith y = ly[(rank % (m*n))/m]; 315d9c9ebe5SPeter Brune if ((y < s) && ((n > 1) || (by == DMDA_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); 31630729d88SBarry Smith 31747c6ae99SBarry Smith ys = 0; 3188865f1eaSKarl Rupp for (i=0; i<(rank % (m*n))/m; i++) ys += ly[i]; 31947c6ae99SBarry Smith 32047c6ae99SBarry Smith if (!lz) { 32147c6ae99SBarry Smith ierr = PetscMalloc(p*sizeof(PetscInt), &dd->lz);CHKERRQ(ierr); 32247c6ae99SBarry Smith lz = dd->lz; 3238865f1eaSKarl Rupp for (i=0; i<p; i++) lz[i] = P/p + ((P % p) > (i % p)); 32447c6ae99SBarry Smith } 32547c6ae99SBarry Smith z = lz[rank/(m*n)]; 326bcea557cSEthan Coon 327fdc81ce8SEthan Coon /* note this is different than x- and y-, as we will handle as an important special 328fdc81ce8SEthan Coon case when p=P=1 and DMDA_BOUNDARY_PERIODIC and s > z. This is to deal with 2D problems 329fdc81ce8SEthan Coon in a 3D code. Additional code for this case is noted with "2d case" comments */ 3306f951b95Secoon twod = PETSC_FALSE; 3318865f1eaSKarl Rupp if (P == 1) twod = PETSC_TRUE; 3328865f1eaSKarl Rupp else if ((z < s) && ((p > 1) || (bz == DMDA_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); 33347c6ae99SBarry Smith zs = 0; 3348865f1eaSKarl Rupp for (i=0; i<(rank/(m*n)); i++) zs += lz[i]; 33547c6ae99SBarry Smith ye = ys + y; 33647c6ae99SBarry Smith xe = xs + x; 33747c6ae99SBarry Smith ze = zs + z; 33847c6ae99SBarry Smith 33988661749SPeter Brune /* determine ghost region (Xs) and region scattered into (IXs) */ 340d9c9ebe5SPeter Brune if (xs-s > 0) { 341d9c9ebe5SPeter Brune Xs = xs - s; IXs = xs - s; 34288661749SPeter Brune } else { 3438865f1eaSKarl Rupp if (bx) Xs = xs - s; 3448865f1eaSKarl Rupp else Xs = 0; 34588661749SPeter Brune IXs = 0; 34688661749SPeter Brune } 347d9c9ebe5SPeter Brune if (xe+s <= M) { 348d9c9ebe5SPeter Brune Xe = xe + s; IXe = xe + s; 34988661749SPeter Brune } else { 35088661749SPeter Brune if (bx) { 351d9c9ebe5SPeter Brune Xs = xs - s; Xe = xe + s; 3528865f1eaSKarl Rupp } else Xe = M; 35388661749SPeter Brune IXe = M; 35488661749SPeter Brune } 35547c6ae99SBarry Smith 356c2859e5eSBarry Smith if (bx == DMDA_BOUNDARY_PERIODIC || bx == DMDA_BOUNDARY_MIRROR) { 357d9c9ebe5SPeter Brune IXs = xs - s; 358d9c9ebe5SPeter Brune IXe = xe + s; 359d9c9ebe5SPeter Brune Xs = xs - s; 360d9c9ebe5SPeter Brune Xe = xe + s; 36188661749SPeter Brune } 36288661749SPeter Brune 363d9c9ebe5SPeter Brune if (ys-s > 0) { 364d9c9ebe5SPeter Brune Ys = ys - s; IYs = ys - s; 36588661749SPeter Brune } else { 3668865f1eaSKarl Rupp if (by) Ys = ys - s; 3678865f1eaSKarl Rupp else Ys = 0; 36888661749SPeter Brune IYs = 0; 36988661749SPeter Brune } 370d9c9ebe5SPeter Brune if (ye+s <= N) { 371d9c9ebe5SPeter Brune Ye = ye + s; IYe = ye + s; 37288661749SPeter Brune } else { 3738865f1eaSKarl Rupp if (by) Ye = ye + s; 3748865f1eaSKarl Rupp else Ye = N; 37588661749SPeter Brune IYe = N; 37688661749SPeter Brune } 37788661749SPeter Brune 378c2859e5eSBarry Smith if (by == DMDA_BOUNDARY_PERIODIC || by == DMDA_BOUNDARY_MIRROR) { 379d9c9ebe5SPeter Brune IYs = ys - s; 380d9c9ebe5SPeter Brune IYe = ye + s; 381d9c9ebe5SPeter Brune Ys = ys - s; 382d9c9ebe5SPeter Brune Ye = ye + s; 38388661749SPeter Brune } 38488661749SPeter Brune 385d9c9ebe5SPeter Brune if (zs-s > 0) { 386d9c9ebe5SPeter Brune Zs = zs - s; IZs = zs - s; 38788661749SPeter Brune } else { 3888865f1eaSKarl Rupp if (bz) Zs = zs - s; 3898865f1eaSKarl Rupp else Zs = 0; 39088661749SPeter Brune IZs = 0; 39188661749SPeter Brune } 392d9c9ebe5SPeter Brune if (ze+s <= P) { 393d9c9ebe5SPeter Brune Ze = ze + s; IZe = ze + s; 39488661749SPeter Brune } else { 3958865f1eaSKarl Rupp if (bz) Ze = ze + s; 3968865f1eaSKarl Rupp else Ze = P; 39788661749SPeter Brune IZe = P; 39888661749SPeter Brune } 39988661749SPeter Brune 400c2859e5eSBarry Smith if (bz == DMDA_BOUNDARY_PERIODIC || bz == DMDA_BOUNDARY_MIRROR) { 401d9c9ebe5SPeter Brune IZs = zs - s; 402d9c9ebe5SPeter Brune IZe = ze + s; 403d9c9ebe5SPeter Brune Zs = zs - s; 404d9c9ebe5SPeter Brune Ze = ze + s; 40588661749SPeter Brune } 40647c6ae99SBarry Smith 40747c6ae99SBarry Smith /* Resize all X parameters to reflect w */ 408d9c9ebe5SPeter Brune s_x = s; 409d9c9ebe5SPeter Brune s_y = s; 410d9c9ebe5SPeter Brune s_z = s; 41147c6ae99SBarry Smith 41247c6ae99SBarry Smith /* determine starting point of each processor */ 41347c6ae99SBarry Smith nn = x*y*z; 41447c6ae99SBarry Smith ierr = PetscMalloc2(size+1,PetscInt,&bases,size,PetscInt,&ldims);CHKERRQ(ierr); 41547c6ae99SBarry Smith ierr = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr); 41647c6ae99SBarry Smith bases[0] = 0; 4178865f1eaSKarl Rupp for (i=1; i<=size; i++) bases[i] = ldims[i-1]; 4188865f1eaSKarl Rupp for (i=1; i<=size; i++) bases[i] += bases[i-1]; 419ce00eea3SSatish Balay base = bases[rank]*dof; 42047c6ae99SBarry Smith 42147c6ae99SBarry Smith /* allocate the base parallel and sequential vectors */ 422ce00eea3SSatish Balay dd->Nlocal = x*y*z*dof; 423778a2246SBarry Smith ierr = VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,0,&global);CHKERRQ(ierr); 424ce00eea3SSatish Balay dd->nlocal = (Xe-Xs)*(Ye-Ys)*(Ze-Zs)*dof; 425778a2246SBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,0,&local);CHKERRQ(ierr); 42647c6ae99SBarry Smith 42747c6ae99SBarry Smith /* generate appropriate vector scatters */ 42847c6ae99SBarry Smith /* local to global inserts non-ghost point region into global */ 42947c6ae99SBarry Smith ierr = VecGetOwnershipRange(global,&start,&end);CHKERRQ(ierr); 430ce00eea3SSatish Balay ierr = ISCreateStride(comm,x*y*z*dof,start,1,&to);CHKERRQ(ierr); 43147c6ae99SBarry Smith 432ce00eea3SSatish Balay ierr = PetscMalloc(x*y*z*sizeof(PetscInt),&idx);CHKERRQ(ierr); 433ce00eea3SSatish Balay left = xs - Xs; right = left + x; 43447c6ae99SBarry Smith bottom = ys - Ys; top = bottom + y; 43547c6ae99SBarry Smith down = zs - Zs; up = down + z; 43647c6ae99SBarry Smith count = 0; 43747c6ae99SBarry Smith for (i=down; i<up; i++) { 43847c6ae99SBarry Smith for (j=bottom; j<top; j++) { 439ce00eea3SSatish Balay for (k=left; k<right; k++) { 440ce00eea3SSatish Balay idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 44147c6ae99SBarry Smith } 44247c6ae99SBarry Smith } 44347c6ae99SBarry Smith } 44447c6ae99SBarry Smith 44547c6ae99SBarry Smith ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&from);CHKERRQ(ierr); 44647c6ae99SBarry Smith ierr = VecScatterCreate(local,from,global,to,<og);CHKERRQ(ierr); 447*3bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)ltog);CHKERRQ(ierr); 448fcfd50ebSBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 449fcfd50ebSBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 45047c6ae99SBarry Smith 451ce00eea3SSatish Balay /* global to local must include ghost points within the domain, 452ce00eea3SSatish Balay but not ghost points outside the domain that aren't periodic */ 453d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_BOX) { 454db87c5ecSEthan Coon count = (IXe-IXs)*(IYe-IYs)*(IZe-IZs); 455db87c5ecSEthan Coon ierr = PetscMalloc(count*sizeof(PetscInt),&idx);CHKERRQ(ierr); 456ce00eea3SSatish Balay 457ce00eea3SSatish Balay left = IXs - Xs; right = left + (IXe-IXs); 458ce00eea3SSatish Balay bottom = IYs - Ys; top = bottom + (IYe-IYs); 459ce00eea3SSatish Balay down = IZs - Zs; up = down + (IZe-IZs); 460ce00eea3SSatish Balay count = 0; 461ce00eea3SSatish Balay for (i=down; i<up; i++) { 462ce00eea3SSatish Balay for (j=bottom; j<top; j++) { 463ce00eea3SSatish Balay for (k=left; k<right; k++) { 464ce00eea3SSatish Balay idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 465ce00eea3SSatish Balay } 466ce00eea3SSatish Balay } 467ce00eea3SSatish Balay } 468ce00eea3SSatish Balay ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr); 469ce00eea3SSatish Balay 47047c6ae99SBarry Smith } else { 47147c6ae99SBarry Smith /* This is way ugly! We need to list the funny cross type region */ 472db87c5ecSEthan Coon count = ((ys-IYs) + (IYe-ye))*x*z + ((xs-IXs) + (IXe-xe))*y*z + ((zs-IZs) + (IZe-ze))*x*y + x*y*z; 473db87c5ecSEthan Coon ierr = PetscMalloc(count*sizeof(PetscInt),&idx);CHKERRQ(ierr); 474ce00eea3SSatish Balay 475ce00eea3SSatish Balay left = xs - Xs; right = left + x; 47647c6ae99SBarry Smith bottom = ys - Ys; top = bottom + y; 47747c6ae99SBarry Smith down = zs - Zs; up = down + z; 47847c6ae99SBarry Smith count = 0; 479ce00eea3SSatish Balay /* the bottom chunck */ 480ce00eea3SSatish Balay for (i=(IZs-Zs); i<down; i++) { 48147c6ae99SBarry Smith for (j=bottom; j<top; j++) { 482ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 48347c6ae99SBarry Smith } 48447c6ae99SBarry Smith } 48547c6ae99SBarry Smith /* the middle piece */ 48647c6ae99SBarry Smith for (i=down; i<up; i++) { 48747c6ae99SBarry Smith /* front */ 488ce00eea3SSatish Balay for (j=(IYs-Ys); j<bottom; j++) { 489ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 49047c6ae99SBarry Smith } 49147c6ae99SBarry Smith /* middle */ 49247c6ae99SBarry Smith for (j=bottom; j<top; j++) { 493ce00eea3SSatish Balay for (k=IXs-Xs; k<IXe-Xs; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 49447c6ae99SBarry Smith } 49547c6ae99SBarry Smith /* back */ 496ce00eea3SSatish Balay for (j=top; j<top+IYe-ye; j++) { 497ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 49847c6ae99SBarry Smith } 49947c6ae99SBarry Smith } 50047c6ae99SBarry Smith /* the top piece */ 501ce00eea3SSatish Balay for (i=up; i<up+IZe-ze; i++) { 50247c6ae99SBarry Smith for (j=bottom; j<top; j++) { 503ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 50447c6ae99SBarry Smith } 50547c6ae99SBarry Smith } 50647c6ae99SBarry Smith ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr); 50747c6ae99SBarry Smith } 50847c6ae99SBarry Smith 50947c6ae99SBarry Smith /* determine who lies on each side of use stored in n24 n25 n26 51047c6ae99SBarry Smith n21 n22 n23 51147c6ae99SBarry Smith n18 n19 n20 51247c6ae99SBarry Smith 51347c6ae99SBarry Smith n15 n16 n17 51447c6ae99SBarry Smith n12 n14 51547c6ae99SBarry Smith n9 n10 n11 51647c6ae99SBarry Smith 51747c6ae99SBarry Smith n6 n7 n8 51847c6ae99SBarry Smith n3 n4 n5 51947c6ae99SBarry Smith n0 n1 n2 52047c6ae99SBarry Smith */ 52147c6ae99SBarry Smith 52247c6ae99SBarry Smith /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */ 52347c6ae99SBarry Smith /* Assume Nodes are Internal to the Cube */ 52447c6ae99SBarry Smith n0 = rank - m*n - m - 1; 52547c6ae99SBarry Smith n1 = rank - m*n - m; 52647c6ae99SBarry Smith n2 = rank - m*n - m + 1; 52747c6ae99SBarry Smith n3 = rank - m*n -1; 52847c6ae99SBarry Smith n4 = rank - m*n; 52947c6ae99SBarry Smith n5 = rank - m*n + 1; 53047c6ae99SBarry Smith n6 = rank - m*n + m - 1; 53147c6ae99SBarry Smith n7 = rank - m*n + m; 53247c6ae99SBarry Smith n8 = rank - m*n + m + 1; 53347c6ae99SBarry Smith 53447c6ae99SBarry Smith n9 = rank - m - 1; 53547c6ae99SBarry Smith n10 = rank - m; 53647c6ae99SBarry Smith n11 = rank - m + 1; 53747c6ae99SBarry Smith n12 = rank - 1; 53847c6ae99SBarry Smith n14 = rank + 1; 53947c6ae99SBarry Smith n15 = rank + m - 1; 54047c6ae99SBarry Smith n16 = rank + m; 54147c6ae99SBarry Smith n17 = rank + m + 1; 54247c6ae99SBarry Smith 54347c6ae99SBarry Smith n18 = rank + m*n - m - 1; 54447c6ae99SBarry Smith n19 = rank + m*n - m; 54547c6ae99SBarry Smith n20 = rank + m*n - m + 1; 54647c6ae99SBarry Smith n21 = rank + m*n - 1; 54747c6ae99SBarry Smith n22 = rank + m*n; 54847c6ae99SBarry Smith n23 = rank + m*n + 1; 54947c6ae99SBarry Smith n24 = rank + m*n + m - 1; 55047c6ae99SBarry Smith n25 = rank + m*n + m; 55147c6ae99SBarry Smith n26 = rank + m*n + m + 1; 55247c6ae99SBarry Smith 55347c6ae99SBarry Smith /* Assume Pieces are on Faces of Cube */ 55447c6ae99SBarry Smith 55547c6ae99SBarry Smith if (xs == 0) { /* First assume not corner or edge */ 55647c6ae99SBarry Smith n0 = rank -1 - (m*n); 55747c6ae99SBarry Smith n3 = rank + m -1 - (m*n); 55847c6ae99SBarry Smith n6 = rank + 2*m -1 - (m*n); 55947c6ae99SBarry Smith n9 = rank -1; 56047c6ae99SBarry Smith n12 = rank + m -1; 56147c6ae99SBarry Smith n15 = rank + 2*m -1; 56247c6ae99SBarry Smith n18 = rank -1 + (m*n); 56347c6ae99SBarry Smith n21 = rank + m -1 + (m*n); 56447c6ae99SBarry Smith n24 = rank + 2*m -1 + (m*n); 56547c6ae99SBarry Smith } 56647c6ae99SBarry Smith 567ce00eea3SSatish Balay if (xe == M) { /* First assume not corner or edge */ 56847c6ae99SBarry Smith n2 = rank -2*m +1 - (m*n); 56947c6ae99SBarry Smith n5 = rank - m +1 - (m*n); 57047c6ae99SBarry Smith n8 = rank +1 - (m*n); 57147c6ae99SBarry Smith n11 = rank -2*m +1; 57247c6ae99SBarry Smith n14 = rank - m +1; 57347c6ae99SBarry Smith n17 = rank +1; 57447c6ae99SBarry Smith n20 = rank -2*m +1 + (m*n); 57547c6ae99SBarry Smith n23 = rank - m +1 + (m*n); 57647c6ae99SBarry Smith n26 = rank +1 + (m*n); 57747c6ae99SBarry Smith } 57847c6ae99SBarry Smith 57947c6ae99SBarry Smith if (ys==0) { /* First assume not corner or edge */ 58047c6ae99SBarry Smith n0 = rank + m * (n-1) -1 - (m*n); 58147c6ae99SBarry Smith n1 = rank + m * (n-1) - (m*n); 58247c6ae99SBarry Smith n2 = rank + m * (n-1) +1 - (m*n); 58347c6ae99SBarry Smith n9 = rank + m * (n-1) -1; 58447c6ae99SBarry Smith n10 = rank + m * (n-1); 58547c6ae99SBarry Smith n11 = rank + m * (n-1) +1; 58647c6ae99SBarry Smith n18 = rank + m * (n-1) -1 + (m*n); 58747c6ae99SBarry Smith n19 = rank + m * (n-1) + (m*n); 58847c6ae99SBarry Smith n20 = rank + m * (n-1) +1 + (m*n); 58947c6ae99SBarry Smith } 59047c6ae99SBarry Smith 59147c6ae99SBarry Smith if (ye == N) { /* First assume not corner or edge */ 59247c6ae99SBarry Smith n6 = rank - m * (n-1) -1 - (m*n); 59347c6ae99SBarry Smith n7 = rank - m * (n-1) - (m*n); 59447c6ae99SBarry Smith n8 = rank - m * (n-1) +1 - (m*n); 59547c6ae99SBarry Smith n15 = rank - m * (n-1) -1; 59647c6ae99SBarry Smith n16 = rank - m * (n-1); 59747c6ae99SBarry Smith n17 = rank - m * (n-1) +1; 59847c6ae99SBarry Smith n24 = rank - m * (n-1) -1 + (m*n); 59947c6ae99SBarry Smith n25 = rank - m * (n-1) + (m*n); 60047c6ae99SBarry Smith n26 = rank - m * (n-1) +1 + (m*n); 60147c6ae99SBarry Smith } 60247c6ae99SBarry Smith 60347c6ae99SBarry Smith if (zs == 0) { /* First assume not corner or edge */ 60447c6ae99SBarry Smith n0 = size - (m*n) + rank - m - 1; 60547c6ae99SBarry Smith n1 = size - (m*n) + rank - m; 60647c6ae99SBarry Smith n2 = size - (m*n) + rank - m + 1; 60747c6ae99SBarry Smith n3 = size - (m*n) + rank - 1; 60847c6ae99SBarry Smith n4 = size - (m*n) + rank; 60947c6ae99SBarry Smith n5 = size - (m*n) + rank + 1; 61047c6ae99SBarry Smith n6 = size - (m*n) + rank + m - 1; 61147c6ae99SBarry Smith n7 = size - (m*n) + rank + m; 61247c6ae99SBarry Smith n8 = size - (m*n) + rank + m + 1; 61347c6ae99SBarry Smith } 61447c6ae99SBarry Smith 61547c6ae99SBarry Smith if (ze == P) { /* First assume not corner or edge */ 61647c6ae99SBarry Smith n18 = (m*n) - (size-rank) - m - 1; 61747c6ae99SBarry Smith n19 = (m*n) - (size-rank) - m; 61847c6ae99SBarry Smith n20 = (m*n) - (size-rank) - m + 1; 61947c6ae99SBarry Smith n21 = (m*n) - (size-rank) - 1; 62047c6ae99SBarry Smith n22 = (m*n) - (size-rank); 62147c6ae99SBarry Smith n23 = (m*n) - (size-rank) + 1; 62247c6ae99SBarry Smith n24 = (m*n) - (size-rank) + m - 1; 62347c6ae99SBarry Smith n25 = (m*n) - (size-rank) + m; 62447c6ae99SBarry Smith n26 = (m*n) - (size-rank) + m + 1; 62547c6ae99SBarry Smith } 62647c6ae99SBarry Smith 62747c6ae99SBarry Smith if ((xs==0) && (zs==0)) { /* Assume an edge, not corner */ 62847c6ae99SBarry Smith n0 = size - m*n + rank + m-1 - m; 62947c6ae99SBarry Smith n3 = size - m*n + rank + m-1; 63047c6ae99SBarry Smith n6 = size - m*n + rank + m-1 + m; 63147c6ae99SBarry Smith } 63247c6ae99SBarry Smith 63347c6ae99SBarry Smith if ((xs==0) && (ze==P)) { /* Assume an edge, not corner */ 63447c6ae99SBarry Smith n18 = m*n - (size - rank) + m-1 - m; 63547c6ae99SBarry Smith n21 = m*n - (size - rank) + m-1; 63647c6ae99SBarry Smith n24 = m*n - (size - rank) + m-1 + m; 63747c6ae99SBarry Smith } 63847c6ae99SBarry Smith 63947c6ae99SBarry Smith if ((xs==0) && (ys==0)) { /* Assume an edge, not corner */ 64047c6ae99SBarry Smith n0 = rank + m*n -1 - m*n; 64147c6ae99SBarry Smith n9 = rank + m*n -1; 64247c6ae99SBarry Smith n18 = rank + m*n -1 + m*n; 64347c6ae99SBarry Smith } 64447c6ae99SBarry Smith 64547c6ae99SBarry Smith if ((xs==0) && (ye==N)) { /* Assume an edge, not corner */ 64647c6ae99SBarry Smith n6 = rank - m*(n-1) + m-1 - m*n; 64747c6ae99SBarry Smith n15 = rank - m*(n-1) + m-1; 64847c6ae99SBarry Smith n24 = rank - m*(n-1) + m-1 + m*n; 64947c6ae99SBarry Smith } 65047c6ae99SBarry Smith 651ce00eea3SSatish Balay if ((xe==M) && (zs==0)) { /* Assume an edge, not corner */ 65247c6ae99SBarry Smith n2 = size - (m*n-rank) - (m-1) - m; 65347c6ae99SBarry Smith n5 = size - (m*n-rank) - (m-1); 65447c6ae99SBarry Smith n8 = size - (m*n-rank) - (m-1) + m; 65547c6ae99SBarry Smith } 65647c6ae99SBarry Smith 657ce00eea3SSatish Balay if ((xe==M) && (ze==P)) { /* Assume an edge, not corner */ 65847c6ae99SBarry Smith n20 = m*n - (size - rank) - (m-1) - m; 65947c6ae99SBarry Smith n23 = m*n - (size - rank) - (m-1); 66047c6ae99SBarry Smith n26 = m*n - (size - rank) - (m-1) + m; 66147c6ae99SBarry Smith } 66247c6ae99SBarry Smith 663ce00eea3SSatish Balay if ((xe==M) && (ys==0)) { /* Assume an edge, not corner */ 66447c6ae99SBarry Smith n2 = rank + m*(n-1) - (m-1) - m*n; 66547c6ae99SBarry Smith n11 = rank + m*(n-1) - (m-1); 66647c6ae99SBarry Smith n20 = rank + m*(n-1) - (m-1) + m*n; 66747c6ae99SBarry Smith } 66847c6ae99SBarry Smith 669ce00eea3SSatish Balay if ((xe==M) && (ye==N)) { /* Assume an edge, not corner */ 67047c6ae99SBarry Smith n8 = rank - m*n +1 - m*n; 67147c6ae99SBarry Smith n17 = rank - m*n +1; 67247c6ae99SBarry Smith n26 = rank - m*n +1 + m*n; 67347c6ae99SBarry Smith } 67447c6ae99SBarry Smith 67547c6ae99SBarry Smith if ((ys==0) && (zs==0)) { /* Assume an edge, not corner */ 67647c6ae99SBarry Smith n0 = size - m + rank -1; 67747c6ae99SBarry Smith n1 = size - m + rank; 67847c6ae99SBarry Smith n2 = size - m + rank +1; 67947c6ae99SBarry Smith } 68047c6ae99SBarry Smith 68147c6ae99SBarry Smith if ((ys==0) && (ze==P)) { /* Assume an edge, not corner */ 68247c6ae99SBarry Smith n18 = m*n - (size - rank) + m*(n-1) -1; 68347c6ae99SBarry Smith n19 = m*n - (size - rank) + m*(n-1); 68447c6ae99SBarry Smith n20 = m*n - (size - rank) + m*(n-1) +1; 68547c6ae99SBarry Smith } 68647c6ae99SBarry Smith 68747c6ae99SBarry Smith if ((ye==N) && (zs==0)) { /* Assume an edge, not corner */ 68847c6ae99SBarry Smith n6 = size - (m*n-rank) - m * (n-1) -1; 68947c6ae99SBarry Smith n7 = size - (m*n-rank) - m * (n-1); 69047c6ae99SBarry Smith n8 = size - (m*n-rank) - m * (n-1) +1; 69147c6ae99SBarry Smith } 69247c6ae99SBarry Smith 69347c6ae99SBarry Smith if ((ye==N) && (ze==P)) { /* Assume an edge, not corner */ 69447c6ae99SBarry Smith n24 = rank - (size-m) -1; 69547c6ae99SBarry Smith n25 = rank - (size-m); 69647c6ae99SBarry Smith n26 = rank - (size-m) +1; 69747c6ae99SBarry Smith } 69847c6ae99SBarry Smith 69947c6ae99SBarry Smith /* Check for Corners */ 7008865f1eaSKarl Rupp if ((xs==0) && (ys==0) && (zs==0)) n0 = size -1; 7018865f1eaSKarl Rupp if ((xs==0) && (ys==0) && (ze==P)) n18 = m*n-1; 7028865f1eaSKarl Rupp if ((xs==0) && (ye==N) && (zs==0)) n6 = (size-1)-m*(n-1); 7038865f1eaSKarl Rupp if ((xs==0) && (ye==N) && (ze==P)) n24 = m-1; 7048865f1eaSKarl Rupp if ((xe==M) && (ys==0) && (zs==0)) n2 = size-m; 7058865f1eaSKarl Rupp if ((xe==M) && (ys==0) && (ze==P)) n20 = m*n-m; 7068865f1eaSKarl Rupp if ((xe==M) && (ye==N) && (zs==0)) n8 = size-m*n; 7078865f1eaSKarl Rupp if ((xe==M) && (ye==N) && (ze==P)) n26 = 0; 70847c6ae99SBarry Smith 70947c6ae99SBarry Smith /* Check for when not X,Y, and Z Periodic */ 71047c6ae99SBarry Smith 71147c6ae99SBarry Smith /* If not X periodic */ 7121321219cSEthan Coon if (bx != DMDA_BOUNDARY_PERIODIC) { 7138865f1eaSKarl Rupp if (xs==0) n0 = n3 = n6 = n9 = n12 = n15 = n18 = n21 = n24 = -2; 7148865f1eaSKarl Rupp if (xe==M) n2 = n5 = n8 = n11 = n14 = n17 = n20 = n23 = n26 = -2; 71547c6ae99SBarry Smith } 71647c6ae99SBarry Smith 71747c6ae99SBarry Smith /* If not Y periodic */ 7181321219cSEthan Coon if (by != DMDA_BOUNDARY_PERIODIC) { 7198865f1eaSKarl Rupp if (ys==0) n0 = n1 = n2 = n9 = n10 = n11 = n18 = n19 = n20 = -2; 7208865f1eaSKarl Rupp if (ye==N) n6 = n7 = n8 = n15 = n16 = n17 = n24 = n25 = n26 = -2; 72147c6ae99SBarry Smith } 72247c6ae99SBarry Smith 72347c6ae99SBarry Smith /* If not Z periodic */ 7241321219cSEthan Coon if (bz != DMDA_BOUNDARY_PERIODIC) { 7258865f1eaSKarl Rupp if (zs==0) n0 = n1 = n2 = n3 = n4 = n5 = n6 = n7 = n8 = -2; 7268865f1eaSKarl Rupp if (ze==P) n18 = n19 = n20 = n21 = n22 = n23 = n24 = n25 = n26 = -2; 72747c6ae99SBarry Smith } 72847c6ae99SBarry Smith 72947c6ae99SBarry Smith ierr = PetscMalloc(27*sizeof(PetscInt),&dd->neighbors);CHKERRQ(ierr); 7308865f1eaSKarl Rupp 73147c6ae99SBarry Smith dd->neighbors[0] = n0; 73247c6ae99SBarry Smith dd->neighbors[1] = n1; 73347c6ae99SBarry Smith dd->neighbors[2] = n2; 73447c6ae99SBarry Smith dd->neighbors[3] = n3; 73547c6ae99SBarry Smith dd->neighbors[4] = n4; 73647c6ae99SBarry Smith dd->neighbors[5] = n5; 73747c6ae99SBarry Smith dd->neighbors[6] = n6; 73847c6ae99SBarry Smith dd->neighbors[7] = n7; 73947c6ae99SBarry Smith dd->neighbors[8] = n8; 74047c6ae99SBarry Smith dd->neighbors[9] = n9; 74147c6ae99SBarry Smith dd->neighbors[10] = n10; 74247c6ae99SBarry Smith dd->neighbors[11] = n11; 74347c6ae99SBarry Smith dd->neighbors[12] = n12; 74447c6ae99SBarry Smith dd->neighbors[13] = rank; 74547c6ae99SBarry Smith dd->neighbors[14] = n14; 74647c6ae99SBarry Smith dd->neighbors[15] = n15; 74747c6ae99SBarry Smith dd->neighbors[16] = n16; 74847c6ae99SBarry Smith dd->neighbors[17] = n17; 74947c6ae99SBarry Smith dd->neighbors[18] = n18; 75047c6ae99SBarry Smith dd->neighbors[19] = n19; 75147c6ae99SBarry Smith dd->neighbors[20] = n20; 75247c6ae99SBarry Smith dd->neighbors[21] = n21; 75347c6ae99SBarry Smith dd->neighbors[22] = n22; 75447c6ae99SBarry Smith dd->neighbors[23] = n23; 75547c6ae99SBarry Smith dd->neighbors[24] = n24; 75647c6ae99SBarry Smith dd->neighbors[25] = n25; 75747c6ae99SBarry Smith dd->neighbors[26] = n26; 75847c6ae99SBarry Smith 75947c6ae99SBarry Smith /* If star stencil then delete the corner neighbors */ 760d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 76147c6ae99SBarry Smith /* save information about corner neighbors */ 76247c6ae99SBarry Smith sn0 = n0; sn1 = n1; sn2 = n2; sn3 = n3; sn5 = n5; sn6 = n6; sn7 = n7; 76347c6ae99SBarry Smith sn8 = n8; sn9 = n9; sn11 = n11; sn15 = n15; sn17 = n17; sn18 = n18; 76447c6ae99SBarry Smith sn19 = n19; sn20 = n20; sn21 = n21; sn23 = n23; sn24 = n24; sn25 = n25; 76547c6ae99SBarry Smith sn26 = n26; 7668865f1eaSKarl Rupp n0 = n1 = n2 = n3 = n5 = n6 = n7 = n8 = n9 = n11 = n15 = n17 = n18 = n19 = n20 = n21 = n23 = n24 = n25 = n26 = -1; 76747c6ae99SBarry Smith } 76847c6ae99SBarry Smith 76947c6ae99SBarry Smith ierr = PetscMalloc((Xe-Xs)*(Ye-Ys)*(Ze-Zs)*sizeof(PetscInt),&idx);CHKERRQ(ierr); 770*3bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)da,(Xe-Xs)*(Ye-Ys)*(Ze-Zs)*sizeof(PetscInt));CHKERRQ(ierr); 77147c6ae99SBarry Smith 77247c6ae99SBarry Smith nn = 0; 77347c6ae99SBarry Smith /* Bottom Level */ 77447c6ae99SBarry Smith for (k=0; k<s_z; k++) { 77547c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 77647c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 777ce00eea3SSatish Balay x_t = lx[n0 % m]; 77847c6ae99SBarry Smith y_t = ly[(n0 % (m*n))/m]; 77947c6ae99SBarry Smith z_t = lz[n0 / (m*n)]; 78047c6ae99SBarry 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; 7818865f1eaSKarl 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 */ 7828865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 78347c6ae99SBarry Smith } 78447c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 78547c6ae99SBarry Smith x_t = x; 78647c6ae99SBarry Smith y_t = ly[(n1 % (m*n))/m]; 78747c6ae99SBarry Smith z_t = lz[n1 / (m*n)]; 78847c6ae99SBarry 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; 7898865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t; /* 2D case */ 7908865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 79147c6ae99SBarry Smith } 79247c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 793ce00eea3SSatish Balay x_t = lx[n2 % m]; 79447c6ae99SBarry Smith y_t = ly[(n2 % (m*n))/m]; 79547c6ae99SBarry Smith z_t = lz[n2 / (m*n)]; 79647c6ae99SBarry 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; 7978865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t; /* 2D case */ 7988865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 79947c6ae99SBarry Smith } 80047c6ae99SBarry Smith } 80147c6ae99SBarry Smith 80247c6ae99SBarry Smith for (i=0; i<y; i++) { 80347c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 804ce00eea3SSatish Balay x_t = lx[n3 % m]; 80547c6ae99SBarry Smith y_t = y; 80647c6ae99SBarry Smith z_t = lz[n3 / (m*n)]; 80747c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8088865f1eaSKarl 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 */ 8098865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 81047c6ae99SBarry Smith } 81147c6ae99SBarry Smith 81247c6ae99SBarry Smith if (n4 >= 0) { /* middle */ 81347c6ae99SBarry Smith x_t = x; 81447c6ae99SBarry Smith y_t = y; 81547c6ae99SBarry Smith z_t = lz[n4 / (m*n)]; 81647c6ae99SBarry Smith s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8178865f1eaSKarl 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 */ 8188865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 819c2859e5eSBarry Smith } else if (bz == DMDA_BOUNDARY_MIRROR) { 8208865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 82147c6ae99SBarry Smith } 82247c6ae99SBarry Smith 82347c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 824ce00eea3SSatish Balay x_t = lx[n5 % m]; 82547c6ae99SBarry Smith y_t = y; 82647c6ae99SBarry Smith z_t = lz[n5 / (m*n)]; 82747c6ae99SBarry Smith s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8288865f1eaSKarl 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 */ 8298865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 83047c6ae99SBarry Smith } 83147c6ae99SBarry Smith } 83247c6ae99SBarry Smith 83347c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 83447c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 835ce00eea3SSatish Balay x_t = lx[n6 % m]; 83647c6ae99SBarry Smith y_t = ly[(n6 % (m*n))/m]; 83747c6ae99SBarry Smith z_t = lz[n6 / (m*n)]; 83847c6ae99SBarry Smith s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8398865f1eaSKarl 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 */ 8408865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 84147c6ae99SBarry Smith } 84247c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 84347c6ae99SBarry Smith x_t = x; 84447c6ae99SBarry Smith y_t = ly[(n7 % (m*n))/m]; 84547c6ae99SBarry Smith z_t = lz[n7 / (m*n)]; 84647c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8478865f1eaSKarl 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 */ 8488865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 84947c6ae99SBarry Smith } 85047c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 851ce00eea3SSatish Balay x_t = lx[n8 % m]; 85247c6ae99SBarry Smith y_t = ly[(n8 % (m*n))/m]; 85347c6ae99SBarry Smith z_t = lz[n8 / (m*n)]; 85447c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8558865f1eaSKarl 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 */ 8568865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 85747c6ae99SBarry Smith } 85847c6ae99SBarry Smith } 85947c6ae99SBarry Smith } 86047c6ae99SBarry Smith 86147c6ae99SBarry Smith /* Middle Level */ 86247c6ae99SBarry Smith for (k=0; k<z; k++) { 86347c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 86447c6ae99SBarry Smith if (n9 >= 0) { /* left below */ 865ce00eea3SSatish Balay x_t = lx[n9 % m]; 86647c6ae99SBarry Smith y_t = ly[(n9 % (m*n))/m]; 86747c6ae99SBarry Smith /* z_t = z; */ 86847c6ae99SBarry Smith s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 8698865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 87047c6ae99SBarry Smith } 87147c6ae99SBarry Smith if (n10 >= 0) { /* directly below */ 87247c6ae99SBarry Smith x_t = x; 87347c6ae99SBarry Smith y_t = ly[(n10 % (m*n))/m]; 87447c6ae99SBarry Smith /* z_t = z; */ 87547c6ae99SBarry Smith s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 8768865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 877c2859e5eSBarry Smith } else if (by == DMDA_BOUNDARY_MIRROR) { 8788865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 87947c6ae99SBarry Smith } 88047c6ae99SBarry Smith if (n11 >= 0) { /* right below */ 881ce00eea3SSatish Balay x_t = lx[n11 % m]; 88247c6ae99SBarry Smith y_t = ly[(n11 % (m*n))/m]; 88347c6ae99SBarry Smith /* z_t = z; */ 88447c6ae99SBarry Smith s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 8858865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 88647c6ae99SBarry Smith } 88747c6ae99SBarry Smith } 88847c6ae99SBarry Smith 88947c6ae99SBarry Smith for (i=0; i<y; i++) { 89047c6ae99SBarry Smith if (n12 >= 0) { /* directly left */ 891ce00eea3SSatish Balay x_t = lx[n12 % m]; 89247c6ae99SBarry Smith y_t = y; 89347c6ae99SBarry Smith /* z_t = z; */ 89447c6ae99SBarry Smith s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t; 8958865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 896c2859e5eSBarry Smith } else if (bx == DMDA_BOUNDARY_MIRROR) { 8978865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = 0; 89847c6ae99SBarry Smith } 89947c6ae99SBarry Smith 90047c6ae99SBarry Smith /* Interior */ 90147c6ae99SBarry Smith s_t = bases[rank] + i*x + k*x*y; 9028865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = s_t++; 90347c6ae99SBarry Smith 90447c6ae99SBarry Smith if (n14 >= 0) { /* directly right */ 905ce00eea3SSatish Balay x_t = lx[n14 % m]; 90647c6ae99SBarry Smith y_t = y; 90747c6ae99SBarry Smith /* z_t = z; */ 90847c6ae99SBarry Smith s_t = bases[n14] + i*x_t + k*x_t*y_t; 9098865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 910c2859e5eSBarry Smith } else if (bx == DMDA_BOUNDARY_MIRROR) { 9118865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = 0; 91247c6ae99SBarry Smith } 91347c6ae99SBarry Smith } 91447c6ae99SBarry Smith 91547c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 91647c6ae99SBarry Smith if (n15 >= 0) { /* left above */ 917ce00eea3SSatish Balay x_t = lx[n15 % m]; 91847c6ae99SBarry Smith y_t = ly[(n15 % (m*n))/m]; 91947c6ae99SBarry Smith /* z_t = z; */ 92047c6ae99SBarry Smith s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t; 9218865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 92247c6ae99SBarry Smith } 92347c6ae99SBarry Smith if (n16 >= 0) { /* directly above */ 92447c6ae99SBarry Smith x_t = x; 92547c6ae99SBarry Smith y_t = ly[(n16 % (m*n))/m]; 92647c6ae99SBarry Smith /* z_t = z; */ 92747c6ae99SBarry Smith s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t; 9288865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 929c2859e5eSBarry Smith } else if (by == DMDA_BOUNDARY_MIRROR) { 9308865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 93147c6ae99SBarry Smith } 93247c6ae99SBarry Smith if (n17 >= 0) { /* right above */ 933ce00eea3SSatish Balay x_t = lx[n17 % m]; 93447c6ae99SBarry Smith y_t = ly[(n17 % (m*n))/m]; 93547c6ae99SBarry Smith /* z_t = z; */ 93647c6ae99SBarry Smith s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t; 9378865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 93847c6ae99SBarry Smith } 93947c6ae99SBarry Smith } 94047c6ae99SBarry Smith } 94147c6ae99SBarry Smith 94247c6ae99SBarry Smith /* Upper Level */ 94347c6ae99SBarry Smith for (k=0; k<s_z; k++) { 94447c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 94547c6ae99SBarry Smith if (n18 >= 0) { /* left below */ 946ce00eea3SSatish Balay x_t = lx[n18 % m]; 94747c6ae99SBarry Smith y_t = ly[(n18 % (m*n))/m]; 94847c6ae99SBarry Smith /* z_t = lz[n18 / (m*n)]; */ 94947c6ae99SBarry Smith s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 9508865f1eaSKarl 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 */ 9518865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 95247c6ae99SBarry Smith } 95347c6ae99SBarry Smith if (n19 >= 0) { /* directly below */ 95447c6ae99SBarry Smith x_t = x; 95547c6ae99SBarry Smith y_t = ly[(n19 % (m*n))/m]; 95647c6ae99SBarry Smith /* z_t = lz[n19 / (m*n)]; */ 95747c6ae99SBarry Smith s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 9588865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n19] - (s_y+1-i)*x_t + x_t*y_t; /* 2d case */ 9598865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 96047c6ae99SBarry Smith } 96147c6ae99SBarry Smith if (n20 >= 0) { /* right below */ 962ce00eea3SSatish Balay x_t = lx[n20 % m]; 96347c6ae99SBarry Smith y_t = ly[(n20 % (m*n))/m]; 96447c6ae99SBarry Smith /* z_t = lz[n20 / (m*n)]; */ 96547c6ae99SBarry Smith s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 9668865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n20] - (s_y+1-i)*x_t + x_t*y_t; /* 2d case */ 9678865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 96847c6ae99SBarry Smith } 96947c6ae99SBarry Smith } 97047c6ae99SBarry Smith 97147c6ae99SBarry Smith for (i=0; i<y; i++) { 97247c6ae99SBarry Smith if (n21 >= 0) { /* directly left */ 973ce00eea3SSatish Balay x_t = lx[n21 % m]; 97447c6ae99SBarry Smith y_t = y; 97547c6ae99SBarry Smith /* z_t = lz[n21 / (m*n)]; */ 97647c6ae99SBarry Smith s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t; 9778865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n21] + (i+1)*x_t - s_x; /* 2d case */ 9788865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 97947c6ae99SBarry Smith } 98047c6ae99SBarry Smith 98147c6ae99SBarry Smith if (n22 >= 0) { /* middle */ 98247c6ae99SBarry Smith x_t = x; 98347c6ae99SBarry Smith y_t = y; 98447c6ae99SBarry Smith /* z_t = lz[n22 / (m*n)]; */ 98547c6ae99SBarry Smith s_t = bases[n22] + i*x_t + k*x_t*y_t; 9868865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n22] + i*x_t; /* 2d case */ 9878865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 988c2859e5eSBarry Smith } else if (bz == DMDA_BOUNDARY_MIRROR) { 9898865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 99047c6ae99SBarry Smith } 99147c6ae99SBarry Smith 99247c6ae99SBarry Smith if (n23 >= 0) { /* directly right */ 993ce00eea3SSatish Balay x_t = lx[n23 % m]; 99447c6ae99SBarry Smith y_t = y; 99547c6ae99SBarry Smith /* z_t = lz[n23 / (m*n)]; */ 99647c6ae99SBarry Smith s_t = bases[n23] + i*x_t + k*x_t*y_t; 9978865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n23] + i*x_t; /* 2d case */ 9988865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 99947c6ae99SBarry Smith } 100047c6ae99SBarry Smith } 100147c6ae99SBarry Smith 100247c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 100347c6ae99SBarry Smith if (n24 >= 0) { /* left above */ 1004ce00eea3SSatish Balay x_t = lx[n24 % m]; 100547c6ae99SBarry Smith y_t = ly[(n24 % (m*n))/m]; 100647c6ae99SBarry Smith /* z_t = lz[n24 / (m*n)]; */ 100747c6ae99SBarry Smith s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t; 10088865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n24] + i*x_t - s_x; /* 2d case */ 10098865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 101047c6ae99SBarry Smith } 101147c6ae99SBarry Smith if (n25 >= 0) { /* directly above */ 101247c6ae99SBarry Smith x_t = x; 101347c6ae99SBarry Smith y_t = ly[(n25 % (m*n))/m]; 101447c6ae99SBarry Smith /* z_t = lz[n25 / (m*n)]; */ 101547c6ae99SBarry Smith s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t; 10168865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n25] + (i-1)*x_t; /* 2d case */ 10178865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 101847c6ae99SBarry Smith } 101947c6ae99SBarry Smith if (n26 >= 0) { /* right above */ 1020ce00eea3SSatish Balay x_t = lx[n26 % m]; 102147c6ae99SBarry Smith y_t = ly[(n26 % (m*n))/m]; 102247c6ae99SBarry Smith /* z_t = lz[n26 / (m*n)]; */ 102347c6ae99SBarry Smith s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t; 10248865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n26] + (i-1)*x_t; /* 2d case */ 10258865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 102647c6ae99SBarry Smith } 102747c6ae99SBarry Smith } 102847c6ae99SBarry Smith } 1029ce00eea3SSatish Balay 1030ce00eea3SSatish Balay ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_COPY_VALUES,&from);CHKERRQ(ierr); 103147c6ae99SBarry Smith ierr = VecScatterCreate(global,from,local,to,>ol);CHKERRQ(ierr); 1032*3bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)gtol);CHKERRQ(ierr); 1033fcfd50ebSBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 1034fcfd50ebSBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 103547c6ae99SBarry Smith 1036d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 103747c6ae99SBarry Smith n0 = sn0; n1 = sn1; n2 = sn2; n3 = sn3; n5 = sn5; n6 = sn6; n7 = sn7; 103847c6ae99SBarry Smith n8 = sn8; n9 = sn9; n11 = sn11; n15 = sn15; n17 = sn17; n18 = sn18; 103947c6ae99SBarry Smith n19 = sn19; n20 = sn20; n21 = sn21; n23 = sn23; n24 = sn24; n25 = sn25; 104047c6ae99SBarry Smith n26 = sn26; 1041ce00eea3SSatish Balay } 104247c6ae99SBarry Smith 104388661749SPeter Brune if (((stencil_type == DMDA_STENCIL_STAR) || 10441321219cSEthan Coon (bx != DMDA_BOUNDARY_PERIODIC && bx) || 10451321219cSEthan Coon (by != DMDA_BOUNDARY_PERIODIC && by) || 1046d9c9ebe5SPeter Brune (bz != DMDA_BOUNDARY_PERIODIC && bz))) { 1047ce00eea3SSatish Balay /* 1048ce00eea3SSatish Balay Recompute the local to global mappings, this time keeping the 1049ce00eea3SSatish Balay information about the cross corner processor numbers. 1050ce00eea3SSatish Balay */ 105147c6ae99SBarry Smith nn = 0; 105247c6ae99SBarry Smith /* Bottom Level */ 105347c6ae99SBarry Smith for (k=0; k<s_z; k++) { 105447c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 105547c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 1056ce00eea3SSatish Balay x_t = lx[n0 % m]; 105747c6ae99SBarry Smith y_t = ly[(n0 % (m*n))/m]; 105847c6ae99SBarry Smith z_t = lz[n0 / (m*n)]; 105947c6ae99SBarry 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; 10608865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1061ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0 && Zs-zs < 0) { 10628865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 106347c6ae99SBarry Smith } 106447c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 106547c6ae99SBarry Smith x_t = x; 106647c6ae99SBarry Smith y_t = ly[(n1 % (m*n))/m]; 106747c6ae99SBarry Smith z_t = lz[n1 / (m*n)]; 106847c6ae99SBarry 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; 10698865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1070ce00eea3SSatish Balay } else if (Ys-ys < 0 && Zs-zs < 0) { 10718865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 107247c6ae99SBarry Smith } 107347c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 1074ce00eea3SSatish Balay x_t = lx[n2 % m]; 107547c6ae99SBarry Smith y_t = ly[(n2 % (m*n))/m]; 107647c6ae99SBarry Smith z_t = lz[n2 / (m*n)]; 107747c6ae99SBarry 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; 10788865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1079ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0 && Zs-zs < 0) { 10808865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 108147c6ae99SBarry Smith } 108247c6ae99SBarry Smith } 108347c6ae99SBarry Smith 108447c6ae99SBarry Smith for (i=0; i<y; i++) { 108547c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 1086ce00eea3SSatish Balay x_t = lx[n3 % m]; 108747c6ae99SBarry Smith y_t = y; 108847c6ae99SBarry Smith z_t = lz[n3 / (m*n)]; 108947c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 10908865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1091ce00eea3SSatish Balay } else if (Xs-xs < 0 && Zs-zs < 0) { 10928865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 109347c6ae99SBarry Smith } 109447c6ae99SBarry Smith 109547c6ae99SBarry Smith if (n4 >= 0) { /* middle */ 109647c6ae99SBarry Smith x_t = x; 109747c6ae99SBarry Smith y_t = y; 109847c6ae99SBarry Smith z_t = lz[n4 / (m*n)]; 109947c6ae99SBarry Smith s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11008865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1101ce00eea3SSatish Balay } else if (Zs-zs < 0) { 1102c2859e5eSBarry Smith if (bz == DMDA_BOUNDARY_MIRROR) { 11038865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 1104c2859e5eSBarry Smith } else { 11058865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 110647c6ae99SBarry Smith } 1107c2859e5eSBarry Smith } 110847c6ae99SBarry Smith 110947c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 1110ce00eea3SSatish Balay x_t = lx[n5 % m]; 111147c6ae99SBarry Smith y_t = y; 111247c6ae99SBarry Smith z_t = lz[n5 / (m*n)]; 111347c6ae99SBarry Smith s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11148865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1115ce00eea3SSatish Balay } else if (xe-Xe < 0 && Zs-zs < 0) { 11168865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 111747c6ae99SBarry Smith } 111847c6ae99SBarry Smith } 111947c6ae99SBarry Smith 112047c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 112147c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 1122ce00eea3SSatish Balay x_t = lx[n6 % m]; 112347c6ae99SBarry Smith y_t = ly[(n6 % (m*n))/m]; 112447c6ae99SBarry Smith z_t = lz[n6 / (m*n)]; 112547c6ae99SBarry Smith s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11268865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1127ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0 && Zs-zs < 0) { 11288865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 112947c6ae99SBarry Smith } 113047c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 113147c6ae99SBarry Smith x_t = x; 113247c6ae99SBarry Smith y_t = ly[(n7 % (m*n))/m]; 113347c6ae99SBarry Smith z_t = lz[n7 / (m*n)]; 113447c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11358865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1136ce00eea3SSatish Balay } else if (ye-Ye < 0 && Zs-zs < 0) { 11378865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 113847c6ae99SBarry Smith } 113947c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 1140ce00eea3SSatish Balay x_t = lx[n8 % m]; 114147c6ae99SBarry Smith y_t = ly[(n8 % (m*n))/m]; 114247c6ae99SBarry Smith z_t = lz[n8 / (m*n)]; 114347c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11448865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1145ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0 && Zs-zs < 0) { 11468865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 114747c6ae99SBarry Smith } 114847c6ae99SBarry Smith } 114947c6ae99SBarry Smith } 115047c6ae99SBarry Smith 115147c6ae99SBarry Smith /* Middle Level */ 115247c6ae99SBarry Smith for (k=0; k<z; k++) { 115347c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 115447c6ae99SBarry Smith if (n9 >= 0) { /* left below */ 1155ce00eea3SSatish Balay x_t = lx[n9 % m]; 115647c6ae99SBarry Smith y_t = ly[(n9 % (m*n))/m]; 115747c6ae99SBarry Smith /* z_t = z; */ 115847c6ae99SBarry Smith s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 11598865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1160ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0) { 11618865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 116247c6ae99SBarry Smith } 116347c6ae99SBarry Smith if (n10 >= 0) { /* directly below */ 116447c6ae99SBarry Smith x_t = x; 116547c6ae99SBarry Smith y_t = ly[(n10 % (m*n))/m]; 116647c6ae99SBarry Smith /* z_t = z; */ 116747c6ae99SBarry Smith s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 11688865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1169ce00eea3SSatish Balay } else if (Ys-ys < 0) { 1170c2859e5eSBarry Smith if (by == DMDA_BOUNDARY_MIRROR) { 11718865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 1172c2859e5eSBarry Smith } else { 11738865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 1174c2859e5eSBarry Smith } 117547c6ae99SBarry Smith } 117647c6ae99SBarry Smith if (n11 >= 0) { /* right below */ 1177ce00eea3SSatish Balay x_t = lx[n11 % m]; 117847c6ae99SBarry Smith y_t = ly[(n11 % (m*n))/m]; 117947c6ae99SBarry Smith /* z_t = z; */ 118047c6ae99SBarry Smith s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 11818865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1182ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0) { 11838865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 118447c6ae99SBarry Smith } 118547c6ae99SBarry Smith } 118647c6ae99SBarry Smith 118747c6ae99SBarry Smith for (i=0; i<y; i++) { 118847c6ae99SBarry Smith if (n12 >= 0) { /* directly left */ 1189ce00eea3SSatish Balay x_t = lx[n12 % m]; 119047c6ae99SBarry Smith y_t = y; 119147c6ae99SBarry Smith /* z_t = z; */ 119247c6ae99SBarry Smith s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t; 11938865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1194ce00eea3SSatish Balay } else if (Xs-xs < 0) { 1195c2859e5eSBarry Smith if (bx == DMDA_BOUNDARY_MIRROR) { 11968865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = 0; 1197c2859e5eSBarry Smith } else { 11988865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 119947c6ae99SBarry Smith } 1200c2859e5eSBarry Smith } 120147c6ae99SBarry Smith 120247c6ae99SBarry Smith /* Interior */ 120347c6ae99SBarry Smith s_t = bases[rank] + i*x + k*x*y; 12048865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = s_t++; 120547c6ae99SBarry Smith 120647c6ae99SBarry Smith if (n14 >= 0) { /* directly right */ 1207ce00eea3SSatish Balay x_t = lx[n14 % m]; 120847c6ae99SBarry Smith y_t = y; 120947c6ae99SBarry Smith /* z_t = z; */ 121047c6ae99SBarry Smith s_t = bases[n14] + i*x_t + k*x_t*y_t; 12118865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1212ce00eea3SSatish Balay } else if (xe-Xe < 0) { 1213c2859e5eSBarry Smith if (bx == DMDA_BOUNDARY_MIRROR) { 12148865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = 0; 1215c2859e5eSBarry Smith } else { 12168865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 121747c6ae99SBarry Smith } 121847c6ae99SBarry Smith } 1219c2859e5eSBarry Smith } 122047c6ae99SBarry Smith 122147c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 122247c6ae99SBarry Smith if (n15 >= 0) { /* left above */ 1223ce00eea3SSatish Balay x_t = lx[n15 % m]; 122447c6ae99SBarry Smith y_t = ly[(n15 % (m*n))/m]; 122547c6ae99SBarry Smith /* z_t = z; */ 122647c6ae99SBarry Smith s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t; 12278865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1228ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0) { 12298865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 123047c6ae99SBarry Smith } 123147c6ae99SBarry Smith if (n16 >= 0) { /* directly above */ 123247c6ae99SBarry Smith x_t = x; 123347c6ae99SBarry Smith y_t = ly[(n16 % (m*n))/m]; 123447c6ae99SBarry Smith /* z_t = z; */ 123547c6ae99SBarry Smith s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t; 12368865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1237ce00eea3SSatish Balay } else if (ye-Ye < 0) { 1238c2859e5eSBarry Smith if (by == DMDA_BOUNDARY_MIRROR) { 12398865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 1240c2859e5eSBarry Smith } else { 12418865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 124247c6ae99SBarry Smith } 1243c2859e5eSBarry Smith } 124447c6ae99SBarry Smith if (n17 >= 0) { /* right above */ 1245ce00eea3SSatish Balay x_t = lx[n17 % m]; 124647c6ae99SBarry Smith y_t = ly[(n17 % (m*n))/m]; 124747c6ae99SBarry Smith /* z_t = z; */ 124847c6ae99SBarry Smith s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t; 12498865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1250ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0) { 12518865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 125247c6ae99SBarry Smith } 125347c6ae99SBarry Smith } 125447c6ae99SBarry Smith } 125547c6ae99SBarry Smith 125647c6ae99SBarry Smith /* Upper Level */ 125747c6ae99SBarry Smith for (k=0; k<s_z; k++) { 125847c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 125947c6ae99SBarry Smith if (n18 >= 0) { /* left below */ 1260ce00eea3SSatish Balay x_t = lx[n18 % m]; 126147c6ae99SBarry Smith y_t = ly[(n18 % (m*n))/m]; 126247c6ae99SBarry Smith /* z_t = lz[n18 / (m*n)]; */ 126347c6ae99SBarry Smith s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 12648865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1265ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0 && ze-Ze < 0) { 12668865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 126747c6ae99SBarry Smith } 126847c6ae99SBarry Smith if (n19 >= 0) { /* directly below */ 126947c6ae99SBarry Smith x_t = x; 127047c6ae99SBarry Smith y_t = ly[(n19 % (m*n))/m]; 127147c6ae99SBarry Smith /* z_t = lz[n19 / (m*n)]; */ 127247c6ae99SBarry Smith s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 12738865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1274ce00eea3SSatish Balay } else if (Ys-ys < 0 && ze-Ze < 0) { 12758865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 127647c6ae99SBarry Smith } 127747c6ae99SBarry Smith if (n20 >= 0) { /* right below */ 1278ce00eea3SSatish Balay x_t = lx[n20 % m]; 127947c6ae99SBarry Smith y_t = ly[(n20 % (m*n))/m]; 128047c6ae99SBarry Smith /* z_t = lz[n20 / (m*n)]; */ 128147c6ae99SBarry Smith s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 12828865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1283ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0 && ze-Ze < 0) { 12848865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 128547c6ae99SBarry Smith } 128647c6ae99SBarry Smith } 128747c6ae99SBarry Smith 128847c6ae99SBarry Smith for (i=0; i<y; i++) { 128947c6ae99SBarry Smith if (n21 >= 0) { /* directly left */ 1290ce00eea3SSatish Balay x_t = lx[n21 % m]; 129147c6ae99SBarry Smith y_t = y; 129247c6ae99SBarry Smith /* z_t = lz[n21 / (m*n)]; */ 129347c6ae99SBarry Smith s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t; 12948865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1295ce00eea3SSatish Balay } else if (Xs-xs < 0 && ze-Ze < 0) { 12968865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 129747c6ae99SBarry Smith } 129847c6ae99SBarry Smith 129947c6ae99SBarry Smith if (n22 >= 0) { /* middle */ 130047c6ae99SBarry Smith x_t = x; 130147c6ae99SBarry Smith y_t = y; 130247c6ae99SBarry Smith /* z_t = lz[n22 / (m*n)]; */ 130347c6ae99SBarry Smith s_t = bases[n22] + i*x_t + k*x_t*y_t; 13048865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1305ce00eea3SSatish Balay } else if (ze-Ze < 0) { 1306c2859e5eSBarry Smith if (bz == DMDA_BOUNDARY_MIRROR) { 13078865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 1308c2859e5eSBarry Smith } else { 13098865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 131047c6ae99SBarry Smith } 1311c2859e5eSBarry Smith } 131247c6ae99SBarry Smith 131347c6ae99SBarry Smith if (n23 >= 0) { /* directly right */ 1314ce00eea3SSatish Balay x_t = lx[n23 % m]; 131547c6ae99SBarry Smith y_t = y; 131647c6ae99SBarry Smith /* z_t = lz[n23 / (m*n)]; */ 131747c6ae99SBarry Smith s_t = bases[n23] + i*x_t + k*x_t*y_t; 13188865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1319ce00eea3SSatish Balay } else if (xe-Xe < 0 && ze-Ze < 0) { 13208865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 132147c6ae99SBarry Smith } 132247c6ae99SBarry Smith } 132347c6ae99SBarry Smith 132447c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 132547c6ae99SBarry Smith if (n24 >= 0) { /* left above */ 1326ce00eea3SSatish Balay x_t = lx[n24 % m]; 132747c6ae99SBarry Smith y_t = ly[(n24 % (m*n))/m]; 132847c6ae99SBarry Smith /* z_t = lz[n24 / (m*n)]; */ 132947c6ae99SBarry Smith s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t; 13308865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1331ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0 && ze-Ze < 0) { 13328865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 133347c6ae99SBarry Smith } 133447c6ae99SBarry Smith if (n25 >= 0) { /* directly above */ 133547c6ae99SBarry Smith x_t = x; 133647c6ae99SBarry Smith y_t = ly[(n25 % (m*n))/m]; 133747c6ae99SBarry Smith /* z_t = lz[n25 / (m*n)]; */ 133847c6ae99SBarry Smith s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t; 13398865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1340ce00eea3SSatish Balay } else if (ye-Ye < 0 && ze-Ze < 0) { 13418865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 134247c6ae99SBarry Smith } 134347c6ae99SBarry Smith if (n26 >= 0) { /* right above */ 1344ce00eea3SSatish Balay x_t = lx[n26 % m]; 134547c6ae99SBarry Smith y_t = ly[(n26 % (m*n))/m]; 134647c6ae99SBarry Smith /* z_t = lz[n26 / (m*n)]; */ 134747c6ae99SBarry Smith s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t; 13488865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1349ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0 && ze-Ze < 0) { 13508865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 135147c6ae99SBarry Smith } 135247c6ae99SBarry Smith } 135347c6ae99SBarry Smith } 135447c6ae99SBarry Smith } 135547c6ae99SBarry Smith /* 135647c6ae99SBarry Smith Set the local to global ordering in the global vector, this allows use 135747c6ae99SBarry Smith of VecSetValuesLocal(). 135847c6ae99SBarry Smith */ 1359ce00eea3SSatish Balay ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_OWN_POINTER,<ogis);CHKERRQ(ierr); 1360ce00eea3SSatish Balay ierr = PetscMalloc(nn*dof*sizeof(PetscInt),&idx_cpy);CHKERRQ(ierr); 1361*3bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)da,nn*dof*sizeof(PetscInt));CHKERRQ(ierr); 13623bf036e2SBarry Smith ierr = ISGetIndices(ltogis, &idx_full);CHKERRQ(ierr); 1363ce00eea3SSatish Balay ierr = PetscMemcpy(idx_cpy,idx_full,nn*dof*sizeof(PetscInt));CHKERRQ(ierr); 13643bf036e2SBarry Smith ierr = ISRestoreIndices(ltogis, &idx_full);CHKERRQ(ierr); 1365ce00eea3SSatish Balay ierr = ISLocalToGlobalMappingCreateIS(ltogis,&da->ltogmap);CHKERRQ(ierr); 1366*3bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap);CHKERRQ(ierr); 1367fcfd50ebSBarry Smith ierr = ISDestroy(<ogis);CHKERRQ(ierr); 13681411c6eeSJed Brown ierr = ISLocalToGlobalMappingBlock(da->ltogmap,dd->w,&da->ltogmapb);CHKERRQ(ierr); 1369*3bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap);CHKERRQ(ierr); 137047c6ae99SBarry Smith 1371ce00eea3SSatish Balay ierr = PetscFree2(bases,ldims);CHKERRQ(ierr); 1372ce00eea3SSatish Balay dd->m = m; dd->n = n; dd->p = p; 1373ce00eea3SSatish Balay /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */ 1374ce00eea3SSatish Balay dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = zs; dd->ze = ze; 1375ce00eea3SSatish Balay dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = Zs; dd->Ze = Ze; 1376ce00eea3SSatish Balay 1377fcfd50ebSBarry Smith ierr = VecDestroy(&local);CHKERRQ(ierr); 1378fcfd50ebSBarry Smith ierr = VecDestroy(&global);CHKERRQ(ierr); 1379ce00eea3SSatish Balay 1380ce00eea3SSatish Balay dd->gtol = gtol; 1381ce00eea3SSatish Balay dd->ltog = ltog; 1382ce00eea3SSatish Balay dd->idx = idx_cpy; 1383ce00eea3SSatish Balay dd->Nl = nn*dof; 1384ce00eea3SSatish Balay dd->base = base; 1385ce00eea3SSatish Balay da->ops->view = DMView_DA_3d; 13860298fd71SBarry Smith dd->ltol = NULL; 13870298fd71SBarry Smith dd->ao = NULL; 138847c6ae99SBarry Smith PetscFunctionReturn(0); 138947c6ae99SBarry Smith } 1390564755cdSBarry Smith 139147c6ae99SBarry Smith 139247c6ae99SBarry Smith #undef __FUNCT__ 1393aa219208SBarry Smith #define __FUNCT__ "DMDACreate3d" 139447c6ae99SBarry Smith /*@C 1395aa219208SBarry Smith DMDACreate3d - Creates an object that will manage the communication of three-dimensional 139647c6ae99SBarry Smith regular array data that is distributed across some processors. 139747c6ae99SBarry Smith 139847c6ae99SBarry Smith Collective on MPI_Comm 139947c6ae99SBarry Smith 140047c6ae99SBarry Smith Input Parameters: 140147c6ae99SBarry Smith + comm - MPI communicator 14021321219cSEthan Coon . bx,by,bz - type of ghost nodes the array have. 14031321219cSEthan Coon Use one of DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_PERIODIC. 1404aa219208SBarry Smith . stencil_type - Type of stencil (DMDA_STENCIL_STAR or DMDA_STENCIL_BOX) 140547c6ae99SBarry 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 140647c6ae99SBarry Smith from the command line with -da_grid_x <M> -da_grid_y <N> -da_grid_z <P>) 140747c6ae99SBarry Smith . m,n,p - corresponding number of processors in each dimension 140847c6ae99SBarry Smith (or PETSC_DECIDE to have calculated) 140947c6ae99SBarry Smith . dof - number of degrees of freedom per node 141010d7c030SMatthew G Knepley . s - stencil width 141110d7c030SMatthew G Knepley - lx, ly, lz - arrays containing the number of nodes in each cell along 14120298fd71SBarry Smith the x, y, and z coordinates, or NULL. If non-null, these 141347c6ae99SBarry Smith must be of length as m,n,p and the corresponding 141447c6ae99SBarry Smith m,n, or p cannot be PETSC_DECIDE. Sum of the lx[] entries must be M, sum of 141547c6ae99SBarry Smith the ly[] must N, sum of the lz[] must be P 141647c6ae99SBarry Smith 141747c6ae99SBarry Smith Output Parameter: 141847c6ae99SBarry Smith . da - the resulting distributed array object 141947c6ae99SBarry Smith 142047c6ae99SBarry Smith Options Database Key: 1421706a11cbSBarry Smith + -dm_view - Calls DMView() at the conclusion of DMDACreate3d() 142247c6ae99SBarry Smith . -da_grid_x <nx> - number of grid points in x direction, if M < 0 142347c6ae99SBarry Smith . -da_grid_y <ny> - number of grid points in y direction, if N < 0 142447c6ae99SBarry Smith . -da_grid_z <nz> - number of grid points in z direction, if P < 0 142547c6ae99SBarry Smith . -da_processors_x <MX> - number of processors in x direction 142647c6ae99SBarry Smith . -da_processors_y <MY> - number of processors in y direction 142747c6ae99SBarry Smith . -da_processors_z <MZ> - number of processors in z direction 1428e0f5d30fSBarry Smith . -da_refine_x <rx> - refinement ratio in x direction 1429e0f5d30fSBarry Smith . -da_refine_y <ry> - refinement ratio in y direction 1430e0f5d30fSBarry Smith . -da_refine_z <rz>- refinement ratio in z directio 1431e0f5d30fSBarry Smith - -da_refine <n> - refine the DMDA n times before creating it, , if M, N, or P < 0 143247c6ae99SBarry Smith 143347c6ae99SBarry Smith Level: beginner 143447c6ae99SBarry Smith 143547c6ae99SBarry Smith Notes: 1436aa219208SBarry Smith The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the 1437aa219208SBarry Smith standard 7-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes 143847c6ae99SBarry Smith the standard 27-pt stencil. 143947c6ae99SBarry Smith 1440aa219208SBarry Smith The array data itself is NOT stored in the DMDA, it is stored in Vec objects; 1441564755cdSBarry Smith The appropriate vector objects can be obtained with calls to DMCreateGlobalVector() 1442564755cdSBarry Smith and DMCreateLocalVector() and calls to VecDuplicate() if more are needed. 144347c6ae99SBarry Smith 144447c6ae99SBarry Smith .keywords: distributed array, create, three-dimensional 144547c6ae99SBarry Smith 1446aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate2d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(), 1447aa219208SBarry Smith DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMDALocalToLocalBegin(), DMDALocalToLocalEnd(), DMDASetRefinementFactor(), 1448d461ba97SBarry Smith DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges() 144947c6ae99SBarry Smith 145047c6ae99SBarry Smith @*/ 14511321219cSEthan Coon PetscErrorCode DMDACreate3d(MPI_Comm comm,DMDABoundaryType bx,DMDABoundaryType by,DMDABoundaryType bz,DMDAStencilType stencil_type,PetscInt M, 14529a42bb27SBarry 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) 145347c6ae99SBarry Smith { 145447c6ae99SBarry Smith PetscErrorCode ierr; 145547c6ae99SBarry Smith 145647c6ae99SBarry Smith PetscFunctionBegin; 1457aa219208SBarry Smith ierr = DMDACreate(comm, da);CHKERRQ(ierr); 1458aa219208SBarry Smith ierr = DMDASetDim(*da, 3);CHKERRQ(ierr); 1459aa219208SBarry Smith ierr = DMDASetSizes(*da, M, N, P);CHKERRQ(ierr); 1460aa219208SBarry Smith ierr = DMDASetNumProcs(*da, m, n, p);CHKERRQ(ierr); 14611321219cSEthan Coon ierr = DMDASetBoundaryType(*da, bx, by, bz);CHKERRQ(ierr); 1462aa219208SBarry Smith ierr = DMDASetDof(*da, dof);CHKERRQ(ierr); 1463aa219208SBarry Smith ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr); 1464aa219208SBarry Smith ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr); 1465aa219208SBarry Smith ierr = DMDASetOwnershipRanges(*da, lx, ly, lz);CHKERRQ(ierr); 146647c6ae99SBarry Smith /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */ 14679a42bb27SBarry Smith ierr = DMSetFromOptions(*da);CHKERRQ(ierr); 14689a42bb27SBarry Smith ierr = DMSetUp(*da);CHKERRQ(ierr); 1469146574abSBarry Smith ierr = DMViewFromOptions(*da,NULL,"-dm_view");CHKERRQ(ierr); 147047c6ae99SBarry Smith PetscFunctionReturn(0); 147147c6ae99SBarry Smith } 1472