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 7b45d2f2cSJed Brown #include <petsc-private/daimpl.h> /*I "petscdmda.h" I*/ 847c6ae99SBarry Smith 947c6ae99SBarry Smith #undef __FUNCT__ 109a42bb27SBarry Smith #define __FUNCT__ "DMView_DA_3d" 119a42bb27SBarry Smith PetscErrorCode DMView_DA_3d(DM da,PetscViewer viewer) 1247c6ae99SBarry Smith { 1347c6ae99SBarry Smith PetscErrorCode ierr; 1447c6ae99SBarry Smith PetscMPIInt rank; 159a42bb27SBarry Smith PetscBool iascii,isdraw,isbinary; 1647c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 179a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 189a42bb27SBarry Smith PetscBool ismatlab; 199a42bb27SBarry Smith #endif 2047c6ae99SBarry Smith 2147c6ae99SBarry Smith PetscFunctionBegin; 2247c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)da)->comm,&rank);CHKERRQ(ierr); 2347c6ae99SBarry Smith 24251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 25251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);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 337b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 3447c6ae99SBarry Smith ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 3547c6ae99SBarry Smith if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL) { 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; 4847c6ae99SBarry 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); 496636e97aSMatthew G Knepley ierr = VecRestoreArrayRead(da->coordinates,&coors);CHKERRQ(ierr); 5047c6ae99SBarry Smith } 5147c6ae99SBarry Smith #endif 5247c6ae99SBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 537b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 54616157d6SJed Brown } else { 55616157d6SJed Brown ierr = DMView_DA_VTK(da,viewer);CHKERRQ(ierr); 5647c6ae99SBarry Smith } 5747c6ae99SBarry Smith } else if (isdraw) { 5847c6ae99SBarry Smith PetscDraw draw; 5947c6ae99SBarry Smith PetscReal ymin = -1.0,ymax = (PetscReal)dd->N; 6047c6ae99SBarry Smith PetscReal xmin = -1.0,xmax = (PetscReal)((dd->M+2)*dd->P),x,y,ycoord,xcoord; 6147c6ae99SBarry Smith PetscInt k,plane,base,*idx; 6247c6ae99SBarry Smith char node[10]; 6347c6ae99SBarry Smith PetscBool isnull; 6447c6ae99SBarry Smith 6547c6ae99SBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 6647c6ae99SBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 6747c6ae99SBarry Smith ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr); 6847c6ae99SBarry Smith ierr = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr); 6947c6ae99SBarry Smith 7047c6ae99SBarry Smith /* first processor draw all node lines */ 7147c6ae99SBarry Smith if (!rank) { 7247c6ae99SBarry Smith for (k=0; k<dd->P; k++) { 7347c6ae99SBarry Smith ymin = 0.0; ymax = (PetscReal)(dd->N - 1); 7447c6ae99SBarry Smith for (xmin=(PetscReal)(k*(dd->M+1)); xmin<(PetscReal)(dd->M+(k*(dd->M+1))); xmin++) { 7547c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr); 7647c6ae99SBarry Smith } 7747c6ae99SBarry Smith 7847c6ae99SBarry Smith xmin = (PetscReal)(k*(dd->M+1)); xmax = xmin + (PetscReal)(dd->M - 1); 7947c6ae99SBarry Smith for (ymin=0; ymin<(PetscReal)dd->N; ymin++) { 8047c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr); 8147c6ae99SBarry Smith } 8247c6ae99SBarry Smith } 8347c6ae99SBarry Smith } 8447c6ae99SBarry Smith ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr); 8547c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 8647c6ae99SBarry Smith 8747c6ae99SBarry Smith for (k=0; k<dd->P; k++) { /*Go through and draw for each plane*/ 8847c6ae99SBarry Smith if ((k >= dd->zs) && (k < dd->ze)) { 8947c6ae99SBarry Smith /* draw my box */ 9047c6ae99SBarry Smith ymin = dd->ys; 9147c6ae99SBarry Smith ymax = dd->ye - 1; 9247c6ae99SBarry Smith xmin = dd->xs/dd->w + (dd->M+1)*k; 9347c6ae99SBarry Smith xmax =(dd->xe-1)/dd->w + (dd->M+1)*k; 9447c6ae99SBarry Smith 9547c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr); 9647c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 9747c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 9847c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 9947c6ae99SBarry Smith 10047c6ae99SBarry Smith xmin = dd->xs/dd->w; 10147c6ae99SBarry Smith xmax =(dd->xe-1)/dd->w; 10247c6ae99SBarry Smith 10347c6ae99SBarry Smith /* put in numbers*/ 10447c6ae99SBarry Smith base = (dd->base+(dd->xe-dd->xs)*(dd->ye-dd->ys)*(k-dd->zs))/dd->w; 10547c6ae99SBarry Smith 10647c6ae99SBarry Smith /* Identify which processor owns the box */ 10747c6ae99SBarry Smith sprintf(node,"%d",rank); 10847c6ae99SBarry Smith ierr = PetscDrawString(draw,xmin+(dd->M+1)*k+.2,ymin+.3,PETSC_DRAW_RED,node);CHKERRQ(ierr); 10947c6ae99SBarry Smith 11047c6ae99SBarry Smith for (y=ymin; y<=ymax; y++) { 11147c6ae99SBarry Smith for (x=xmin+(dd->M+1)*k; x<=xmax+(dd->M+1)*k; x++) { 11247c6ae99SBarry Smith sprintf(node,"%d",(int)base++); 11347c6ae99SBarry Smith ierr = PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);CHKERRQ(ierr); 11447c6ae99SBarry Smith } 11547c6ae99SBarry Smith } 11647c6ae99SBarry Smith 11747c6ae99SBarry Smith } 11847c6ae99SBarry Smith } 11947c6ae99SBarry Smith ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr); 12047c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 12147c6ae99SBarry Smith 12247c6ae99SBarry Smith for (k=0-dd->s; k<dd->P+dd->s; k++) { 12347c6ae99SBarry Smith /* Go through and draw for each plane */ 12447c6ae99SBarry Smith if ((k >= dd->Zs) && (k < dd->Ze)) { 12547c6ae99SBarry Smith 12647c6ae99SBarry Smith /* overlay ghost numbers, useful for error checking */ 12747c6ae99SBarry Smith base = (dd->Xe-dd->Xs)*(dd->Ye-dd->Ys)*(k-dd->Zs); idx = dd->idx; 12847c6ae99SBarry Smith plane=k; 12947c6ae99SBarry Smith /* Keep z wrap around points on the dradrawg */ 130*8865f1eaSKarl Rupp if (k<0) plane=dd->P+k; 131*8865f1eaSKarl Rupp if (k>=dd->P) plane=k-dd->P; 13247c6ae99SBarry Smith ymin = dd->Ys; ymax = dd->Ye; 13347c6ae99SBarry Smith xmin = (dd->M+1)*plane*dd->w; 13447c6ae99SBarry Smith xmax = (dd->M+1)*plane*dd->w+dd->M*dd->w; 13547c6ae99SBarry Smith for (y=ymin; y<ymax; y++) { 13647c6ae99SBarry Smith for (x=xmin+dd->Xs; x<xmin+dd->Xe; x+=dd->w) { 13747c6ae99SBarry Smith sprintf(node,"%d",(int)(idx[base]/dd->w)); 13847c6ae99SBarry Smith ycoord = y; 13947c6ae99SBarry Smith /*Keep y wrap around points on drawing */ 140*8865f1eaSKarl Rupp if (y<0) ycoord = dd->N+y; 14147c6ae99SBarry Smith 142*8865f1eaSKarl Rupp if (y>=dd->N) ycoord = y-dd->N; 14347c6ae99SBarry Smith xcoord = x; /* Keep x wrap points on drawing */ 14447c6ae99SBarry Smith 145*8865f1eaSKarl Rupp if (x<xmin) xcoord = xmax - (xmin-x); 146*8865f1eaSKarl Rupp if (x>=xmax) xcoord = xmin + (x-xmax); 14747c6ae99SBarry Smith ierr = PetscDrawString(draw,xcoord/dd->w,ycoord,PETSC_DRAW_BLUE,node);CHKERRQ(ierr); 14847c6ae99SBarry Smith base+=dd->w; 14947c6ae99SBarry Smith } 15047c6ae99SBarry Smith } 15147c6ae99SBarry Smith } 15247c6ae99SBarry Smith } 15347c6ae99SBarry Smith ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr); 15447c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 1559a42bb27SBarry Smith } else if (isbinary) { 1569a42bb27SBarry Smith ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr); 1579a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 1589a42bb27SBarry Smith } else if (ismatlab) { 1599a42bb27SBarry Smith ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr); 1609a42bb27SBarry Smith #endif 16111aeaf0aSBarry Smith } 16247c6ae99SBarry Smith PetscFunctionReturn(0); 16347c6ae99SBarry Smith } 16447c6ae99SBarry Smith 16547c6ae99SBarry Smith #undef __FUNCT__ 1669a42bb27SBarry Smith #define __FUNCT__ "DMSetUp_DA_3D" 1677087cfbeSBarry Smith PetscErrorCode DMSetUp_DA_3D(DM da) 16847c6ae99SBarry Smith { 16947c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 17047c6ae99SBarry Smith const PetscInt M = dd->M; 17147c6ae99SBarry Smith const PetscInt N = dd->N; 17247c6ae99SBarry Smith const PetscInt P = dd->P; 17347c6ae99SBarry Smith PetscInt m = dd->m; 17447c6ae99SBarry Smith PetscInt n = dd->n; 17547c6ae99SBarry Smith PetscInt p = dd->p; 17647c6ae99SBarry Smith const PetscInt dof = dd->w; 17747c6ae99SBarry Smith const PetscInt s = dd->s; 17819fd82e9SBarry Smith DMDABoundaryType bx = dd->bx; 17919fd82e9SBarry Smith DMDABoundaryType by = dd->by; 18019fd82e9SBarry Smith DMDABoundaryType bz = dd->bz; 18119fd82e9SBarry Smith DMDAStencilType stencil_type = dd->stencil_type; 18247c6ae99SBarry Smith PetscInt *lx = dd->lx; 18347c6ae99SBarry Smith PetscInt *ly = dd->ly; 18447c6ae99SBarry Smith PetscInt *lz = dd->lz; 18547c6ae99SBarry Smith MPI_Comm comm; 18647c6ae99SBarry Smith PetscMPIInt rank,size; 187ce00eea3SSatish Balay PetscInt xs = 0,xe,ys = 0,ye,zs = 0,ze,x = 0,y = 0,z = 0; 188ce00eea3SSatish Balay PetscInt Xs,Xe,Ys,Ye,Zs,Ze,IXs,IXe,IYs,IYe,IZs,IZe,start,end,pm; 189ce00eea3SSatish Balay PetscInt left,right,up,down,bottom,top,i,j,k,*idx,*idx_cpy,nn; 190ce00eea3SSatish Balay const PetscInt *idx_full; 19147c6ae99SBarry Smith PetscInt n0,n1,n2,n3,n4,n5,n6,n7,n8,n9,n10,n11,n12,n14; 19247c6ae99SBarry Smith PetscInt n15,n16,n17,n18,n19,n20,n21,n22,n23,n24,n25,n26; 193db87c5ecSEthan Coon PetscInt *bases,*ldims,base,x_t,y_t,z_t,s_t,count,s_x,s_y,s_z; 19447c6ae99SBarry Smith PetscInt sn0 = 0,sn1 = 0,sn2 = 0,sn3 = 0,sn5 = 0,sn6 = 0,sn7 = 0; 19547c6ae99SBarry Smith PetscInt sn8 = 0,sn9 = 0,sn11 = 0,sn15 = 0,sn24 = 0,sn25 = 0,sn26 = 0; 19647c6ae99SBarry Smith PetscInt sn17 = 0,sn18 = 0,sn19 = 0,sn20 = 0,sn21 = 0,sn23 = 0; 19747c6ae99SBarry Smith Vec local,global; 19847c6ae99SBarry Smith VecScatter ltog,gtol; 199ce00eea3SSatish Balay IS to,from,ltogis; 2006f951b95Secoon PetscBool twod; 20147c6ae99SBarry Smith PetscErrorCode ierr; 20247c6ae99SBarry Smith 2036f951b95Secoon 20447c6ae99SBarry Smith PetscFunctionBegin; 205c2859e5eSBarry Smith if (stencil_type == DMDA_STENCIL_BOX && (bx == DMDA_BOUNDARY_MIRROR || by == DMDA_BOUNDARY_MIRROR || bz == DMDA_BOUNDARY_MIRROR)) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Mirror boundary and box stencil"); 20647c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject) da, &comm);CHKERRQ(ierr); 2073855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES) 2083855c12bSBarry 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); 2093855c12bSBarry Smith #endif 2103855c12bSBarry Smith 21147c6ae99SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 21247c6ae99SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 21347c6ae99SBarry Smith 21447c6ae99SBarry Smith if (m != PETSC_DECIDE) { 21547c6ae99SBarry Smith if (m < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m); 21647c6ae99SBarry Smith else if (m > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size); 21747c6ae99SBarry Smith } 21847c6ae99SBarry Smith if (n != PETSC_DECIDE) { 21947c6ae99SBarry Smith if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n); 22047c6ae99SBarry Smith else if (n > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size); 22147c6ae99SBarry Smith } 22247c6ae99SBarry Smith if (p != PETSC_DECIDE) { 22347c6ae99SBarry Smith if (p < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Z direction: %D",p); 22447c6ae99SBarry Smith else if (p > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Z direction: %D %d",p,size); 22547c6ae99SBarry Smith } 2260716a85fSBarry 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); 22747c6ae99SBarry Smith 22847c6ae99SBarry Smith /* Partition the array among the processors */ 22947c6ae99SBarry Smith if (m == PETSC_DECIDE && n != PETSC_DECIDE && p != PETSC_DECIDE) { 23047c6ae99SBarry Smith m = size/(n*p); 23147c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) { 23247c6ae99SBarry Smith n = size/(m*p); 23347c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) { 23447c6ae99SBarry Smith p = size/(m*n); 23547c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) { 23647c6ae99SBarry Smith /* try for squarish distribution */ 2378f1a2a5eSBarry Smith m = (int)(0.5 + sqrt(((double)M)*((double)size)/((double)N*p))); 23847c6ae99SBarry Smith if (!m) m = 1; 23947c6ae99SBarry Smith while (m > 0) { 24047c6ae99SBarry Smith n = size/(m*p); 24147c6ae99SBarry Smith if (m*n*p == size) break; 24247c6ae99SBarry Smith m--; 24347c6ae99SBarry Smith } 24447c6ae99SBarry Smith if (!m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad p value: p = %D",p); 24547c6ae99SBarry Smith if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;} 24647c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) { 24747c6ae99SBarry Smith /* try for squarish distribution */ 2488f1a2a5eSBarry Smith m = (int)(0.5 + sqrt(((double)M)*((double)size)/((double)P*n))); 24947c6ae99SBarry Smith if (!m) m = 1; 25047c6ae99SBarry Smith while (m > 0) { 25147c6ae99SBarry Smith p = size/(m*n); 25247c6ae99SBarry Smith if (m*n*p == size) break; 25347c6ae99SBarry Smith m--; 25447c6ae99SBarry Smith } 25547c6ae99SBarry Smith if (!m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad n value: n = %D",n); 25647c6ae99SBarry Smith if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;} 25747c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) { 25847c6ae99SBarry Smith /* try for squarish distribution */ 2598f1a2a5eSBarry Smith n = (int)(0.5 + sqrt(((double)N)*((double)size)/((double)P*m))); 26047c6ae99SBarry Smith if (!n) n = 1; 26147c6ae99SBarry Smith while (n > 0) { 26247c6ae99SBarry Smith p = size/(m*n); 26347c6ae99SBarry Smith if (m*n*p == size) break; 26447c6ae99SBarry Smith n--; 26547c6ae99SBarry Smith } 26647c6ae99SBarry Smith if (!n) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad m value: m = %D",n); 26747c6ae99SBarry Smith if (N > P && n < p) {PetscInt _n = n; n = p; p = _n;} 26847c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) { 26947c6ae99SBarry Smith /* try for squarish distribution */ 2708f1a2a5eSBarry Smith n = (PetscInt)(0.5 + pow(((double)N*N)*((double)size)/((double)P*M),(double)(1./3.))); 27147c6ae99SBarry Smith if (!n) n = 1; 27247c6ae99SBarry Smith while (n > 0) { 27347c6ae99SBarry Smith pm = size/n; 27447c6ae99SBarry Smith if (n*pm == size) break; 27547c6ae99SBarry Smith n--; 27647c6ae99SBarry Smith } 27747c6ae99SBarry Smith if (!n) n = 1; 2788f1a2a5eSBarry Smith m = (PetscInt)(0.5 + sqrt(((double)M)*((double)size)/((double)P*n))); 27947c6ae99SBarry Smith if (!m) m = 1; 28047c6ae99SBarry Smith while (m > 0) { 28147c6ae99SBarry Smith p = size/(m*n); 28247c6ae99SBarry Smith if (m*n*p == size) break; 28347c6ae99SBarry Smith m--; 28447c6ae99SBarry Smith } 28547c6ae99SBarry Smith if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;} 28647c6ae99SBarry Smith } else if (m*n*p != size) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition"); 28747c6ae99SBarry Smith 28847c6ae99SBarry Smith if (m*n*p != size) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_PLIB,"Could not find good partition"); 28947c6ae99SBarry Smith if (M < m) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m); 29047c6ae99SBarry Smith if (N < n) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n); 29147c6ae99SBarry Smith if (P < p) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in z direction is too fine! %D %D",P,p); 29247c6ae99SBarry Smith 29347c6ae99SBarry Smith /* 29447c6ae99SBarry Smith Determine locally owned region 29547c6ae99SBarry Smith [x, y, or z]s is the first local node number, [x, y, z] is the number of local nodes 29647c6ae99SBarry Smith */ 29747c6ae99SBarry Smith 29847c6ae99SBarry Smith if (!lx) { 29947c6ae99SBarry Smith ierr = PetscMalloc(m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr); 30047c6ae99SBarry Smith lx = dd->lx; 301*8865f1eaSKarl Rupp for (i=0; i<m; i++) lx[i] = M/m + ((M % m) > (i % m)); 30247c6ae99SBarry Smith } 30347c6ae99SBarry Smith x = lx[rank % m]; 30447c6ae99SBarry Smith xs = 0; 305*8865f1eaSKarl Rupp for (i=0; i<(rank%m); i++) xs += lx[i]; 306d9c9ebe5SPeter 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); 30747c6ae99SBarry Smith 30847c6ae99SBarry Smith if (!ly) { 30947c6ae99SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr); 31047c6ae99SBarry Smith ly = dd->ly; 311*8865f1eaSKarl Rupp for (i=0; i<n; i++) ly[i] = N/n + ((N % n) > (i % n)); 31247c6ae99SBarry Smith } 31347c6ae99SBarry Smith y = ly[(rank % (m*n))/m]; 314d9c9ebe5SPeter 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); 31530729d88SBarry Smith 31647c6ae99SBarry Smith ys = 0; 317*8865f1eaSKarl Rupp for (i=0; i<(rank % (m*n))/m; i++) ys += ly[i]; 31847c6ae99SBarry Smith 31947c6ae99SBarry Smith if (!lz) { 32047c6ae99SBarry Smith ierr = PetscMalloc(p*sizeof(PetscInt), &dd->lz);CHKERRQ(ierr); 32147c6ae99SBarry Smith lz = dd->lz; 322*8865f1eaSKarl Rupp for (i=0; i<p; i++) lz[i] = P/p + ((P % p) > (i % p)); 32347c6ae99SBarry Smith } 32447c6ae99SBarry Smith z = lz[rank/(m*n)]; 325bcea557cSEthan Coon 326fdc81ce8SEthan Coon /* note this is different than x- and y-, as we will handle as an important special 327fdc81ce8SEthan Coon case when p=P=1 and DMDA_BOUNDARY_PERIODIC and s > z. This is to deal with 2D problems 328fdc81ce8SEthan Coon in a 3D code. Additional code for this case is noted with "2d case" comments */ 3296f951b95Secoon twod = PETSC_FALSE; 330*8865f1eaSKarl Rupp if (P == 1) twod = PETSC_TRUE; 331*8865f1eaSKarl 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); 33247c6ae99SBarry Smith zs = 0; 333*8865f1eaSKarl Rupp for (i=0; i<(rank/(m*n)); i++) zs += lz[i]; 33447c6ae99SBarry Smith ye = ys + y; 33547c6ae99SBarry Smith xe = xs + x; 33647c6ae99SBarry Smith ze = zs + z; 33747c6ae99SBarry Smith 33888661749SPeter Brune /* determine ghost region (Xs) and region scattered into (IXs) */ 339d9c9ebe5SPeter Brune if (xs-s > 0) { 340d9c9ebe5SPeter Brune Xs = xs - s; IXs = xs - s; 34188661749SPeter Brune } else { 342*8865f1eaSKarl Rupp if (bx) Xs = xs - s; 343*8865f1eaSKarl Rupp else Xs = 0; 34488661749SPeter Brune IXs = 0; 34588661749SPeter Brune } 346d9c9ebe5SPeter Brune if (xe+s <= M) { 347d9c9ebe5SPeter Brune Xe = xe + s; IXe = xe + s; 34888661749SPeter Brune } else { 34988661749SPeter Brune if (bx) { 350d9c9ebe5SPeter Brune Xs = xs - s; Xe = xe + s; 351*8865f1eaSKarl Rupp } else Xe = M; 35288661749SPeter Brune IXe = M; 35388661749SPeter Brune } 35447c6ae99SBarry Smith 355c2859e5eSBarry Smith if (bx == DMDA_BOUNDARY_PERIODIC || bx == DMDA_BOUNDARY_MIRROR) { 356d9c9ebe5SPeter Brune IXs = xs - s; 357d9c9ebe5SPeter Brune IXe = xe + s; 358d9c9ebe5SPeter Brune Xs = xs - s; 359d9c9ebe5SPeter Brune Xe = xe + s; 36088661749SPeter Brune } 36188661749SPeter Brune 362d9c9ebe5SPeter Brune if (ys-s > 0) { 363d9c9ebe5SPeter Brune Ys = ys - s; IYs = ys - s; 36488661749SPeter Brune } else { 365*8865f1eaSKarl Rupp if (by) Ys = ys - s; 366*8865f1eaSKarl Rupp else Ys = 0; 36788661749SPeter Brune IYs = 0; 36888661749SPeter Brune } 369d9c9ebe5SPeter Brune if (ye+s <= N) { 370d9c9ebe5SPeter Brune Ye = ye + s; IYe = ye + s; 37188661749SPeter Brune } else { 372*8865f1eaSKarl Rupp if (by) Ye = ye + s; 373*8865f1eaSKarl Rupp else Ye = N; 37488661749SPeter Brune IYe = N; 37588661749SPeter Brune } 37688661749SPeter Brune 377c2859e5eSBarry Smith if (by == DMDA_BOUNDARY_PERIODIC || by == DMDA_BOUNDARY_MIRROR) { 378d9c9ebe5SPeter Brune IYs = ys - s; 379d9c9ebe5SPeter Brune IYe = ye + s; 380d9c9ebe5SPeter Brune Ys = ys - s; 381d9c9ebe5SPeter Brune Ye = ye + s; 38288661749SPeter Brune } 38388661749SPeter Brune 384d9c9ebe5SPeter Brune if (zs-s > 0) { 385d9c9ebe5SPeter Brune Zs = zs - s; IZs = zs - s; 38688661749SPeter Brune } else { 387*8865f1eaSKarl Rupp if (bz) Zs = zs - s; 388*8865f1eaSKarl Rupp else Zs = 0; 38988661749SPeter Brune IZs = 0; 39088661749SPeter Brune } 391d9c9ebe5SPeter Brune if (ze+s <= P) { 392d9c9ebe5SPeter Brune Ze = ze + s; IZe = ze + s; 39388661749SPeter Brune } else { 394*8865f1eaSKarl Rupp if (bz) Ze = ze + s; 395*8865f1eaSKarl Rupp else Ze = P; 39688661749SPeter Brune IZe = P; 39788661749SPeter Brune } 39888661749SPeter Brune 399c2859e5eSBarry Smith if (bz == DMDA_BOUNDARY_PERIODIC || bz == DMDA_BOUNDARY_MIRROR) { 400d9c9ebe5SPeter Brune IZs = zs - s; 401d9c9ebe5SPeter Brune IZe = ze + s; 402d9c9ebe5SPeter Brune Zs = zs - s; 403d9c9ebe5SPeter Brune Ze = ze + s; 40488661749SPeter Brune } 40547c6ae99SBarry Smith 40647c6ae99SBarry Smith /* Resize all X parameters to reflect w */ 407d9c9ebe5SPeter Brune s_x = s; 408d9c9ebe5SPeter Brune s_y = s; 409d9c9ebe5SPeter Brune s_z = s; 41047c6ae99SBarry Smith 41147c6ae99SBarry Smith /* determine starting point of each processor */ 41247c6ae99SBarry Smith nn = x*y*z; 41347c6ae99SBarry Smith ierr = PetscMalloc2(size+1,PetscInt,&bases,size,PetscInt,&ldims);CHKERRQ(ierr); 41447c6ae99SBarry Smith ierr = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr); 41547c6ae99SBarry Smith bases[0] = 0; 416*8865f1eaSKarl Rupp for (i=1; i<=size; i++) bases[i] = ldims[i-1]; 417*8865f1eaSKarl Rupp for (i=1; i<=size; i++) bases[i] += bases[i-1]; 418ce00eea3SSatish Balay base = bases[rank]*dof; 41947c6ae99SBarry Smith 42047c6ae99SBarry Smith /* allocate the base parallel and sequential vectors */ 421ce00eea3SSatish Balay dd->Nlocal = x*y*z*dof; 422778a2246SBarry Smith ierr = VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,0,&global);CHKERRQ(ierr); 423ce00eea3SSatish Balay dd->nlocal = (Xe-Xs)*(Ye-Ys)*(Ze-Zs)*dof; 424778a2246SBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,0,&local);CHKERRQ(ierr); 42547c6ae99SBarry Smith 42647c6ae99SBarry Smith /* generate appropriate vector scatters */ 42747c6ae99SBarry Smith /* local to global inserts non-ghost point region into global */ 42847c6ae99SBarry Smith ierr = VecGetOwnershipRange(global,&start,&end);CHKERRQ(ierr); 429ce00eea3SSatish Balay ierr = ISCreateStride(comm,x*y*z*dof,start,1,&to);CHKERRQ(ierr); 43047c6ae99SBarry Smith 431ce00eea3SSatish Balay ierr = PetscMalloc(x*y*z*sizeof(PetscInt),&idx);CHKERRQ(ierr); 432ce00eea3SSatish Balay left = xs - Xs; right = left + x; 43347c6ae99SBarry Smith bottom = ys - Ys; top = bottom + y; 43447c6ae99SBarry Smith down = zs - Zs; up = down + z; 43547c6ae99SBarry Smith count = 0; 43647c6ae99SBarry Smith for (i=down; i<up; i++) { 43747c6ae99SBarry Smith for (j=bottom; j<top; j++) { 438ce00eea3SSatish Balay for (k=left; k<right; k++) { 439ce00eea3SSatish Balay idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 44047c6ae99SBarry Smith } 44147c6ae99SBarry Smith } 44247c6ae99SBarry Smith } 44347c6ae99SBarry Smith 44447c6ae99SBarry Smith ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&from);CHKERRQ(ierr); 44547c6ae99SBarry Smith ierr = VecScatterCreate(local,from,global,to,<og);CHKERRQ(ierr); 44647c6ae99SBarry Smith ierr = PetscLogObjectParent(da,ltog);CHKERRQ(ierr); 447fcfd50ebSBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 448fcfd50ebSBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 44947c6ae99SBarry Smith 450ce00eea3SSatish Balay /* global to local must include ghost points within the domain, 451ce00eea3SSatish Balay but not ghost points outside the domain that aren't periodic */ 452d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_BOX) { 453db87c5ecSEthan Coon count = (IXe-IXs)*(IYe-IYs)*(IZe-IZs); 454db87c5ecSEthan Coon ierr = PetscMalloc(count*sizeof(PetscInt),&idx);CHKERRQ(ierr); 455ce00eea3SSatish Balay 456ce00eea3SSatish Balay left = IXs - Xs; right = left + (IXe-IXs); 457ce00eea3SSatish Balay bottom = IYs - Ys; top = bottom + (IYe-IYs); 458ce00eea3SSatish Balay down = IZs - Zs; up = down + (IZe-IZs); 459ce00eea3SSatish Balay count = 0; 460ce00eea3SSatish Balay for (i=down; i<up; i++) { 461ce00eea3SSatish Balay for (j=bottom; j<top; j++) { 462ce00eea3SSatish Balay for (k=left; k<right; k++) { 463ce00eea3SSatish Balay idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 464ce00eea3SSatish Balay } 465ce00eea3SSatish Balay } 466ce00eea3SSatish Balay } 467ce00eea3SSatish Balay ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr); 468ce00eea3SSatish Balay 46947c6ae99SBarry Smith } else { 47047c6ae99SBarry Smith /* This is way ugly! We need to list the funny cross type region */ 471db87c5ecSEthan Coon count = ((ys-IYs) + (IYe-ye))*x*z + ((xs-IXs) + (IXe-xe))*y*z + ((zs-IZs) + (IZe-ze))*x*y + x*y*z; 472db87c5ecSEthan Coon ierr = PetscMalloc(count*sizeof(PetscInt),&idx);CHKERRQ(ierr); 473ce00eea3SSatish Balay 474ce00eea3SSatish Balay left = xs - Xs; right = left + x; 47547c6ae99SBarry Smith bottom = ys - Ys; top = bottom + y; 47647c6ae99SBarry Smith down = zs - Zs; up = down + z; 47747c6ae99SBarry Smith count = 0; 478ce00eea3SSatish Balay /* the bottom chunck */ 479ce00eea3SSatish Balay for (i=(IZs-Zs); i<down; i++) { 48047c6ae99SBarry Smith for (j=bottom; j<top; j++) { 481ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 48247c6ae99SBarry Smith } 48347c6ae99SBarry Smith } 48447c6ae99SBarry Smith /* the middle piece */ 48547c6ae99SBarry Smith for (i=down; i<up; i++) { 48647c6ae99SBarry Smith /* front */ 487ce00eea3SSatish Balay for (j=(IYs-Ys); j<bottom; j++) { 488ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 48947c6ae99SBarry Smith } 49047c6ae99SBarry Smith /* middle */ 49147c6ae99SBarry Smith for (j=bottom; j<top; j++) { 492ce00eea3SSatish Balay for (k=IXs-Xs; k<IXe-Xs; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 49347c6ae99SBarry Smith } 49447c6ae99SBarry Smith /* back */ 495ce00eea3SSatish Balay for (j=top; j<top+IYe-ye; j++) { 496ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 49747c6ae99SBarry Smith } 49847c6ae99SBarry Smith } 49947c6ae99SBarry Smith /* the top piece */ 500ce00eea3SSatish Balay for (i=up; i<up+IZe-ze; i++) { 50147c6ae99SBarry Smith for (j=bottom; j<top; j++) { 502ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 50347c6ae99SBarry Smith } 50447c6ae99SBarry Smith } 50547c6ae99SBarry Smith ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr); 50647c6ae99SBarry Smith } 50747c6ae99SBarry Smith 50847c6ae99SBarry Smith /* determine who lies on each side of use stored in n24 n25 n26 50947c6ae99SBarry Smith n21 n22 n23 51047c6ae99SBarry Smith n18 n19 n20 51147c6ae99SBarry Smith 51247c6ae99SBarry Smith n15 n16 n17 51347c6ae99SBarry Smith n12 n14 51447c6ae99SBarry Smith n9 n10 n11 51547c6ae99SBarry Smith 51647c6ae99SBarry Smith n6 n7 n8 51747c6ae99SBarry Smith n3 n4 n5 51847c6ae99SBarry Smith n0 n1 n2 51947c6ae99SBarry Smith */ 52047c6ae99SBarry Smith 52147c6ae99SBarry Smith /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */ 52247c6ae99SBarry Smith /* Assume Nodes are Internal to the Cube */ 52347c6ae99SBarry Smith n0 = rank - m*n - m - 1; 52447c6ae99SBarry Smith n1 = rank - m*n - m; 52547c6ae99SBarry Smith n2 = rank - m*n - m + 1; 52647c6ae99SBarry Smith n3 = rank - m*n -1; 52747c6ae99SBarry Smith n4 = rank - m*n; 52847c6ae99SBarry Smith n5 = rank - m*n + 1; 52947c6ae99SBarry Smith n6 = rank - m*n + m - 1; 53047c6ae99SBarry Smith n7 = rank - m*n + m; 53147c6ae99SBarry Smith n8 = rank - m*n + m + 1; 53247c6ae99SBarry Smith 53347c6ae99SBarry Smith n9 = rank - m - 1; 53447c6ae99SBarry Smith n10 = rank - m; 53547c6ae99SBarry Smith n11 = rank - m + 1; 53647c6ae99SBarry Smith n12 = rank - 1; 53747c6ae99SBarry Smith n14 = rank + 1; 53847c6ae99SBarry Smith n15 = rank + m - 1; 53947c6ae99SBarry Smith n16 = rank + m; 54047c6ae99SBarry Smith n17 = rank + m + 1; 54147c6ae99SBarry Smith 54247c6ae99SBarry Smith n18 = rank + m*n - m - 1; 54347c6ae99SBarry Smith n19 = rank + m*n - m; 54447c6ae99SBarry Smith n20 = rank + m*n - m + 1; 54547c6ae99SBarry Smith n21 = rank + m*n - 1; 54647c6ae99SBarry Smith n22 = rank + m*n; 54747c6ae99SBarry Smith n23 = rank + m*n + 1; 54847c6ae99SBarry Smith n24 = rank + m*n + m - 1; 54947c6ae99SBarry Smith n25 = rank + m*n + m; 55047c6ae99SBarry Smith n26 = rank + m*n + m + 1; 55147c6ae99SBarry Smith 55247c6ae99SBarry Smith /* Assume Pieces are on Faces of Cube */ 55347c6ae99SBarry Smith 55447c6ae99SBarry Smith if (xs == 0) { /* First assume not corner or edge */ 55547c6ae99SBarry Smith n0 = rank -1 - (m*n); 55647c6ae99SBarry Smith n3 = rank + m -1 - (m*n); 55747c6ae99SBarry Smith n6 = rank + 2*m -1 - (m*n); 55847c6ae99SBarry Smith n9 = rank -1; 55947c6ae99SBarry Smith n12 = rank + m -1; 56047c6ae99SBarry Smith n15 = rank + 2*m -1; 56147c6ae99SBarry Smith n18 = rank -1 + (m*n); 56247c6ae99SBarry Smith n21 = rank + m -1 + (m*n); 56347c6ae99SBarry Smith n24 = rank + 2*m -1 + (m*n); 56447c6ae99SBarry Smith } 56547c6ae99SBarry Smith 566ce00eea3SSatish Balay if (xe == M) { /* First assume not corner or edge */ 56747c6ae99SBarry Smith n2 = rank -2*m +1 - (m*n); 56847c6ae99SBarry Smith n5 = rank - m +1 - (m*n); 56947c6ae99SBarry Smith n8 = rank +1 - (m*n); 57047c6ae99SBarry Smith n11 = rank -2*m +1; 57147c6ae99SBarry Smith n14 = rank - m +1; 57247c6ae99SBarry Smith n17 = rank +1; 57347c6ae99SBarry Smith n20 = rank -2*m +1 + (m*n); 57447c6ae99SBarry Smith n23 = rank - m +1 + (m*n); 57547c6ae99SBarry Smith n26 = rank +1 + (m*n); 57647c6ae99SBarry Smith } 57747c6ae99SBarry Smith 57847c6ae99SBarry Smith if (ys==0) { /* First assume not corner or edge */ 57947c6ae99SBarry Smith n0 = rank + m * (n-1) -1 - (m*n); 58047c6ae99SBarry Smith n1 = rank + m * (n-1) - (m*n); 58147c6ae99SBarry Smith n2 = rank + m * (n-1) +1 - (m*n); 58247c6ae99SBarry Smith n9 = rank + m * (n-1) -1; 58347c6ae99SBarry Smith n10 = rank + m * (n-1); 58447c6ae99SBarry Smith n11 = rank + m * (n-1) +1; 58547c6ae99SBarry Smith n18 = rank + m * (n-1) -1 + (m*n); 58647c6ae99SBarry Smith n19 = rank + m * (n-1) + (m*n); 58747c6ae99SBarry Smith n20 = rank + m * (n-1) +1 + (m*n); 58847c6ae99SBarry Smith } 58947c6ae99SBarry Smith 59047c6ae99SBarry Smith if (ye == N) { /* First assume not corner or edge */ 59147c6ae99SBarry Smith n6 = rank - m * (n-1) -1 - (m*n); 59247c6ae99SBarry Smith n7 = rank - m * (n-1) - (m*n); 59347c6ae99SBarry Smith n8 = rank - m * (n-1) +1 - (m*n); 59447c6ae99SBarry Smith n15 = rank - m * (n-1) -1; 59547c6ae99SBarry Smith n16 = rank - m * (n-1); 59647c6ae99SBarry Smith n17 = rank - m * (n-1) +1; 59747c6ae99SBarry Smith n24 = rank - m * (n-1) -1 + (m*n); 59847c6ae99SBarry Smith n25 = rank - m * (n-1) + (m*n); 59947c6ae99SBarry Smith n26 = rank - m * (n-1) +1 + (m*n); 60047c6ae99SBarry Smith } 60147c6ae99SBarry Smith 60247c6ae99SBarry Smith if (zs == 0) { /* First assume not corner or edge */ 60347c6ae99SBarry Smith n0 = size - (m*n) + rank - m - 1; 60447c6ae99SBarry Smith n1 = size - (m*n) + rank - m; 60547c6ae99SBarry Smith n2 = size - (m*n) + rank - m + 1; 60647c6ae99SBarry Smith n3 = size - (m*n) + rank - 1; 60747c6ae99SBarry Smith n4 = size - (m*n) + rank; 60847c6ae99SBarry Smith n5 = size - (m*n) + rank + 1; 60947c6ae99SBarry Smith n6 = size - (m*n) + rank + m - 1; 61047c6ae99SBarry Smith n7 = size - (m*n) + rank + m; 61147c6ae99SBarry Smith n8 = size - (m*n) + rank + m + 1; 61247c6ae99SBarry Smith } 61347c6ae99SBarry Smith 61447c6ae99SBarry Smith if (ze == P) { /* First assume not corner or edge */ 61547c6ae99SBarry Smith n18 = (m*n) - (size-rank) - m - 1; 61647c6ae99SBarry Smith n19 = (m*n) - (size-rank) - m; 61747c6ae99SBarry Smith n20 = (m*n) - (size-rank) - m + 1; 61847c6ae99SBarry Smith n21 = (m*n) - (size-rank) - 1; 61947c6ae99SBarry Smith n22 = (m*n) - (size-rank); 62047c6ae99SBarry Smith n23 = (m*n) - (size-rank) + 1; 62147c6ae99SBarry Smith n24 = (m*n) - (size-rank) + m - 1; 62247c6ae99SBarry Smith n25 = (m*n) - (size-rank) + m; 62347c6ae99SBarry Smith n26 = (m*n) - (size-rank) + m + 1; 62447c6ae99SBarry Smith } 62547c6ae99SBarry Smith 62647c6ae99SBarry Smith if ((xs==0) && (zs==0)) { /* Assume an edge, not corner */ 62747c6ae99SBarry Smith n0 = size - m*n + rank + m-1 - m; 62847c6ae99SBarry Smith n3 = size - m*n + rank + m-1; 62947c6ae99SBarry Smith n6 = size - m*n + rank + m-1 + m; 63047c6ae99SBarry Smith } 63147c6ae99SBarry Smith 63247c6ae99SBarry Smith if ((xs==0) && (ze==P)) { /* Assume an edge, not corner */ 63347c6ae99SBarry Smith n18 = m*n - (size - rank) + m-1 - m; 63447c6ae99SBarry Smith n21 = m*n - (size - rank) + m-1; 63547c6ae99SBarry Smith n24 = m*n - (size - rank) + m-1 + m; 63647c6ae99SBarry Smith } 63747c6ae99SBarry Smith 63847c6ae99SBarry Smith if ((xs==0) && (ys==0)) { /* Assume an edge, not corner */ 63947c6ae99SBarry Smith n0 = rank + m*n -1 - m*n; 64047c6ae99SBarry Smith n9 = rank + m*n -1; 64147c6ae99SBarry Smith n18 = rank + m*n -1 + m*n; 64247c6ae99SBarry Smith } 64347c6ae99SBarry Smith 64447c6ae99SBarry Smith if ((xs==0) && (ye==N)) { /* Assume an edge, not corner */ 64547c6ae99SBarry Smith n6 = rank - m*(n-1) + m-1 - m*n; 64647c6ae99SBarry Smith n15 = rank - m*(n-1) + m-1; 64747c6ae99SBarry Smith n24 = rank - m*(n-1) + m-1 + m*n; 64847c6ae99SBarry Smith } 64947c6ae99SBarry Smith 650ce00eea3SSatish Balay if ((xe==M) && (zs==0)) { /* Assume an edge, not corner */ 65147c6ae99SBarry Smith n2 = size - (m*n-rank) - (m-1) - m; 65247c6ae99SBarry Smith n5 = size - (m*n-rank) - (m-1); 65347c6ae99SBarry Smith n8 = size - (m*n-rank) - (m-1) + m; 65447c6ae99SBarry Smith } 65547c6ae99SBarry Smith 656ce00eea3SSatish Balay if ((xe==M) && (ze==P)) { /* Assume an edge, not corner */ 65747c6ae99SBarry Smith n20 = m*n - (size - rank) - (m-1) - m; 65847c6ae99SBarry Smith n23 = m*n - (size - rank) - (m-1); 65947c6ae99SBarry Smith n26 = m*n - (size - rank) - (m-1) + m; 66047c6ae99SBarry Smith } 66147c6ae99SBarry Smith 662ce00eea3SSatish Balay if ((xe==M) && (ys==0)) { /* Assume an edge, not corner */ 66347c6ae99SBarry Smith n2 = rank + m*(n-1) - (m-1) - m*n; 66447c6ae99SBarry Smith n11 = rank + m*(n-1) - (m-1); 66547c6ae99SBarry Smith n20 = rank + m*(n-1) - (m-1) + m*n; 66647c6ae99SBarry Smith } 66747c6ae99SBarry Smith 668ce00eea3SSatish Balay if ((xe==M) && (ye==N)) { /* Assume an edge, not corner */ 66947c6ae99SBarry Smith n8 = rank - m*n +1 - m*n; 67047c6ae99SBarry Smith n17 = rank - m*n +1; 67147c6ae99SBarry Smith n26 = rank - m*n +1 + m*n; 67247c6ae99SBarry Smith } 67347c6ae99SBarry Smith 67447c6ae99SBarry Smith if ((ys==0) && (zs==0)) { /* Assume an edge, not corner */ 67547c6ae99SBarry Smith n0 = size - m + rank -1; 67647c6ae99SBarry Smith n1 = size - m + rank; 67747c6ae99SBarry Smith n2 = size - m + rank +1; 67847c6ae99SBarry Smith } 67947c6ae99SBarry Smith 68047c6ae99SBarry Smith if ((ys==0) && (ze==P)) { /* Assume an edge, not corner */ 68147c6ae99SBarry Smith n18 = m*n - (size - rank) + m*(n-1) -1; 68247c6ae99SBarry Smith n19 = m*n - (size - rank) + m*(n-1); 68347c6ae99SBarry Smith n20 = m*n - (size - rank) + m*(n-1) +1; 68447c6ae99SBarry Smith } 68547c6ae99SBarry Smith 68647c6ae99SBarry Smith if ((ye==N) && (zs==0)) { /* Assume an edge, not corner */ 68747c6ae99SBarry Smith n6 = size - (m*n-rank) - m * (n-1) -1; 68847c6ae99SBarry Smith n7 = size - (m*n-rank) - m * (n-1); 68947c6ae99SBarry Smith n8 = size - (m*n-rank) - m * (n-1) +1; 69047c6ae99SBarry Smith } 69147c6ae99SBarry Smith 69247c6ae99SBarry Smith if ((ye==N) && (ze==P)) { /* Assume an edge, not corner */ 69347c6ae99SBarry Smith n24 = rank - (size-m) -1; 69447c6ae99SBarry Smith n25 = rank - (size-m); 69547c6ae99SBarry Smith n26 = rank - (size-m) +1; 69647c6ae99SBarry Smith } 69747c6ae99SBarry Smith 69847c6ae99SBarry Smith /* Check for Corners */ 699*8865f1eaSKarl Rupp if ((xs==0) && (ys==0) && (zs==0)) n0 = size -1; 700*8865f1eaSKarl Rupp if ((xs==0) && (ys==0) && (ze==P)) n18 = m*n-1; 701*8865f1eaSKarl Rupp if ((xs==0) && (ye==N) && (zs==0)) n6 = (size-1)-m*(n-1); 702*8865f1eaSKarl Rupp if ((xs==0) && (ye==N) && (ze==P)) n24 = m-1; 703*8865f1eaSKarl Rupp if ((xe==M) && (ys==0) && (zs==0)) n2 = size-m; 704*8865f1eaSKarl Rupp if ((xe==M) && (ys==0) && (ze==P)) n20 = m*n-m; 705*8865f1eaSKarl Rupp if ((xe==M) && (ye==N) && (zs==0)) n8 = size-m*n; 706*8865f1eaSKarl Rupp if ((xe==M) && (ye==N) && (ze==P)) n26 = 0; 70747c6ae99SBarry Smith 70847c6ae99SBarry Smith /* Check for when not X,Y, and Z Periodic */ 70947c6ae99SBarry Smith 71047c6ae99SBarry Smith /* If not X periodic */ 7111321219cSEthan Coon if (bx != DMDA_BOUNDARY_PERIODIC) { 712*8865f1eaSKarl Rupp if (xs==0) n0 = n3 = n6 = n9 = n12 = n15 = n18 = n21 = n24 = -2; 713*8865f1eaSKarl Rupp if (xe==M) n2 = n5 = n8 = n11 = n14 = n17 = n20 = n23 = n26 = -2; 71447c6ae99SBarry Smith } 71547c6ae99SBarry Smith 71647c6ae99SBarry Smith /* If not Y periodic */ 7171321219cSEthan Coon if (by != DMDA_BOUNDARY_PERIODIC) { 718*8865f1eaSKarl Rupp if (ys==0) n0 = n1 = n2 = n9 = n10 = n11 = n18 = n19 = n20 = -2; 719*8865f1eaSKarl Rupp if (ye==N) n6 = n7 = n8 = n15 = n16 = n17 = n24 = n25 = n26 = -2; 72047c6ae99SBarry Smith } 72147c6ae99SBarry Smith 72247c6ae99SBarry Smith /* If not Z periodic */ 7231321219cSEthan Coon if (bz != DMDA_BOUNDARY_PERIODIC) { 724*8865f1eaSKarl Rupp if (zs==0) n0 = n1 = n2 = n3 = n4 = n5 = n6 = n7 = n8 = -2; 725*8865f1eaSKarl Rupp if (ze==P) n18 = n19 = n20 = n21 = n22 = n23 = n24 = n25 = n26 = -2; 72647c6ae99SBarry Smith } 72747c6ae99SBarry Smith 72847c6ae99SBarry Smith ierr = PetscMalloc(27*sizeof(PetscInt),&dd->neighbors);CHKERRQ(ierr); 729*8865f1eaSKarl Rupp 73047c6ae99SBarry Smith dd->neighbors[0] = n0; 73147c6ae99SBarry Smith dd->neighbors[1] = n1; 73247c6ae99SBarry Smith dd->neighbors[2] = n2; 73347c6ae99SBarry Smith dd->neighbors[3] = n3; 73447c6ae99SBarry Smith dd->neighbors[4] = n4; 73547c6ae99SBarry Smith dd->neighbors[5] = n5; 73647c6ae99SBarry Smith dd->neighbors[6] = n6; 73747c6ae99SBarry Smith dd->neighbors[7] = n7; 73847c6ae99SBarry Smith dd->neighbors[8] = n8; 73947c6ae99SBarry Smith dd->neighbors[9] = n9; 74047c6ae99SBarry Smith dd->neighbors[10] = n10; 74147c6ae99SBarry Smith dd->neighbors[11] = n11; 74247c6ae99SBarry Smith dd->neighbors[12] = n12; 74347c6ae99SBarry Smith dd->neighbors[13] = rank; 74447c6ae99SBarry Smith dd->neighbors[14] = n14; 74547c6ae99SBarry Smith dd->neighbors[15] = n15; 74647c6ae99SBarry Smith dd->neighbors[16] = n16; 74747c6ae99SBarry Smith dd->neighbors[17] = n17; 74847c6ae99SBarry Smith dd->neighbors[18] = n18; 74947c6ae99SBarry Smith dd->neighbors[19] = n19; 75047c6ae99SBarry Smith dd->neighbors[20] = n20; 75147c6ae99SBarry Smith dd->neighbors[21] = n21; 75247c6ae99SBarry Smith dd->neighbors[22] = n22; 75347c6ae99SBarry Smith dd->neighbors[23] = n23; 75447c6ae99SBarry Smith dd->neighbors[24] = n24; 75547c6ae99SBarry Smith dd->neighbors[25] = n25; 75647c6ae99SBarry Smith dd->neighbors[26] = n26; 75747c6ae99SBarry Smith 75847c6ae99SBarry Smith /* If star stencil then delete the corner neighbors */ 759d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 76047c6ae99SBarry Smith /* save information about corner neighbors */ 76147c6ae99SBarry Smith sn0 = n0; sn1 = n1; sn2 = n2; sn3 = n3; sn5 = n5; sn6 = n6; sn7 = n7; 76247c6ae99SBarry Smith sn8 = n8; sn9 = n9; sn11 = n11; sn15 = n15; sn17 = n17; sn18 = n18; 76347c6ae99SBarry Smith sn19 = n19; sn20 = n20; sn21 = n21; sn23 = n23; sn24 = n24; sn25 = n25; 76447c6ae99SBarry Smith sn26 = n26; 765*8865f1eaSKarl Rupp n0 = n1 = n2 = n3 = n5 = n6 = n7 = n8 = n9 = n11 = n15 = n17 = n18 = n19 = n20 = n21 = n23 = n24 = n25 = n26 = -1; 76647c6ae99SBarry Smith } 76747c6ae99SBarry Smith 76847c6ae99SBarry Smith ierr = PetscMalloc((Xe-Xs)*(Ye-Ys)*(Ze-Zs)*sizeof(PetscInt),&idx);CHKERRQ(ierr); 76947c6ae99SBarry Smith ierr = PetscLogObjectMemory(da,(Xe-Xs)*(Ye-Ys)*(Ze-Zs)*sizeof(PetscInt));CHKERRQ(ierr); 77047c6ae99SBarry Smith 77147c6ae99SBarry Smith nn = 0; 77247c6ae99SBarry Smith /* Bottom Level */ 77347c6ae99SBarry Smith for (k=0; k<s_z; k++) { 77447c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 77547c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 776ce00eea3SSatish Balay x_t = lx[n0 % m]; 77747c6ae99SBarry Smith y_t = ly[(n0 % (m*n))/m]; 77847c6ae99SBarry Smith z_t = lz[n0 / (m*n)]; 77947c6ae99SBarry 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; 780*8865f1eaSKarl 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 */ 781*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 78247c6ae99SBarry Smith } 78347c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 78447c6ae99SBarry Smith x_t = x; 78547c6ae99SBarry Smith y_t = ly[(n1 % (m*n))/m]; 78647c6ae99SBarry Smith z_t = lz[n1 / (m*n)]; 78747c6ae99SBarry 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; 788*8865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t; /* 2D case */ 789*8865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 79047c6ae99SBarry Smith } 79147c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 792ce00eea3SSatish Balay x_t = lx[n2 % m]; 79347c6ae99SBarry Smith y_t = ly[(n2 % (m*n))/m]; 79447c6ae99SBarry Smith z_t = lz[n2 / (m*n)]; 79547c6ae99SBarry 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; 796*8865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t; /* 2D case */ 797*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 79847c6ae99SBarry Smith } 79947c6ae99SBarry Smith } 80047c6ae99SBarry Smith 80147c6ae99SBarry Smith for (i=0; i<y; i++) { 80247c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 803ce00eea3SSatish Balay x_t = lx[n3 % m]; 80447c6ae99SBarry Smith y_t = y; 80547c6ae99SBarry Smith z_t = lz[n3 / (m*n)]; 80647c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 807*8865f1eaSKarl 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 */ 808*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 80947c6ae99SBarry Smith } 81047c6ae99SBarry Smith 81147c6ae99SBarry Smith if (n4 >= 0) { /* middle */ 81247c6ae99SBarry Smith x_t = x; 81347c6ae99SBarry Smith y_t = y; 81447c6ae99SBarry Smith z_t = lz[n4 / (m*n)]; 81547c6ae99SBarry Smith s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 816*8865f1eaSKarl 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 */ 817*8865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 818c2859e5eSBarry Smith } else if (bz == DMDA_BOUNDARY_MIRROR) { 819*8865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 82047c6ae99SBarry Smith } 82147c6ae99SBarry Smith 82247c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 823ce00eea3SSatish Balay x_t = lx[n5 % m]; 82447c6ae99SBarry Smith y_t = y; 82547c6ae99SBarry Smith z_t = lz[n5 / (m*n)]; 82647c6ae99SBarry Smith s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 827*8865f1eaSKarl 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 */ 828*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 82947c6ae99SBarry Smith } 83047c6ae99SBarry Smith } 83147c6ae99SBarry Smith 83247c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 83347c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 834ce00eea3SSatish Balay x_t = lx[n6 % m]; 83547c6ae99SBarry Smith y_t = ly[(n6 % (m*n))/m]; 83647c6ae99SBarry Smith z_t = lz[n6 / (m*n)]; 83747c6ae99SBarry Smith s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 838*8865f1eaSKarl 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 */ 839*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 84047c6ae99SBarry Smith } 84147c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 84247c6ae99SBarry Smith x_t = x; 84347c6ae99SBarry Smith y_t = ly[(n7 % (m*n))/m]; 84447c6ae99SBarry Smith z_t = lz[n7 / (m*n)]; 84547c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 846*8865f1eaSKarl 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 */ 847*8865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 84847c6ae99SBarry Smith } 84947c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 850ce00eea3SSatish Balay x_t = lx[n8 % m]; 85147c6ae99SBarry Smith y_t = ly[(n8 % (m*n))/m]; 85247c6ae99SBarry Smith z_t = lz[n8 / (m*n)]; 85347c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 854*8865f1eaSKarl 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 */ 855*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 85647c6ae99SBarry Smith } 85747c6ae99SBarry Smith } 85847c6ae99SBarry Smith } 85947c6ae99SBarry Smith 86047c6ae99SBarry Smith /* Middle Level */ 86147c6ae99SBarry Smith for (k=0; k<z; k++) { 86247c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 86347c6ae99SBarry Smith if (n9 >= 0) { /* left below */ 864ce00eea3SSatish Balay x_t = lx[n9 % m]; 86547c6ae99SBarry Smith y_t = ly[(n9 % (m*n))/m]; 86647c6ae99SBarry Smith /* z_t = z; */ 86747c6ae99SBarry Smith s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 868*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 86947c6ae99SBarry Smith } 87047c6ae99SBarry Smith if (n10 >= 0) { /* directly below */ 87147c6ae99SBarry Smith x_t = x; 87247c6ae99SBarry Smith y_t = ly[(n10 % (m*n))/m]; 87347c6ae99SBarry Smith /* z_t = z; */ 87447c6ae99SBarry Smith s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 875*8865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 876c2859e5eSBarry Smith } else if (by == DMDA_BOUNDARY_MIRROR) { 877*8865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 87847c6ae99SBarry Smith } 87947c6ae99SBarry Smith if (n11 >= 0) { /* right below */ 880ce00eea3SSatish Balay x_t = lx[n11 % m]; 88147c6ae99SBarry Smith y_t = ly[(n11 % (m*n))/m]; 88247c6ae99SBarry Smith /* z_t = z; */ 88347c6ae99SBarry Smith s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 884*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 88547c6ae99SBarry Smith } 88647c6ae99SBarry Smith } 88747c6ae99SBarry Smith 88847c6ae99SBarry Smith for (i=0; i<y; i++) { 88947c6ae99SBarry Smith if (n12 >= 0) { /* directly left */ 890ce00eea3SSatish Balay x_t = lx[n12 % m]; 89147c6ae99SBarry Smith y_t = y; 89247c6ae99SBarry Smith /* z_t = z; */ 89347c6ae99SBarry Smith s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t; 894*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 895c2859e5eSBarry Smith } else if (bx == DMDA_BOUNDARY_MIRROR) { 896*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = 0; 89747c6ae99SBarry Smith } 89847c6ae99SBarry Smith 89947c6ae99SBarry Smith /* Interior */ 90047c6ae99SBarry Smith s_t = bases[rank] + i*x + k*x*y; 901*8865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = s_t++; 90247c6ae99SBarry Smith 90347c6ae99SBarry Smith if (n14 >= 0) { /* directly right */ 904ce00eea3SSatish Balay x_t = lx[n14 % m]; 90547c6ae99SBarry Smith y_t = y; 90647c6ae99SBarry Smith /* z_t = z; */ 90747c6ae99SBarry Smith s_t = bases[n14] + i*x_t + k*x_t*y_t; 908*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 909c2859e5eSBarry Smith } else if (bx == DMDA_BOUNDARY_MIRROR) { 910*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = 0; 91147c6ae99SBarry Smith } 91247c6ae99SBarry Smith } 91347c6ae99SBarry Smith 91447c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 91547c6ae99SBarry Smith if (n15 >= 0) { /* left above */ 916ce00eea3SSatish Balay x_t = lx[n15 % m]; 91747c6ae99SBarry Smith y_t = ly[(n15 % (m*n))/m]; 91847c6ae99SBarry Smith /* z_t = z; */ 91947c6ae99SBarry Smith s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t; 920*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 92147c6ae99SBarry Smith } 92247c6ae99SBarry Smith if (n16 >= 0) { /* directly above */ 92347c6ae99SBarry Smith x_t = x; 92447c6ae99SBarry Smith y_t = ly[(n16 % (m*n))/m]; 92547c6ae99SBarry Smith /* z_t = z; */ 92647c6ae99SBarry Smith s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t; 927*8865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 928c2859e5eSBarry Smith } else if (by == DMDA_BOUNDARY_MIRROR) { 929*8865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 93047c6ae99SBarry Smith } 93147c6ae99SBarry Smith if (n17 >= 0) { /* right above */ 932ce00eea3SSatish Balay x_t = lx[n17 % m]; 93347c6ae99SBarry Smith y_t = ly[(n17 % (m*n))/m]; 93447c6ae99SBarry Smith /* z_t = z; */ 93547c6ae99SBarry Smith s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t; 936*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 93747c6ae99SBarry Smith } 93847c6ae99SBarry Smith } 93947c6ae99SBarry Smith } 94047c6ae99SBarry Smith 94147c6ae99SBarry Smith /* Upper Level */ 94247c6ae99SBarry Smith for (k=0; k<s_z; k++) { 94347c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 94447c6ae99SBarry Smith if (n18 >= 0) { /* left below */ 945ce00eea3SSatish Balay x_t = lx[n18 % m]; 94647c6ae99SBarry Smith y_t = ly[(n18 % (m*n))/m]; 94747c6ae99SBarry Smith /* z_t = lz[n18 / (m*n)]; */ 94847c6ae99SBarry Smith s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 949*8865f1eaSKarl 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 */ 950*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 95147c6ae99SBarry Smith } 95247c6ae99SBarry Smith if (n19 >= 0) { /* directly below */ 95347c6ae99SBarry Smith x_t = x; 95447c6ae99SBarry Smith y_t = ly[(n19 % (m*n))/m]; 95547c6ae99SBarry Smith /* z_t = lz[n19 / (m*n)]; */ 95647c6ae99SBarry Smith s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 957*8865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n19] - (s_y+1-i)*x_t + x_t*y_t; /* 2d case */ 958*8865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 95947c6ae99SBarry Smith } 96047c6ae99SBarry Smith if (n20 >= 0) { /* right below */ 961ce00eea3SSatish Balay x_t = lx[n20 % m]; 96247c6ae99SBarry Smith y_t = ly[(n20 % (m*n))/m]; 96347c6ae99SBarry Smith /* z_t = lz[n20 / (m*n)]; */ 96447c6ae99SBarry Smith s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 965*8865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n20] - (s_y+1-i)*x_t + x_t*y_t; /* 2d case */ 966*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 96747c6ae99SBarry Smith } 96847c6ae99SBarry Smith } 96947c6ae99SBarry Smith 97047c6ae99SBarry Smith for (i=0; i<y; i++) { 97147c6ae99SBarry Smith if (n21 >= 0) { /* directly left */ 972ce00eea3SSatish Balay x_t = lx[n21 % m]; 97347c6ae99SBarry Smith y_t = y; 97447c6ae99SBarry Smith /* z_t = lz[n21 / (m*n)]; */ 97547c6ae99SBarry Smith s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t; 976*8865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n21] + (i+1)*x_t - s_x; /* 2d case */ 977*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 97847c6ae99SBarry Smith } 97947c6ae99SBarry Smith 98047c6ae99SBarry Smith if (n22 >= 0) { /* middle */ 98147c6ae99SBarry Smith x_t = x; 98247c6ae99SBarry Smith y_t = y; 98347c6ae99SBarry Smith /* z_t = lz[n22 / (m*n)]; */ 98447c6ae99SBarry Smith s_t = bases[n22] + i*x_t + k*x_t*y_t; 985*8865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n22] + i*x_t; /* 2d case */ 986*8865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 987c2859e5eSBarry Smith } else if (bz == DMDA_BOUNDARY_MIRROR) { 988*8865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 98947c6ae99SBarry Smith } 99047c6ae99SBarry Smith 99147c6ae99SBarry Smith if (n23 >= 0) { /* directly right */ 992ce00eea3SSatish Balay x_t = lx[n23 % m]; 99347c6ae99SBarry Smith y_t = y; 99447c6ae99SBarry Smith /* z_t = lz[n23 / (m*n)]; */ 99547c6ae99SBarry Smith s_t = bases[n23] + i*x_t + k*x_t*y_t; 996*8865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n23] + i*x_t; /* 2d case */ 997*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 99847c6ae99SBarry Smith } 99947c6ae99SBarry Smith } 100047c6ae99SBarry Smith 100147c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 100247c6ae99SBarry Smith if (n24 >= 0) { /* left above */ 1003ce00eea3SSatish Balay x_t = lx[n24 % m]; 100447c6ae99SBarry Smith y_t = ly[(n24 % (m*n))/m]; 100547c6ae99SBarry Smith /* z_t = lz[n24 / (m*n)]; */ 100647c6ae99SBarry Smith s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t; 1007*8865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n24] + i*x_t - s_x; /* 2d case */ 1008*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 100947c6ae99SBarry Smith } 101047c6ae99SBarry Smith if (n25 >= 0) { /* directly above */ 101147c6ae99SBarry Smith x_t = x; 101247c6ae99SBarry Smith y_t = ly[(n25 % (m*n))/m]; 101347c6ae99SBarry Smith /* z_t = lz[n25 / (m*n)]; */ 101447c6ae99SBarry Smith s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t; 1015*8865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n25] + (i-1)*x_t; /* 2d case */ 1016*8865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 101747c6ae99SBarry Smith } 101847c6ae99SBarry Smith if (n26 >= 0) { /* right above */ 1019ce00eea3SSatish Balay x_t = lx[n26 % m]; 102047c6ae99SBarry Smith y_t = ly[(n26 % (m*n))/m]; 102147c6ae99SBarry Smith /* z_t = lz[n26 / (m*n)]; */ 102247c6ae99SBarry Smith s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t; 1023*8865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n26] + (i-1)*x_t; /* 2d case */ 1024*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 102547c6ae99SBarry Smith } 102647c6ae99SBarry Smith } 102747c6ae99SBarry Smith } 1028ce00eea3SSatish Balay 1029ce00eea3SSatish Balay ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_COPY_VALUES,&from);CHKERRQ(ierr); 103047c6ae99SBarry Smith ierr = VecScatterCreate(global,from,local,to,>ol);CHKERRQ(ierr); 103147c6ae99SBarry Smith ierr = PetscLogObjectParent(da,gtol);CHKERRQ(ierr); 1032fcfd50ebSBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 1033fcfd50ebSBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 103447c6ae99SBarry Smith 1035d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 103647c6ae99SBarry Smith n0 = sn0; n1 = sn1; n2 = sn2; n3 = sn3; n5 = sn5; n6 = sn6; n7 = sn7; 103747c6ae99SBarry Smith n8 = sn8; n9 = sn9; n11 = sn11; n15 = sn15; n17 = sn17; n18 = sn18; 103847c6ae99SBarry Smith n19 = sn19; n20 = sn20; n21 = sn21; n23 = sn23; n24 = sn24; n25 = sn25; 103947c6ae99SBarry Smith n26 = sn26; 1040ce00eea3SSatish Balay } 104147c6ae99SBarry Smith 104288661749SPeter Brune if (((stencil_type == DMDA_STENCIL_STAR) || 10431321219cSEthan Coon (bx != DMDA_BOUNDARY_PERIODIC && bx) || 10441321219cSEthan Coon (by != DMDA_BOUNDARY_PERIODIC && by) || 1045d9c9ebe5SPeter Brune (bz != DMDA_BOUNDARY_PERIODIC && bz))) { 1046ce00eea3SSatish Balay /* 1047ce00eea3SSatish Balay Recompute the local to global mappings, this time keeping the 1048ce00eea3SSatish Balay information about the cross corner processor numbers. 1049ce00eea3SSatish Balay */ 105047c6ae99SBarry Smith nn = 0; 105147c6ae99SBarry Smith /* Bottom Level */ 105247c6ae99SBarry Smith for (k=0; k<s_z; k++) { 105347c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 105447c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 1055ce00eea3SSatish Balay x_t = lx[n0 % m]; 105647c6ae99SBarry Smith y_t = ly[(n0 % (m*n))/m]; 105747c6ae99SBarry Smith z_t = lz[n0 / (m*n)]; 105847c6ae99SBarry 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; 1059*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1060ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0 && Zs-zs < 0) { 1061*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 106247c6ae99SBarry Smith } 106347c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 106447c6ae99SBarry Smith x_t = x; 106547c6ae99SBarry Smith y_t = ly[(n1 % (m*n))/m]; 106647c6ae99SBarry Smith z_t = lz[n1 / (m*n)]; 106747c6ae99SBarry 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; 1068*8865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1069ce00eea3SSatish Balay } else if (Ys-ys < 0 && Zs-zs < 0) { 1070*8865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 107147c6ae99SBarry Smith } 107247c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 1073ce00eea3SSatish Balay x_t = lx[n2 % m]; 107447c6ae99SBarry Smith y_t = ly[(n2 % (m*n))/m]; 107547c6ae99SBarry Smith z_t = lz[n2 / (m*n)]; 107647c6ae99SBarry 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; 1077*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1078ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0 && Zs-zs < 0) { 1079*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 108047c6ae99SBarry Smith } 108147c6ae99SBarry Smith } 108247c6ae99SBarry Smith 108347c6ae99SBarry Smith for (i=0; i<y; i++) { 108447c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 1085ce00eea3SSatish Balay x_t = lx[n3 % m]; 108647c6ae99SBarry Smith y_t = y; 108747c6ae99SBarry Smith z_t = lz[n3 / (m*n)]; 108847c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 1089*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1090ce00eea3SSatish Balay } else if (Xs-xs < 0 && Zs-zs < 0) { 1091*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 109247c6ae99SBarry Smith } 109347c6ae99SBarry Smith 109447c6ae99SBarry Smith if (n4 >= 0) { /* middle */ 109547c6ae99SBarry Smith x_t = x; 109647c6ae99SBarry Smith y_t = y; 109747c6ae99SBarry Smith z_t = lz[n4 / (m*n)]; 109847c6ae99SBarry Smith s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 1099*8865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1100ce00eea3SSatish Balay } else if (Zs-zs < 0) { 1101c2859e5eSBarry Smith if (bz == DMDA_BOUNDARY_MIRROR) { 1102*8865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 1103c2859e5eSBarry Smith } else { 1104*8865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 110547c6ae99SBarry Smith } 1106c2859e5eSBarry Smith } 110747c6ae99SBarry Smith 110847c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 1109ce00eea3SSatish Balay x_t = lx[n5 % m]; 111047c6ae99SBarry Smith y_t = y; 111147c6ae99SBarry Smith z_t = lz[n5 / (m*n)]; 111247c6ae99SBarry Smith s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 1113*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1114ce00eea3SSatish Balay } else if (xe-Xe < 0 && Zs-zs < 0) { 1115*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 111647c6ae99SBarry Smith } 111747c6ae99SBarry Smith } 111847c6ae99SBarry Smith 111947c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 112047c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 1121ce00eea3SSatish Balay x_t = lx[n6 % m]; 112247c6ae99SBarry Smith y_t = ly[(n6 % (m*n))/m]; 112347c6ae99SBarry Smith z_t = lz[n6 / (m*n)]; 112447c6ae99SBarry Smith s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 1125*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1126ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0 && Zs-zs < 0) { 1127*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 112847c6ae99SBarry Smith } 112947c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 113047c6ae99SBarry Smith x_t = x; 113147c6ae99SBarry Smith y_t = ly[(n7 % (m*n))/m]; 113247c6ae99SBarry Smith z_t = lz[n7 / (m*n)]; 113347c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 1134*8865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1135ce00eea3SSatish Balay } else if (ye-Ye < 0 && Zs-zs < 0) { 1136*8865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 113747c6ae99SBarry Smith } 113847c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 1139ce00eea3SSatish Balay x_t = lx[n8 % m]; 114047c6ae99SBarry Smith y_t = ly[(n8 % (m*n))/m]; 114147c6ae99SBarry Smith z_t = lz[n8 / (m*n)]; 114247c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 1143*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1144ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0 && Zs-zs < 0) { 1145*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 114647c6ae99SBarry Smith } 114747c6ae99SBarry Smith } 114847c6ae99SBarry Smith } 114947c6ae99SBarry Smith 115047c6ae99SBarry Smith /* Middle Level */ 115147c6ae99SBarry Smith for (k=0; k<z; k++) { 115247c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 115347c6ae99SBarry Smith if (n9 >= 0) { /* left below */ 1154ce00eea3SSatish Balay x_t = lx[n9 % m]; 115547c6ae99SBarry Smith y_t = ly[(n9 % (m*n))/m]; 115647c6ae99SBarry Smith /* z_t = z; */ 115747c6ae99SBarry Smith s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 1158*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1159ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0) { 1160*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 116147c6ae99SBarry Smith } 116247c6ae99SBarry Smith if (n10 >= 0) { /* directly below */ 116347c6ae99SBarry Smith x_t = x; 116447c6ae99SBarry Smith y_t = ly[(n10 % (m*n))/m]; 116547c6ae99SBarry Smith /* z_t = z; */ 116647c6ae99SBarry Smith s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 1167*8865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1168ce00eea3SSatish Balay } else if (Ys-ys < 0) { 1169c2859e5eSBarry Smith if (by == DMDA_BOUNDARY_MIRROR) { 1170*8865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 1171c2859e5eSBarry Smith } else { 1172*8865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 1173c2859e5eSBarry Smith } 117447c6ae99SBarry Smith } 117547c6ae99SBarry Smith if (n11 >= 0) { /* right below */ 1176ce00eea3SSatish Balay x_t = lx[n11 % m]; 117747c6ae99SBarry Smith y_t = ly[(n11 % (m*n))/m]; 117847c6ae99SBarry Smith /* z_t = z; */ 117947c6ae99SBarry Smith s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 1180*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1181ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0) { 1182*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 118347c6ae99SBarry Smith } 118447c6ae99SBarry Smith } 118547c6ae99SBarry Smith 118647c6ae99SBarry Smith for (i=0; i<y; i++) { 118747c6ae99SBarry Smith if (n12 >= 0) { /* directly left */ 1188ce00eea3SSatish Balay x_t = lx[n12 % m]; 118947c6ae99SBarry Smith y_t = y; 119047c6ae99SBarry Smith /* z_t = z; */ 119147c6ae99SBarry Smith s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t; 1192*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1193ce00eea3SSatish Balay } else if (Xs-xs < 0) { 1194c2859e5eSBarry Smith if (bx == DMDA_BOUNDARY_MIRROR) { 1195*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = 0; 1196c2859e5eSBarry Smith } else { 1197*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 119847c6ae99SBarry Smith } 1199c2859e5eSBarry Smith } 120047c6ae99SBarry Smith 120147c6ae99SBarry Smith /* Interior */ 120247c6ae99SBarry Smith s_t = bases[rank] + i*x + k*x*y; 1203*8865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = s_t++; 120447c6ae99SBarry Smith 120547c6ae99SBarry Smith if (n14 >= 0) { /* directly right */ 1206ce00eea3SSatish Balay x_t = lx[n14 % m]; 120747c6ae99SBarry Smith y_t = y; 120847c6ae99SBarry Smith /* z_t = z; */ 120947c6ae99SBarry Smith s_t = bases[n14] + i*x_t + k*x_t*y_t; 1210*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1211ce00eea3SSatish Balay } else if (xe-Xe < 0) { 1212c2859e5eSBarry Smith if (bx == DMDA_BOUNDARY_MIRROR) { 1213*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = 0; 1214c2859e5eSBarry Smith } else { 1215*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 121647c6ae99SBarry Smith } 121747c6ae99SBarry Smith } 1218c2859e5eSBarry Smith } 121947c6ae99SBarry Smith 122047c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 122147c6ae99SBarry Smith if (n15 >= 0) { /* left above */ 1222ce00eea3SSatish Balay x_t = lx[n15 % m]; 122347c6ae99SBarry Smith y_t = ly[(n15 % (m*n))/m]; 122447c6ae99SBarry Smith /* z_t = z; */ 122547c6ae99SBarry Smith s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t; 1226*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1227ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0) { 1228*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 122947c6ae99SBarry Smith } 123047c6ae99SBarry Smith if (n16 >= 0) { /* directly above */ 123147c6ae99SBarry Smith x_t = x; 123247c6ae99SBarry Smith y_t = ly[(n16 % (m*n))/m]; 123347c6ae99SBarry Smith /* z_t = z; */ 123447c6ae99SBarry Smith s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t; 1235*8865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1236ce00eea3SSatish Balay } else if (ye-Ye < 0) { 1237c2859e5eSBarry Smith if (by == DMDA_BOUNDARY_MIRROR) { 1238*8865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 1239c2859e5eSBarry Smith } else { 1240*8865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 124147c6ae99SBarry Smith } 1242c2859e5eSBarry Smith } 124347c6ae99SBarry Smith if (n17 >= 0) { /* right above */ 1244ce00eea3SSatish Balay x_t = lx[n17 % m]; 124547c6ae99SBarry Smith y_t = ly[(n17 % (m*n))/m]; 124647c6ae99SBarry Smith /* z_t = z; */ 124747c6ae99SBarry Smith s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t; 1248*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1249ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0) { 1250*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 125147c6ae99SBarry Smith } 125247c6ae99SBarry Smith } 125347c6ae99SBarry Smith } 125447c6ae99SBarry Smith 125547c6ae99SBarry Smith /* Upper Level */ 125647c6ae99SBarry Smith for (k=0; k<s_z; k++) { 125747c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 125847c6ae99SBarry Smith if (n18 >= 0) { /* left below */ 1259ce00eea3SSatish Balay x_t = lx[n18 % m]; 126047c6ae99SBarry Smith y_t = ly[(n18 % (m*n))/m]; 126147c6ae99SBarry Smith /* z_t = lz[n18 / (m*n)]; */ 126247c6ae99SBarry Smith s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 1263*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1264ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0 && ze-Ze < 0) { 1265*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 126647c6ae99SBarry Smith } 126747c6ae99SBarry Smith if (n19 >= 0) { /* directly below */ 126847c6ae99SBarry Smith x_t = x; 126947c6ae99SBarry Smith y_t = ly[(n19 % (m*n))/m]; 127047c6ae99SBarry Smith /* z_t = lz[n19 / (m*n)]; */ 127147c6ae99SBarry Smith s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 1272*8865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1273ce00eea3SSatish Balay } else if (Ys-ys < 0 && ze-Ze < 0) { 1274*8865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 127547c6ae99SBarry Smith } 127647c6ae99SBarry Smith if (n20 >= 0) { /* right below */ 1277ce00eea3SSatish Balay x_t = lx[n20 % m]; 127847c6ae99SBarry Smith y_t = ly[(n20 % (m*n))/m]; 127947c6ae99SBarry Smith /* z_t = lz[n20 / (m*n)]; */ 128047c6ae99SBarry Smith s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 1281*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1282ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0 && ze-Ze < 0) { 1283*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 128447c6ae99SBarry Smith } 128547c6ae99SBarry Smith } 128647c6ae99SBarry Smith 128747c6ae99SBarry Smith for (i=0; i<y; i++) { 128847c6ae99SBarry Smith if (n21 >= 0) { /* directly left */ 1289ce00eea3SSatish Balay x_t = lx[n21 % m]; 129047c6ae99SBarry Smith y_t = y; 129147c6ae99SBarry Smith /* z_t = lz[n21 / (m*n)]; */ 129247c6ae99SBarry Smith s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t; 1293*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1294ce00eea3SSatish Balay } else if (Xs-xs < 0 && ze-Ze < 0) { 1295*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 129647c6ae99SBarry Smith } 129747c6ae99SBarry Smith 129847c6ae99SBarry Smith if (n22 >= 0) { /* middle */ 129947c6ae99SBarry Smith x_t = x; 130047c6ae99SBarry Smith y_t = y; 130147c6ae99SBarry Smith /* z_t = lz[n22 / (m*n)]; */ 130247c6ae99SBarry Smith s_t = bases[n22] + i*x_t + k*x_t*y_t; 1303*8865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1304ce00eea3SSatish Balay } else if (ze-Ze < 0) { 1305c2859e5eSBarry Smith if (bz == DMDA_BOUNDARY_MIRROR) { 1306*8865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 1307c2859e5eSBarry Smith } else { 1308*8865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 130947c6ae99SBarry Smith } 1310c2859e5eSBarry Smith } 131147c6ae99SBarry Smith 131247c6ae99SBarry Smith if (n23 >= 0) { /* directly right */ 1313ce00eea3SSatish Balay x_t = lx[n23 % m]; 131447c6ae99SBarry Smith y_t = y; 131547c6ae99SBarry Smith /* z_t = lz[n23 / (m*n)]; */ 131647c6ae99SBarry Smith s_t = bases[n23] + i*x_t + k*x_t*y_t; 1317*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1318ce00eea3SSatish Balay } else if (xe-Xe < 0 && ze-Ze < 0) { 1319*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 132047c6ae99SBarry Smith } 132147c6ae99SBarry Smith } 132247c6ae99SBarry Smith 132347c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 132447c6ae99SBarry Smith if (n24 >= 0) { /* left above */ 1325ce00eea3SSatish Balay x_t = lx[n24 % m]; 132647c6ae99SBarry Smith y_t = ly[(n24 % (m*n))/m]; 132747c6ae99SBarry Smith /* z_t = lz[n24 / (m*n)]; */ 132847c6ae99SBarry Smith s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t; 1329*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1330ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0 && ze-Ze < 0) { 1331*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 133247c6ae99SBarry Smith } 133347c6ae99SBarry Smith if (n25 >= 0) { /* directly above */ 133447c6ae99SBarry Smith x_t = x; 133547c6ae99SBarry Smith y_t = ly[(n25 % (m*n))/m]; 133647c6ae99SBarry Smith /* z_t = lz[n25 / (m*n)]; */ 133747c6ae99SBarry Smith s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t; 1338*8865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1339ce00eea3SSatish Balay } else if (ye-Ye < 0 && ze-Ze < 0) { 1340*8865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 134147c6ae99SBarry Smith } 134247c6ae99SBarry Smith if (n26 >= 0) { /* right above */ 1343ce00eea3SSatish Balay x_t = lx[n26 % m]; 134447c6ae99SBarry Smith y_t = ly[(n26 % (m*n))/m]; 134547c6ae99SBarry Smith /* z_t = lz[n26 / (m*n)]; */ 134647c6ae99SBarry Smith s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t; 1347*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1348ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0 && ze-Ze < 0) { 1349*8865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 135047c6ae99SBarry Smith } 135147c6ae99SBarry Smith } 135247c6ae99SBarry Smith } 135347c6ae99SBarry Smith } 135447c6ae99SBarry Smith /* 135547c6ae99SBarry Smith Set the local to global ordering in the global vector, this allows use 135647c6ae99SBarry Smith of VecSetValuesLocal(). 135747c6ae99SBarry Smith */ 1358ce00eea3SSatish Balay ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_OWN_POINTER,<ogis);CHKERRQ(ierr); 1359ce00eea3SSatish Balay ierr = PetscMalloc(nn*dof*sizeof(PetscInt),&idx_cpy);CHKERRQ(ierr); 1360db87c5ecSEthan Coon ierr = PetscLogObjectMemory(da,nn*dof*sizeof(PetscInt));CHKERRQ(ierr); 13613bf036e2SBarry Smith ierr = ISGetIndices(ltogis, &idx_full);CHKERRQ(ierr); 1362ce00eea3SSatish Balay ierr = PetscMemcpy(idx_cpy,idx_full,nn*dof*sizeof(PetscInt));CHKERRQ(ierr); 13633bf036e2SBarry Smith ierr = ISRestoreIndices(ltogis, &idx_full);CHKERRQ(ierr); 1364ce00eea3SSatish Balay ierr = ISLocalToGlobalMappingCreateIS(ltogis,&da->ltogmap);CHKERRQ(ierr); 1365ce00eea3SSatish Balay ierr = PetscLogObjectParent(da,da->ltogmap);CHKERRQ(ierr); 1366fcfd50ebSBarry Smith ierr = ISDestroy(<ogis);CHKERRQ(ierr); 13671411c6eeSJed Brown ierr = ISLocalToGlobalMappingBlock(da->ltogmap,dd->w,&da->ltogmapb);CHKERRQ(ierr); 13681411c6eeSJed Brown ierr = PetscLogObjectParent(da,da->ltogmap);CHKERRQ(ierr); 136947c6ae99SBarry Smith 1370ce00eea3SSatish Balay ierr = PetscFree2(bases,ldims);CHKERRQ(ierr); 1371ce00eea3SSatish Balay dd->m = m; dd->n = n; dd->p = p; 1372ce00eea3SSatish Balay /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */ 1373ce00eea3SSatish Balay dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = zs; dd->ze = ze; 1374ce00eea3SSatish Balay dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = Zs; dd->Ze = Ze; 1375ce00eea3SSatish Balay 1376fcfd50ebSBarry Smith ierr = VecDestroy(&local);CHKERRQ(ierr); 1377fcfd50ebSBarry Smith ierr = VecDestroy(&global);CHKERRQ(ierr); 1378ce00eea3SSatish Balay 1379ce00eea3SSatish Balay dd->gtol = gtol; 1380ce00eea3SSatish Balay dd->ltog = ltog; 1381ce00eea3SSatish Balay dd->idx = idx_cpy; 1382ce00eea3SSatish Balay dd->Nl = nn*dof; 1383ce00eea3SSatish Balay dd->base = base; 1384ce00eea3SSatish Balay da->ops->view = DMView_DA_3d; 138547c6ae99SBarry Smith dd->ltol = PETSC_NULL; 138647c6ae99SBarry Smith dd->ao = PETSC_NULL; 138747c6ae99SBarry Smith PetscFunctionReturn(0); 138847c6ae99SBarry Smith } 1389564755cdSBarry Smith 139047c6ae99SBarry Smith 139147c6ae99SBarry Smith #undef __FUNCT__ 1392aa219208SBarry Smith #define __FUNCT__ "DMDACreate3d" 139347c6ae99SBarry Smith /*@C 1394aa219208SBarry Smith DMDACreate3d - Creates an object that will manage the communication of three-dimensional 139547c6ae99SBarry Smith regular array data that is distributed across some processors. 139647c6ae99SBarry Smith 139747c6ae99SBarry Smith Collective on MPI_Comm 139847c6ae99SBarry Smith 139947c6ae99SBarry Smith Input Parameters: 140047c6ae99SBarry Smith + comm - MPI communicator 14011321219cSEthan Coon . bx,by,bz - type of ghost nodes the array have. 14021321219cSEthan Coon Use one of DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_PERIODIC. 1403aa219208SBarry Smith . stencil_type - Type of stencil (DMDA_STENCIL_STAR or DMDA_STENCIL_BOX) 140447c6ae99SBarry 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 140547c6ae99SBarry Smith from the command line with -da_grid_x <M> -da_grid_y <N> -da_grid_z <P>) 140647c6ae99SBarry Smith . m,n,p - corresponding number of processors in each dimension 140747c6ae99SBarry Smith (or PETSC_DECIDE to have calculated) 140847c6ae99SBarry Smith . dof - number of degrees of freedom per node 140910d7c030SMatthew G Knepley . s - stencil width 141010d7c030SMatthew G Knepley - lx, ly, lz - arrays containing the number of nodes in each cell along 141147c6ae99SBarry Smith the x, y, and z coordinates, or PETSC_NULL. If non-null, these 141247c6ae99SBarry Smith must be of length as m,n,p and the corresponding 141347c6ae99SBarry Smith m,n, or p cannot be PETSC_DECIDE. Sum of the lx[] entries must be M, sum of 141447c6ae99SBarry Smith the ly[] must N, sum of the lz[] must be P 141547c6ae99SBarry Smith 141647c6ae99SBarry Smith Output Parameter: 141747c6ae99SBarry Smith . da - the resulting distributed array object 141847c6ae99SBarry Smith 141947c6ae99SBarry Smith Options Database Key: 1420706a11cbSBarry Smith + -dm_view - Calls DMView() at the conclusion of DMDACreate3d() 142147c6ae99SBarry Smith . -da_grid_x <nx> - number of grid points in x direction, if M < 0 142247c6ae99SBarry Smith . -da_grid_y <ny> - number of grid points in y direction, if N < 0 142347c6ae99SBarry Smith . -da_grid_z <nz> - number of grid points in z direction, if P < 0 142447c6ae99SBarry Smith . -da_processors_x <MX> - number of processors in x direction 142547c6ae99SBarry Smith . -da_processors_y <MY> - number of processors in y direction 142647c6ae99SBarry Smith . -da_processors_z <MZ> - number of processors in z direction 1427e0f5d30fSBarry Smith . -da_refine_x <rx> - refinement ratio in x direction 1428e0f5d30fSBarry Smith . -da_refine_y <ry> - refinement ratio in y direction 1429e0f5d30fSBarry Smith . -da_refine_z <rz>- refinement ratio in z directio 1430e0f5d30fSBarry Smith - -da_refine <n> - refine the DMDA n times before creating it, , if M, N, or P < 0 143147c6ae99SBarry Smith 143247c6ae99SBarry Smith Level: beginner 143347c6ae99SBarry Smith 143447c6ae99SBarry Smith Notes: 1435aa219208SBarry Smith The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the 1436aa219208SBarry Smith standard 7-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes 143747c6ae99SBarry Smith the standard 27-pt stencil. 143847c6ae99SBarry Smith 1439aa219208SBarry Smith The array data itself is NOT stored in the DMDA, it is stored in Vec objects; 1440564755cdSBarry Smith The appropriate vector objects can be obtained with calls to DMCreateGlobalVector() 1441564755cdSBarry Smith and DMCreateLocalVector() and calls to VecDuplicate() if more are needed. 144247c6ae99SBarry Smith 144347c6ae99SBarry Smith .keywords: distributed array, create, three-dimensional 144447c6ae99SBarry Smith 1445aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate2d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(), 1446aa219208SBarry Smith DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMDALocalToLocalBegin(), DMDALocalToLocalEnd(), DMDASetRefinementFactor(), 1447d461ba97SBarry Smith DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges() 144847c6ae99SBarry Smith 144947c6ae99SBarry Smith @*/ 14501321219cSEthan Coon PetscErrorCode DMDACreate3d(MPI_Comm comm,DMDABoundaryType bx,DMDABoundaryType by,DMDABoundaryType bz,DMDAStencilType stencil_type,PetscInt M, 14519a42bb27SBarry 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) 145247c6ae99SBarry Smith { 145347c6ae99SBarry Smith PetscErrorCode ierr; 145447c6ae99SBarry Smith 145547c6ae99SBarry Smith PetscFunctionBegin; 1456aa219208SBarry Smith ierr = DMDACreate(comm, da);CHKERRQ(ierr); 1457aa219208SBarry Smith ierr = DMDASetDim(*da, 3);CHKERRQ(ierr); 1458aa219208SBarry Smith ierr = DMDASetSizes(*da, M, N, P);CHKERRQ(ierr); 1459aa219208SBarry Smith ierr = DMDASetNumProcs(*da, m, n, p);CHKERRQ(ierr); 14601321219cSEthan Coon ierr = DMDASetBoundaryType(*da, bx, by, bz);CHKERRQ(ierr); 1461aa219208SBarry Smith ierr = DMDASetDof(*da, dof);CHKERRQ(ierr); 1462aa219208SBarry Smith ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr); 1463aa219208SBarry Smith ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr); 1464aa219208SBarry Smith ierr = DMDASetOwnershipRanges(*da, lx, ly, lz);CHKERRQ(ierr); 146547c6ae99SBarry Smith /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */ 14669a42bb27SBarry Smith ierr = DMSetFromOptions(*da);CHKERRQ(ierr); 14679a42bb27SBarry Smith ierr = DMSetUp(*da);CHKERRQ(ierr); 1468ca266f36SBarry Smith ierr = DMViewFromOptions(*da,"-dm_view");CHKERRQ(ierr); 146947c6ae99SBarry Smith PetscFunctionReturn(0); 147047c6ae99SBarry Smith } 1471