147c6ae99SBarry Smith #define PETSCDM_DLL 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 7e1589f56SBarry Smith #include "private/daimpl.h" /*I "petscdm.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 2447c6ae99SBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 2547c6ae99SBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 269a42bb27SBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 279a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 289a42bb27SBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab);CHKERRQ(ierr); 299a42bb27SBarry Smith #endif 3047c6ae99SBarry Smith if (iascii) { 3147c6ae99SBarry Smith PetscViewerFormat format; 3247c6ae99SBarry Smith 3347c6ae99SBarry Smith ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 3447c6ae99SBarry Smith if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL) { 35aa219208SBarry Smith DMDALocalInfo info; 36aa219208SBarry Smith ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr); 3747c6ae99SBarry 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); 3847c6ae99SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"X range of indices: %D %D, Y range of indices: %D %D, Z range of indices: %D %D\n", 3947c6ae99SBarry Smith info.xs,info.xs+info.xm,info.ys,info.ys+info.ym,info.zs,info.zs+info.zm);CHKERRQ(ierr); 4047c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX) 4147c6ae99SBarry Smith if (dd->coordinates) { 4247c6ae99SBarry Smith PetscInt last; 4347c6ae99SBarry Smith const PetscReal *coors; 4447c6ae99SBarry Smith ierr = VecGetArrayRead(dd->coordinates,&coors);CHKERRQ(ierr); 4547c6ae99SBarry Smith ierr = VecGetLocalSize(dd->coordinates,&last);CHKERRQ(ierr); 4647c6ae99SBarry Smith last = last - 3; 4747c6ae99SBarry 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); 4847c6ae99SBarry Smith ierr = VecRestoreArrayRead(dd->coordinates,&coors);CHKERRQ(ierr); 4947c6ae99SBarry Smith } 5047c6ae99SBarry Smith #endif 5147c6ae99SBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 52616157d6SJed Brown } else { 53616157d6SJed Brown ierr = DMView_DA_VTK(da,viewer);CHKERRQ(ierr); 5447c6ae99SBarry Smith } 5547c6ae99SBarry Smith } else if (isdraw) { 5647c6ae99SBarry Smith PetscDraw draw; 5747c6ae99SBarry Smith PetscReal ymin = -1.0,ymax = (PetscReal)dd->N; 5847c6ae99SBarry Smith PetscReal xmin = -1.0,xmax = (PetscReal)((dd->M+2)*dd->P),x,y,ycoord,xcoord; 5947c6ae99SBarry Smith PetscInt k,plane,base,*idx; 6047c6ae99SBarry Smith char node[10]; 6147c6ae99SBarry Smith PetscBool isnull; 6247c6ae99SBarry Smith 6347c6ae99SBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 6447c6ae99SBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 6547c6ae99SBarry Smith ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr); 6647c6ae99SBarry Smith ierr = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr); 6747c6ae99SBarry Smith 6847c6ae99SBarry Smith /* first processor draw all node lines */ 6947c6ae99SBarry Smith if (!rank) { 7047c6ae99SBarry Smith for (k=0; k<dd->P; k++) { 7147c6ae99SBarry Smith ymin = 0.0; ymax = (PetscReal)(dd->N - 1); 7247c6ae99SBarry Smith for (xmin=(PetscReal)(k*(dd->M+1)); xmin<(PetscReal)(dd->M+(k*(dd->M+1))); xmin++) { 7347c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr); 7447c6ae99SBarry Smith } 7547c6ae99SBarry Smith 7647c6ae99SBarry Smith xmin = (PetscReal)(k*(dd->M+1)); xmax = xmin + (PetscReal)(dd->M - 1); 7747c6ae99SBarry Smith for (ymin=0; ymin<(PetscReal)dd->N; ymin++) { 7847c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr); 7947c6ae99SBarry Smith } 8047c6ae99SBarry Smith } 8147c6ae99SBarry Smith } 8247c6ae99SBarry Smith ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr); 8347c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 8447c6ae99SBarry Smith 8547c6ae99SBarry Smith for (k=0; k<dd->P; k++) { /*Go through and draw for each plane*/ 8647c6ae99SBarry Smith if ((k >= dd->zs) && (k < dd->ze)) { 8747c6ae99SBarry Smith /* draw my box */ 8847c6ae99SBarry Smith ymin = dd->ys; 8947c6ae99SBarry Smith ymax = dd->ye - 1; 9047c6ae99SBarry Smith xmin = dd->xs/dd->w + (dd->M+1)*k; 9147c6ae99SBarry Smith xmax =(dd->xe-1)/dd->w + (dd->M+1)*k; 9247c6ae99SBarry Smith 9347c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr); 9447c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 9547c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 9647c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 9747c6ae99SBarry Smith 9847c6ae99SBarry Smith xmin = dd->xs/dd->w; 9947c6ae99SBarry Smith xmax =(dd->xe-1)/dd->w; 10047c6ae99SBarry Smith 10147c6ae99SBarry Smith /* put in numbers*/ 10247c6ae99SBarry Smith base = (dd->base+(dd->xe-dd->xs)*(dd->ye-dd->ys)*(k-dd->zs))/dd->w; 10347c6ae99SBarry Smith 10447c6ae99SBarry Smith /* Identify which processor owns the box */ 10547c6ae99SBarry Smith sprintf(node,"%d",rank); 10647c6ae99SBarry Smith ierr = PetscDrawString(draw,xmin+(dd->M+1)*k+.2,ymin+.3,PETSC_DRAW_RED,node);CHKERRQ(ierr); 10747c6ae99SBarry Smith 10847c6ae99SBarry Smith for (y=ymin; y<=ymax; y++) { 10947c6ae99SBarry Smith for (x=xmin+(dd->M+1)*k; x<=xmax+(dd->M+1)*k; x++) { 11047c6ae99SBarry Smith sprintf(node,"%d",(int)base++); 11147c6ae99SBarry Smith ierr = PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);CHKERRQ(ierr); 11247c6ae99SBarry Smith } 11347c6ae99SBarry Smith } 11447c6ae99SBarry Smith 11547c6ae99SBarry Smith } 11647c6ae99SBarry Smith } 11747c6ae99SBarry Smith ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr); 11847c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 11947c6ae99SBarry Smith 12047c6ae99SBarry Smith for (k=0-dd->s; k<dd->P+dd->s; k++) { 12147c6ae99SBarry Smith /* Go through and draw for each plane */ 12247c6ae99SBarry Smith if ((k >= dd->Zs) && (k < dd->Ze)) { 12347c6ae99SBarry Smith 12447c6ae99SBarry Smith /* overlay ghost numbers, useful for error checking */ 12547c6ae99SBarry Smith base = (dd->Xe-dd->Xs)*(dd->Ye-dd->Ys)*(k-dd->Zs); idx = dd->idx; 12647c6ae99SBarry Smith plane=k; 12747c6ae99SBarry Smith /* Keep z wrap around points on the dradrawg */ 12847c6ae99SBarry Smith if (k<0) { plane=dd->P+k; } 12947c6ae99SBarry Smith if (k>=dd->P) { plane=k-dd->P; } 13047c6ae99SBarry Smith ymin = dd->Ys; ymax = dd->Ye; 13147c6ae99SBarry Smith xmin = (dd->M+1)*plane*dd->w; 13247c6ae99SBarry Smith xmax = (dd->M+1)*plane*dd->w+dd->M*dd->w; 13347c6ae99SBarry Smith for (y=ymin; y<ymax; y++) { 13447c6ae99SBarry Smith for (x=xmin+dd->Xs; x<xmin+dd->Xe; x+=dd->w) { 13547c6ae99SBarry Smith sprintf(node,"%d",(int)(idx[base]/dd->w)); 13647c6ae99SBarry Smith ycoord = y; 13747c6ae99SBarry Smith /*Keep y wrap around points on drawing */ 13847c6ae99SBarry Smith if (y<0) { ycoord = dd->N+y; } 13947c6ae99SBarry Smith 14047c6ae99SBarry Smith if (y>=dd->N) { ycoord = y-dd->N; } 14147c6ae99SBarry Smith xcoord = x; /* Keep x wrap points on drawing */ 14247c6ae99SBarry Smith 14347c6ae99SBarry Smith if (x<xmin) { xcoord = xmax - (xmin-x); } 14447c6ae99SBarry Smith if (x>=xmax) { xcoord = xmin + (x-xmax); } 14547c6ae99SBarry Smith ierr = PetscDrawString(draw,xcoord/dd->w,ycoord,PETSC_DRAW_BLUE,node);CHKERRQ(ierr); 14647c6ae99SBarry Smith base+=dd->w; 14747c6ae99SBarry Smith } 14847c6ae99SBarry Smith } 14947c6ae99SBarry Smith } 15047c6ae99SBarry Smith } 15147c6ae99SBarry Smith ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr); 15247c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 1539a42bb27SBarry Smith } else if (isbinary){ 1549a42bb27SBarry Smith ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr); 1559a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 1569a42bb27SBarry Smith } else if (ismatlab) { 1579a42bb27SBarry Smith ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr); 1589a42bb27SBarry Smith #endif 159aa219208SBarry Smith } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for DMDA 1d",((PetscObject)viewer)->type_name); 16047c6ae99SBarry Smith PetscFunctionReturn(0); 16147c6ae99SBarry Smith } 16247c6ae99SBarry Smith 16347c6ae99SBarry Smith #undef __FUNCT__ 1649a42bb27SBarry Smith #define __FUNCT__ "DMSetUp_DA_3D" 1657087cfbeSBarry Smith PetscErrorCode DMSetUp_DA_3D(DM da) 16647c6ae99SBarry Smith { 16747c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 16847c6ae99SBarry Smith const PetscInt M = dd->M; 16947c6ae99SBarry Smith const PetscInt N = dd->N; 17047c6ae99SBarry Smith const PetscInt P = dd->P; 17147c6ae99SBarry Smith PetscInt m = dd->m; 17247c6ae99SBarry Smith PetscInt n = dd->n; 17347c6ae99SBarry Smith PetscInt p = dd->p; 17447c6ae99SBarry Smith const PetscInt dof = dd->w; 17547c6ae99SBarry Smith const PetscInt s = dd->s; 176aa219208SBarry Smith const DMDAPeriodicType wrap = dd->wrap; 177aa219208SBarry Smith const DMDAStencilType stencil_type = dd->stencil_type; 17847c6ae99SBarry Smith PetscInt *lx = dd->lx; 17947c6ae99SBarry Smith PetscInt *ly = dd->ly; 18047c6ae99SBarry Smith PetscInt *lz = dd->lz; 18147c6ae99SBarry Smith MPI_Comm comm; 18247c6ae99SBarry Smith PetscMPIInt rank,size; 183*ce00eea3SSatish Balay PetscInt xs = 0,xe,ys = 0,ye,zs = 0,ze,x = 0,y = 0,z = 0; 184*ce00eea3SSatish Balay PetscInt Xs,Xe,Ys,Ye,Zs,Ze,IXs,IXe,IYs,IYe,IZs,IZe,start,end,pm; 185*ce00eea3SSatish Balay PetscInt left,right,up,down,bottom,top,i,j,k,*idx,*idx_cpy,nn; 186*ce00eea3SSatish Balay const PetscInt *idx_full; 18747c6ae99SBarry Smith PetscInt n0,n1,n2,n3,n4,n5,n6,n7,n8,n9,n10,n11,n12,n14; 18847c6ae99SBarry Smith PetscInt n15,n16,n17,n18,n19,n20,n21,n22,n23,n24,n25,n26; 189*ce00eea3SSatish Balay PetscInt *bases,*ldims,base,x_t,y_t,z_t,s_t,count,count_dbg,s_x,s_y,s_z; 19047c6ae99SBarry Smith PetscInt sn0 = 0,sn1 = 0,sn2 = 0,sn3 = 0,sn5 = 0,sn6 = 0,sn7 = 0; 19147c6ae99SBarry Smith PetscInt sn8 = 0,sn9 = 0,sn11 = 0,sn15 = 0,sn24 = 0,sn25 = 0,sn26 = 0; 19247c6ae99SBarry Smith PetscInt sn17 = 0,sn18 = 0,sn19 = 0,sn20 = 0,sn21 = 0,sn23 = 0; 19347c6ae99SBarry Smith Vec local,global; 19447c6ae99SBarry Smith VecScatter ltog,gtol; 195*ce00eea3SSatish Balay IS to,from,ltogis; 19647c6ae99SBarry Smith PetscErrorCode ierr; 19747c6ae99SBarry Smith 19847c6ae99SBarry Smith PetscFunctionBegin; 19947c6ae99SBarry Smith if (dof < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof); 20047c6ae99SBarry Smith if (s < 0) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s); 20147c6ae99SBarry Smith 20247c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject) da, &comm);CHKERRQ(ierr); 20347c6ae99SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 20447c6ae99SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 20547c6ae99SBarry Smith 20647c6ae99SBarry Smith dd->dim = 3; 20747c6ae99SBarry Smith ierr = PetscMalloc(dof*sizeof(char*),&dd->fieldname);CHKERRQ(ierr); 20847c6ae99SBarry Smith ierr = PetscMemzero(dd->fieldname,dof*sizeof(char*));CHKERRQ(ierr); 20947c6ae99SBarry Smith 21047c6ae99SBarry Smith if (m != PETSC_DECIDE) { 21147c6ae99SBarry Smith if (m < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m); 21247c6ae99SBarry Smith else if (m > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size); 21347c6ae99SBarry Smith } 21447c6ae99SBarry Smith if (n != PETSC_DECIDE) { 21547c6ae99SBarry Smith if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n); 21647c6ae99SBarry Smith else if (n > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size); 21747c6ae99SBarry Smith } 21847c6ae99SBarry Smith if (p != PETSC_DECIDE) { 21947c6ae99SBarry Smith if (p < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Z direction: %D",p); 22047c6ae99SBarry Smith else if (p > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Z direction: %D %d",p,size); 22147c6ae99SBarry Smith } 22247c6ae99SBarry Smith 22347c6ae99SBarry Smith /* Partition the array among the processors */ 22447c6ae99SBarry Smith if (m == PETSC_DECIDE && n != PETSC_DECIDE && p != PETSC_DECIDE) { 22547c6ae99SBarry Smith m = size/(n*p); 22647c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) { 22747c6ae99SBarry Smith n = size/(m*p); 22847c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) { 22947c6ae99SBarry Smith p = size/(m*n); 23047c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) { 23147c6ae99SBarry Smith /* try for squarish distribution */ 23247c6ae99SBarry Smith m = (int)(0.5 + sqrt(((PetscReal)M)*((PetscReal)size)/((PetscReal)N*p))); 23347c6ae99SBarry Smith if (!m) m = 1; 23447c6ae99SBarry Smith while (m > 0) { 23547c6ae99SBarry Smith n = size/(m*p); 23647c6ae99SBarry Smith if (m*n*p == size) break; 23747c6ae99SBarry Smith m--; 23847c6ae99SBarry Smith } 23947c6ae99SBarry Smith if (!m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad p value: p = %D",p); 24047c6ae99SBarry Smith if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;} 24147c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) { 24247c6ae99SBarry Smith /* try for squarish distribution */ 24347c6ae99SBarry Smith m = (int)(0.5 + sqrt(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n))); 24447c6ae99SBarry Smith if (!m) m = 1; 24547c6ae99SBarry Smith while (m > 0) { 24647c6ae99SBarry Smith p = size/(m*n); 24747c6ae99SBarry Smith if (m*n*p == size) break; 24847c6ae99SBarry Smith m--; 24947c6ae99SBarry Smith } 25047c6ae99SBarry Smith if (!m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad n value: n = %D",n); 25147c6ae99SBarry Smith if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;} 25247c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) { 25347c6ae99SBarry Smith /* try for squarish distribution */ 25447c6ae99SBarry Smith n = (int)(0.5 + sqrt(((PetscReal)N)*((PetscReal)size)/((PetscReal)P*m))); 25547c6ae99SBarry Smith if (!n) n = 1; 25647c6ae99SBarry Smith while (n > 0) { 25747c6ae99SBarry Smith p = size/(m*n); 25847c6ae99SBarry Smith if (m*n*p == size) break; 25947c6ae99SBarry Smith n--; 26047c6ae99SBarry Smith } 26147c6ae99SBarry Smith if (!n) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad m value: m = %D",n); 26247c6ae99SBarry Smith if (N > P && n < p) {PetscInt _n = n; n = p; p = _n;} 26347c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) { 26447c6ae99SBarry Smith /* try for squarish distribution */ 26547c6ae99SBarry Smith n = (PetscInt)(0.5 + pow(((PetscReal)N*N)*((PetscReal)size)/((PetscReal)P*M),(PetscReal)(1./3.))); 26647c6ae99SBarry Smith if (!n) n = 1; 26747c6ae99SBarry Smith while (n > 0) { 26847c6ae99SBarry Smith pm = size/n; 26947c6ae99SBarry Smith if (n*pm == size) break; 27047c6ae99SBarry Smith n--; 27147c6ae99SBarry Smith } 27247c6ae99SBarry Smith if (!n) n = 1; 27347c6ae99SBarry Smith m = (PetscInt)(0.5 + sqrt(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n))); 27447c6ae99SBarry Smith if (!m) m = 1; 27547c6ae99SBarry Smith while (m > 0) { 27647c6ae99SBarry Smith p = size/(m*n); 27747c6ae99SBarry Smith if (m*n*p == size) break; 27847c6ae99SBarry Smith m--; 27947c6ae99SBarry Smith } 28047c6ae99SBarry Smith if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;} 28147c6ae99SBarry Smith } else if (m*n*p != size) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition"); 28247c6ae99SBarry Smith 28347c6ae99SBarry Smith if (m*n*p != size) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_PLIB,"Could not find good partition"); 28447c6ae99SBarry Smith if (M < m) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m); 28547c6ae99SBarry Smith if (N < n) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n); 28647c6ae99SBarry Smith if (P < p) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in z direction is too fine! %D %D",P,p); 28747c6ae99SBarry Smith 28847c6ae99SBarry Smith /* 28947c6ae99SBarry Smith Determine locally owned region 29047c6ae99SBarry Smith [x, y, or z]s is the first local node number, [x, y, z] is the number of local nodes 29147c6ae99SBarry Smith */ 29247c6ae99SBarry Smith 29347c6ae99SBarry Smith if (!lx) { 29447c6ae99SBarry Smith ierr = PetscMalloc(m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr); 29547c6ae99SBarry Smith lx = dd->lx; 29647c6ae99SBarry Smith for (i=0; i<m; i++) { 29747c6ae99SBarry Smith lx[i] = M/m + ((M % m) > (i % m)); 29847c6ae99SBarry Smith } 29947c6ae99SBarry Smith } 30047c6ae99SBarry Smith x = lx[rank % m]; 30147c6ae99SBarry Smith xs = 0; 30247c6ae99SBarry Smith for (i=0; i<(rank%m); i++) { xs += lx[i];} 30347c6ae99SBarry Smith if (m > 1 && x < s) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column width is too thin for stencil! %D %D",x,s); 30447c6ae99SBarry Smith 30547c6ae99SBarry Smith if (!ly) { 30647c6ae99SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr); 30747c6ae99SBarry Smith ly = dd->ly; 30847c6ae99SBarry Smith for (i=0; i<n; i++) { 30947c6ae99SBarry Smith ly[i] = N/n + ((N % n) > (i % n)); 31047c6ae99SBarry Smith } 31147c6ae99SBarry Smith } 31247c6ae99SBarry Smith y = ly[(rank % (m*n))/m]; 31347c6ae99SBarry Smith if (n > 1 && y < s) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row width is too thin for stencil! %D %D",y,s); 31447c6ae99SBarry Smith ys = 0; 31547c6ae99SBarry Smith for (i=0; i<(rank % (m*n))/m; i++) { ys += ly[i];} 31647c6ae99SBarry Smith 31747c6ae99SBarry Smith if (!lz) { 31847c6ae99SBarry Smith ierr = PetscMalloc(p*sizeof(PetscInt), &dd->lz);CHKERRQ(ierr); 31947c6ae99SBarry Smith lz = dd->lz; 32047c6ae99SBarry Smith for (i=0; i<p; i++) { 32147c6ae99SBarry Smith lz[i] = P/p + ((P % p) > (i % p)); 32247c6ae99SBarry Smith } 32347c6ae99SBarry Smith } 32447c6ae99SBarry Smith z = lz[rank/(m*n)]; 32547c6ae99SBarry Smith if (p > 1 && z < s) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Plane width is too thin for stencil! %D %D",z,s); 32647c6ae99SBarry Smith zs = 0; 32747c6ae99SBarry Smith for (i=0; i<(rank/(m*n)); i++) { zs += lz[i];} 32847c6ae99SBarry Smith ye = ys + y; 32947c6ae99SBarry Smith xe = xs + x; 33047c6ae99SBarry Smith ze = zs + z; 33147c6ae99SBarry Smith 33247c6ae99SBarry Smith /* determine ghost region */ 33347c6ae99SBarry Smith /* Assume No Periodicity */ 334*ce00eea3SSatish Balay if (xs-s > 0) { Xs = xs - s; IXs = xs - s; } else { Xs = 0; IXs = 0; } 335*ce00eea3SSatish Balay if (xe+s <= M) { Xe = xe + s; IXe = xe + s; } else { Xe = M; IXe = M; } 336*ce00eea3SSatish Balay if (ys-s > 0) { Ys = ys - s; IYs = ys - s; } else { Ys = 0; IYs = 0; } 337*ce00eea3SSatish Balay if (ye+s <= N) { Ye = ye + s; IYe = ye + s; } else { Ye = N; IYe = N; } 338*ce00eea3SSatish Balay if (zs-s > 0) { Zs = zs - s; IZs = zs - s; } else { Zs = 0; IZs = 0; } 339*ce00eea3SSatish Balay if (ze+s <= P) { Ze = ze + s; IZe = ze + s; } else { Ze = P; IZe = P; } 34047c6ae99SBarry Smith 341*ce00eea3SSatish Balay /* fix for periodicity/ghosted */ 342*ce00eea3SSatish Balay if (DMDAXGhosted(wrap)) { Xs = xs - s; Xe = xe + s; } 343*ce00eea3SSatish Balay if (DMDAXPeriodic(wrap)) { IXs = xs - s; IXe = xe + s; } 344*ce00eea3SSatish Balay if (DMDAYGhosted(wrap)) { Ys = ys - s; Ye = ye + s; } 345*ce00eea3SSatish Balay if (DMDAYPeriodic(wrap)) { IYs = ys - s; IYe = ye + s; } 346*ce00eea3SSatish Balay if (DMDAZGhosted(wrap)) { Zs = zs - s; Ze = ze + s; } 347*ce00eea3SSatish Balay if (DMDAZPeriodic(wrap)) { IZs = zs - s; IZe = ze + s; } 34847c6ae99SBarry Smith 34947c6ae99SBarry Smith /* Resize all X parameters to reflect w */ 350*ce00eea3SSatish Balay s_x = s; 35147c6ae99SBarry Smith s_y = s; 35247c6ae99SBarry Smith s_z = s; 35347c6ae99SBarry Smith 35447c6ae99SBarry Smith /* determine starting point of each processor */ 35547c6ae99SBarry Smith nn = x*y*z; 35647c6ae99SBarry Smith ierr = PetscMalloc2(size+1,PetscInt,&bases,size,PetscInt,&ldims);CHKERRQ(ierr); 35747c6ae99SBarry Smith ierr = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr); 35847c6ae99SBarry Smith bases[0] = 0; 35947c6ae99SBarry Smith for (i=1; i<=size; i++) { 36047c6ae99SBarry Smith bases[i] = ldims[i-1]; 36147c6ae99SBarry Smith } 36247c6ae99SBarry Smith for (i=1; i<=size; i++) { 36347c6ae99SBarry Smith bases[i] += bases[i-1]; 36447c6ae99SBarry Smith } 365*ce00eea3SSatish Balay base = bases[rank]*dof; 36647c6ae99SBarry Smith 36747c6ae99SBarry Smith /* allocate the base parallel and sequential vectors */ 368*ce00eea3SSatish Balay dd->Nlocal = x*y*z*dof; 36947c6ae99SBarry Smith ierr = VecCreateMPIWithArray(comm,dd->Nlocal,PETSC_DECIDE,0,&global);CHKERRQ(ierr); 37047c6ae99SBarry Smith ierr = VecSetBlockSize(global,dof);CHKERRQ(ierr); 371*ce00eea3SSatish Balay dd->nlocal = (Xe-Xs)*(Ye-Ys)*(Ze-Zs)*dof; 372*ce00eea3SSatish Balay ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dd->nlocal,0,&local);CHKERRQ(ierr); 37347c6ae99SBarry Smith ierr = VecSetBlockSize(local,dof);CHKERRQ(ierr); 37447c6ae99SBarry Smith 37547c6ae99SBarry Smith /* generate appropriate vector scatters */ 37647c6ae99SBarry Smith /* local to global inserts non-ghost point region into global */ 37747c6ae99SBarry Smith ierr = VecGetOwnershipRange(global,&start,&end);CHKERRQ(ierr); 378*ce00eea3SSatish Balay ierr = ISCreateStride(comm,x*y*z*dof,start,1,&to);CHKERRQ(ierr); 37947c6ae99SBarry Smith 380*ce00eea3SSatish Balay count_dbg = x*y*z; 381*ce00eea3SSatish Balay ierr = PetscMalloc(x*y*z*sizeof(PetscInt),&idx);CHKERRQ(ierr); 382*ce00eea3SSatish Balay left = xs - Xs; right = left + x; 38347c6ae99SBarry Smith bottom = ys - Ys; top = bottom + y; 38447c6ae99SBarry Smith down = zs - Zs; up = down + z; 38547c6ae99SBarry Smith count = 0; 38647c6ae99SBarry Smith for (i=down; i<up; i++) { 38747c6ae99SBarry Smith for (j=bottom; j<top; j++) { 388*ce00eea3SSatish Balay for (k=left; k<right; k++) { 389*ce00eea3SSatish Balay idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 39047c6ae99SBarry Smith } 39147c6ae99SBarry Smith } 39247c6ae99SBarry Smith } 39347c6ae99SBarry Smith 394*ce00eea3SSatish Balay if (count != count_dbg) { 395*ce00eea3SSatish Balay SETERRQ(((PetscObject)da)->comm, PETSC_ERR_SUP,"count != count_dbg"); 396*ce00eea3SSatish Balay PetscFunctionReturn(1); 397*ce00eea3SSatish Balay } 39847c6ae99SBarry Smith ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&from);CHKERRQ(ierr); 39947c6ae99SBarry Smith ierr = VecScatterCreate(local,from,global,to,<og);CHKERRQ(ierr); 40047c6ae99SBarry Smith ierr = PetscLogObjectParent(da,ltog);CHKERRQ(ierr); 40147c6ae99SBarry Smith ierr = ISDestroy(from);CHKERRQ(ierr); 40247c6ae99SBarry Smith ierr = ISDestroy(to);CHKERRQ(ierr); 40347c6ae99SBarry Smith 404*ce00eea3SSatish Balay /* global to local must include ghost points within the domain, 405*ce00eea3SSatish Balay but not ghost points outside the domain that aren't periodic */ 406aa219208SBarry Smith if (stencil_type == DMDA_STENCIL_BOX) { 407*ce00eea3SSatish Balay count_dbg = (IXe-IXs)*(IYe-IYs)*(IZe-IZs); 408*ce00eea3SSatish Balay ierr = PetscMalloc(count_dbg*sizeof(PetscInt),&idx);CHKERRQ(ierr); 409*ce00eea3SSatish Balay 410*ce00eea3SSatish Balay left = IXs - Xs; right = left + (IXe-IXs); 411*ce00eea3SSatish Balay bottom = IYs - Ys; top = bottom + (IYe-IYs); 412*ce00eea3SSatish Balay down = IZs - Zs; up = down + (IZe-IZs); 413*ce00eea3SSatish Balay count = 0; 414*ce00eea3SSatish Balay for (i=down; i<up; i++) { 415*ce00eea3SSatish Balay for (j=bottom; j<top; j++) { 416*ce00eea3SSatish Balay for (k=left; k<right; k++) { 417*ce00eea3SSatish Balay idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 418*ce00eea3SSatish Balay } 419*ce00eea3SSatish Balay } 420*ce00eea3SSatish Balay } 421*ce00eea3SSatish Balay if (count != count_dbg) { 422*ce00eea3SSatish Balay SETERRQ(((PetscObject)da)->comm, PETSC_ERR_SUP,"count != count_dbg"); 423*ce00eea3SSatish Balay PetscFunctionReturn(1); 424*ce00eea3SSatish Balay } 425*ce00eea3SSatish Balay ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr); 426*ce00eea3SSatish Balay 42747c6ae99SBarry Smith } else { 42847c6ae99SBarry Smith /* This is way ugly! We need to list the funny cross type region */ 429*ce00eea3SSatish Balay count_dbg = ((ys-IYs) + (IYe-ye))*x*z + ((xs-IXs) + (IXe-xe))*y*z + ((zs-IZs) + (IZe-ze))*x*y + x*y*z; 430*ce00eea3SSatish Balay ierr = PetscMalloc(count_dbg*sizeof(PetscInt),&idx);CHKERRQ(ierr); 431*ce00eea3SSatish Balay 432*ce00eea3SSatish 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; 436*ce00eea3SSatish Balay /* the bottom chunck */ 437*ce00eea3SSatish Balay for (i=(IZs-Zs); i<down; i++) { 43847c6ae99SBarry Smith for (j=bottom; j<top; j++) { 439*ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 44047c6ae99SBarry Smith } 44147c6ae99SBarry Smith } 44247c6ae99SBarry Smith /* the middle piece */ 44347c6ae99SBarry Smith for (i=down; i<up; i++) { 44447c6ae99SBarry Smith /* front */ 445*ce00eea3SSatish Balay for (j=(IYs-Ys); j<bottom; j++) { 446*ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 44747c6ae99SBarry Smith } 44847c6ae99SBarry Smith /* middle */ 44947c6ae99SBarry Smith for (j=bottom; j<top; j++) { 450*ce00eea3SSatish Balay for (k=IXs-Xs; k<IXe-Xs; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 45147c6ae99SBarry Smith } 45247c6ae99SBarry Smith /* back */ 453*ce00eea3SSatish Balay for (j=top; j<top+IYe-ye; j++) { 454*ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 45547c6ae99SBarry Smith } 45647c6ae99SBarry Smith } 45747c6ae99SBarry Smith /* the top piece */ 458*ce00eea3SSatish Balay for (i=up; i<up+IZe-ze; i++) { 45947c6ae99SBarry Smith for (j=bottom; j<top; j++) { 460*ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 46147c6ae99SBarry Smith } 46247c6ae99SBarry Smith } 463*ce00eea3SSatish Balay if (count != count_dbg) { 464*ce00eea3SSatish Balay SETERRQ(((PetscObject)da)->comm, PETSC_ERR_SUP,"count != count_dbg"); 465*ce00eea3SSatish Balay PetscFunctionReturn(1); 466*ce00eea3SSatish Balay } 46747c6ae99SBarry Smith ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr); 46847c6ae99SBarry Smith } 46947c6ae99SBarry Smith 47047c6ae99SBarry Smith /* determine who lies on each side of use stored in n24 n25 n26 47147c6ae99SBarry Smith n21 n22 n23 47247c6ae99SBarry Smith n18 n19 n20 47347c6ae99SBarry Smith 47447c6ae99SBarry Smith n15 n16 n17 47547c6ae99SBarry Smith n12 n14 47647c6ae99SBarry Smith n9 n10 n11 47747c6ae99SBarry Smith 47847c6ae99SBarry Smith n6 n7 n8 47947c6ae99SBarry Smith n3 n4 n5 48047c6ae99SBarry Smith n0 n1 n2 48147c6ae99SBarry Smith */ 48247c6ae99SBarry Smith 48347c6ae99SBarry Smith /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */ 48447c6ae99SBarry Smith /* Assume Nodes are Internal to the Cube */ 48547c6ae99SBarry Smith n0 = rank - m*n - m - 1; 48647c6ae99SBarry Smith n1 = rank - m*n - m; 48747c6ae99SBarry Smith n2 = rank - m*n - m + 1; 48847c6ae99SBarry Smith n3 = rank - m*n -1; 48947c6ae99SBarry Smith n4 = rank - m*n; 49047c6ae99SBarry Smith n5 = rank - m*n + 1; 49147c6ae99SBarry Smith n6 = rank - m*n + m - 1; 49247c6ae99SBarry Smith n7 = rank - m*n + m; 49347c6ae99SBarry Smith n8 = rank - m*n + m + 1; 49447c6ae99SBarry Smith 49547c6ae99SBarry Smith n9 = rank - m - 1; 49647c6ae99SBarry Smith n10 = rank - m; 49747c6ae99SBarry Smith n11 = rank - m + 1; 49847c6ae99SBarry Smith n12 = rank - 1; 49947c6ae99SBarry Smith n14 = rank + 1; 50047c6ae99SBarry Smith n15 = rank + m - 1; 50147c6ae99SBarry Smith n16 = rank + m; 50247c6ae99SBarry Smith n17 = rank + m + 1; 50347c6ae99SBarry Smith 50447c6ae99SBarry Smith n18 = rank + m*n - m - 1; 50547c6ae99SBarry Smith n19 = rank + m*n - m; 50647c6ae99SBarry Smith n20 = rank + m*n - m + 1; 50747c6ae99SBarry Smith n21 = rank + m*n - 1; 50847c6ae99SBarry Smith n22 = rank + m*n; 50947c6ae99SBarry Smith n23 = rank + m*n + 1; 51047c6ae99SBarry Smith n24 = rank + m*n + m - 1; 51147c6ae99SBarry Smith n25 = rank + m*n + m; 51247c6ae99SBarry Smith n26 = rank + m*n + m + 1; 51347c6ae99SBarry Smith 51447c6ae99SBarry Smith /* Assume Pieces are on Faces of Cube */ 51547c6ae99SBarry Smith 51647c6ae99SBarry Smith if (xs == 0) { /* First assume not corner or edge */ 51747c6ae99SBarry Smith n0 = rank -1 - (m*n); 51847c6ae99SBarry Smith n3 = rank + m -1 - (m*n); 51947c6ae99SBarry Smith n6 = rank + 2*m -1 - (m*n); 52047c6ae99SBarry Smith n9 = rank -1; 52147c6ae99SBarry Smith n12 = rank + m -1; 52247c6ae99SBarry Smith n15 = rank + 2*m -1; 52347c6ae99SBarry Smith n18 = rank -1 + (m*n); 52447c6ae99SBarry Smith n21 = rank + m -1 + (m*n); 52547c6ae99SBarry Smith n24 = rank + 2*m -1 + (m*n); 52647c6ae99SBarry Smith } 52747c6ae99SBarry Smith 528*ce00eea3SSatish Balay if (xe == M) { /* First assume not corner or edge */ 52947c6ae99SBarry Smith n2 = rank -2*m +1 - (m*n); 53047c6ae99SBarry Smith n5 = rank - m +1 - (m*n); 53147c6ae99SBarry Smith n8 = rank +1 - (m*n); 53247c6ae99SBarry Smith n11 = rank -2*m +1; 53347c6ae99SBarry Smith n14 = rank - m +1; 53447c6ae99SBarry Smith n17 = rank +1; 53547c6ae99SBarry Smith n20 = rank -2*m +1 + (m*n); 53647c6ae99SBarry Smith n23 = rank - m +1 + (m*n); 53747c6ae99SBarry Smith n26 = rank +1 + (m*n); 53847c6ae99SBarry Smith } 53947c6ae99SBarry Smith 54047c6ae99SBarry Smith if (ys==0) { /* First assume not corner or edge */ 54147c6ae99SBarry Smith n0 = rank + m * (n-1) -1 - (m*n); 54247c6ae99SBarry Smith n1 = rank + m * (n-1) - (m*n); 54347c6ae99SBarry Smith n2 = rank + m * (n-1) +1 - (m*n); 54447c6ae99SBarry Smith n9 = rank + m * (n-1) -1; 54547c6ae99SBarry Smith n10 = rank + m * (n-1); 54647c6ae99SBarry Smith n11 = rank + m * (n-1) +1; 54747c6ae99SBarry Smith n18 = rank + m * (n-1) -1 + (m*n); 54847c6ae99SBarry Smith n19 = rank + m * (n-1) + (m*n); 54947c6ae99SBarry Smith n20 = rank + m * (n-1) +1 + (m*n); 55047c6ae99SBarry Smith } 55147c6ae99SBarry Smith 55247c6ae99SBarry Smith if (ye == N) { /* First assume not corner or edge */ 55347c6ae99SBarry Smith n6 = rank - m * (n-1) -1 - (m*n); 55447c6ae99SBarry Smith n7 = rank - m * (n-1) - (m*n); 55547c6ae99SBarry Smith n8 = rank - m * (n-1) +1 - (m*n); 55647c6ae99SBarry Smith n15 = rank - m * (n-1) -1; 55747c6ae99SBarry Smith n16 = rank - m * (n-1); 55847c6ae99SBarry Smith n17 = rank - m * (n-1) +1; 55947c6ae99SBarry Smith n24 = rank - m * (n-1) -1 + (m*n); 56047c6ae99SBarry Smith n25 = rank - m * (n-1) + (m*n); 56147c6ae99SBarry Smith n26 = rank - m * (n-1) +1 + (m*n); 56247c6ae99SBarry Smith } 56347c6ae99SBarry Smith 56447c6ae99SBarry Smith if (zs == 0) { /* First assume not corner or edge */ 56547c6ae99SBarry Smith n0 = size - (m*n) + rank - m - 1; 56647c6ae99SBarry Smith n1 = size - (m*n) + rank - m; 56747c6ae99SBarry Smith n2 = size - (m*n) + rank - m + 1; 56847c6ae99SBarry Smith n3 = size - (m*n) + rank - 1; 56947c6ae99SBarry Smith n4 = size - (m*n) + rank; 57047c6ae99SBarry Smith n5 = size - (m*n) + rank + 1; 57147c6ae99SBarry Smith n6 = size - (m*n) + rank + m - 1; 57247c6ae99SBarry Smith n7 = size - (m*n) + rank + m ; 57347c6ae99SBarry Smith n8 = size - (m*n) + rank + m + 1; 57447c6ae99SBarry Smith } 57547c6ae99SBarry Smith 57647c6ae99SBarry Smith if (ze == P) { /* First assume not corner or edge */ 57747c6ae99SBarry Smith n18 = (m*n) - (size-rank) - m - 1; 57847c6ae99SBarry Smith n19 = (m*n) - (size-rank) - m; 57947c6ae99SBarry Smith n20 = (m*n) - (size-rank) - m + 1; 58047c6ae99SBarry Smith n21 = (m*n) - (size-rank) - 1; 58147c6ae99SBarry Smith n22 = (m*n) - (size-rank); 58247c6ae99SBarry Smith n23 = (m*n) - (size-rank) + 1; 58347c6ae99SBarry Smith n24 = (m*n) - (size-rank) + m - 1; 58447c6ae99SBarry Smith n25 = (m*n) - (size-rank) + m; 58547c6ae99SBarry Smith n26 = (m*n) - (size-rank) + m + 1; 58647c6ae99SBarry Smith } 58747c6ae99SBarry Smith 58847c6ae99SBarry Smith if ((xs==0) && (zs==0)) { /* Assume an edge, not corner */ 58947c6ae99SBarry Smith n0 = size - m*n + rank + m-1 - m; 59047c6ae99SBarry Smith n3 = size - m*n + rank + m-1; 59147c6ae99SBarry Smith n6 = size - m*n + rank + m-1 + m; 59247c6ae99SBarry Smith } 59347c6ae99SBarry Smith 59447c6ae99SBarry Smith if ((xs==0) && (ze==P)) { /* Assume an edge, not corner */ 59547c6ae99SBarry Smith n18 = m*n - (size - rank) + m-1 - m; 59647c6ae99SBarry Smith n21 = m*n - (size - rank) + m-1; 59747c6ae99SBarry Smith n24 = m*n - (size - rank) + m-1 + m; 59847c6ae99SBarry Smith } 59947c6ae99SBarry Smith 60047c6ae99SBarry Smith if ((xs==0) && (ys==0)) { /* Assume an edge, not corner */ 60147c6ae99SBarry Smith n0 = rank + m*n -1 - m*n; 60247c6ae99SBarry Smith n9 = rank + m*n -1; 60347c6ae99SBarry Smith n18 = rank + m*n -1 + m*n; 60447c6ae99SBarry Smith } 60547c6ae99SBarry Smith 60647c6ae99SBarry Smith if ((xs==0) && (ye==N)) { /* Assume an edge, not corner */ 60747c6ae99SBarry Smith n6 = rank - m*(n-1) + m-1 - m*n; 60847c6ae99SBarry Smith n15 = rank - m*(n-1) + m-1; 60947c6ae99SBarry Smith n24 = rank - m*(n-1) + m-1 + m*n; 61047c6ae99SBarry Smith } 61147c6ae99SBarry Smith 612*ce00eea3SSatish Balay if ((xe==M) && (zs==0)) { /* Assume an edge, not corner */ 61347c6ae99SBarry Smith n2 = size - (m*n-rank) - (m-1) - m; 61447c6ae99SBarry Smith n5 = size - (m*n-rank) - (m-1); 61547c6ae99SBarry Smith n8 = size - (m*n-rank) - (m-1) + m; 61647c6ae99SBarry Smith } 61747c6ae99SBarry Smith 618*ce00eea3SSatish Balay if ((xe==M) && (ze==P)) { /* Assume an edge, not corner */ 61947c6ae99SBarry Smith n20 = m*n - (size - rank) - (m-1) - m; 62047c6ae99SBarry Smith n23 = m*n - (size - rank) - (m-1); 62147c6ae99SBarry Smith n26 = m*n - (size - rank) - (m-1) + m; 62247c6ae99SBarry Smith } 62347c6ae99SBarry Smith 624*ce00eea3SSatish Balay if ((xe==M) && (ys==0)) { /* Assume an edge, not corner */ 62547c6ae99SBarry Smith n2 = rank + m*(n-1) - (m-1) - m*n; 62647c6ae99SBarry Smith n11 = rank + m*(n-1) - (m-1); 62747c6ae99SBarry Smith n20 = rank + m*(n-1) - (m-1) + m*n; 62847c6ae99SBarry Smith } 62947c6ae99SBarry Smith 630*ce00eea3SSatish Balay if ((xe==M) && (ye==N)) { /* Assume an edge, not corner */ 63147c6ae99SBarry Smith n8 = rank - m*n +1 - m*n; 63247c6ae99SBarry Smith n17 = rank - m*n +1; 63347c6ae99SBarry Smith n26 = rank - m*n +1 + m*n; 63447c6ae99SBarry Smith } 63547c6ae99SBarry Smith 63647c6ae99SBarry Smith if ((ys==0) && (zs==0)) { /* Assume an edge, not corner */ 63747c6ae99SBarry Smith n0 = size - m + rank -1; 63847c6ae99SBarry Smith n1 = size - m + rank; 63947c6ae99SBarry Smith n2 = size - m + rank +1; 64047c6ae99SBarry Smith } 64147c6ae99SBarry Smith 64247c6ae99SBarry Smith if ((ys==0) && (ze==P)) { /* Assume an edge, not corner */ 64347c6ae99SBarry Smith n18 = m*n - (size - rank) + m*(n-1) -1; 64447c6ae99SBarry Smith n19 = m*n - (size - rank) + m*(n-1); 64547c6ae99SBarry Smith n20 = m*n - (size - rank) + m*(n-1) +1; 64647c6ae99SBarry Smith } 64747c6ae99SBarry Smith 64847c6ae99SBarry Smith if ((ye==N) && (zs==0)) { /* Assume an edge, not corner */ 64947c6ae99SBarry Smith n6 = size - (m*n-rank) - m * (n-1) -1; 65047c6ae99SBarry Smith n7 = size - (m*n-rank) - m * (n-1); 65147c6ae99SBarry Smith n8 = size - (m*n-rank) - m * (n-1) +1; 65247c6ae99SBarry Smith } 65347c6ae99SBarry Smith 65447c6ae99SBarry Smith if ((ye==N) && (ze==P)) { /* Assume an edge, not corner */ 65547c6ae99SBarry Smith n24 = rank - (size-m) -1; 65647c6ae99SBarry Smith n25 = rank - (size-m); 65747c6ae99SBarry Smith n26 = rank - (size-m) +1; 65847c6ae99SBarry Smith } 65947c6ae99SBarry Smith 66047c6ae99SBarry Smith /* Check for Corners */ 66147c6ae99SBarry Smith if ((xs==0) && (ys==0) && (zs==0)) { n0 = size -1;} 66247c6ae99SBarry Smith if ((xs==0) && (ys==0) && (ze==P)) { n18 = m*n-1;} 66347c6ae99SBarry Smith if ((xs==0) && (ye==N) && (zs==0)) { n6 = (size-1)-m*(n-1);} 66447c6ae99SBarry Smith if ((xs==0) && (ye==N) && (ze==P)) { n24 = m-1;} 665*ce00eea3SSatish Balay if ((xe==M) && (ys==0) && (zs==0)) { n2 = size-m;} 666*ce00eea3SSatish Balay if ((xe==M) && (ys==0) && (ze==P)) { n20 = m*n-m;} 667*ce00eea3SSatish Balay if ((xe==M) && (ye==N) && (zs==0)) { n8 = size-m*n;} 668*ce00eea3SSatish Balay if ((xe==M) && (ye==N) && (ze==P)) { n26 = 0;} 66947c6ae99SBarry Smith 67047c6ae99SBarry Smith /* Check for when not X,Y, and Z Periodic */ 67147c6ae99SBarry Smith 67247c6ae99SBarry Smith /* If not X periodic */ 673*ce00eea3SSatish Balay if (!DMDAXPeriodic(wrap)) { 67447c6ae99SBarry Smith if (xs==0) {n0 = n3 = n6 = n9 = n12 = n15 = n18 = n21 = n24 = -2;} 675*ce00eea3SSatish Balay if (xe==M) {n2 = n5 = n8 = n11 = n14 = n17 = n20 = n23 = n26 = -2;} 67647c6ae99SBarry Smith } 67747c6ae99SBarry Smith 67847c6ae99SBarry Smith /* If not Y periodic */ 679*ce00eea3SSatish Balay if (!DMDAYPeriodic(wrap)) { 68047c6ae99SBarry Smith if (ys==0) {n0 = n1 = n2 = n9 = n10 = n11 = n18 = n19 = n20 = -2;} 68147c6ae99SBarry Smith if (ye==N) {n6 = n7 = n8 = n15 = n16 = n17 = n24 = n25 = n26 = -2;} 68247c6ae99SBarry Smith } 68347c6ae99SBarry Smith 68447c6ae99SBarry Smith /* If not Z periodic */ 685*ce00eea3SSatish Balay if (!DMDAZPeriodic(wrap)) { 68647c6ae99SBarry Smith if (zs==0) {n0 = n1 = n2 = n3 = n4 = n5 = n6 = n7 = n8 = -2;} 68747c6ae99SBarry Smith if (ze==P) {n18 = n19 = n20 = n21 = n22 = n23 = n24 = n25 = n26 = -2;} 68847c6ae99SBarry Smith } 68947c6ae99SBarry Smith 69047c6ae99SBarry Smith ierr = PetscMalloc(27*sizeof(PetscInt),&dd->neighbors);CHKERRQ(ierr); 69147c6ae99SBarry Smith dd->neighbors[0] = n0; 69247c6ae99SBarry Smith dd->neighbors[1] = n1; 69347c6ae99SBarry Smith dd->neighbors[2] = n2; 69447c6ae99SBarry Smith dd->neighbors[3] = n3; 69547c6ae99SBarry Smith dd->neighbors[4] = n4; 69647c6ae99SBarry Smith dd->neighbors[5] = n5; 69747c6ae99SBarry Smith dd->neighbors[6] = n6; 69847c6ae99SBarry Smith dd->neighbors[7] = n7; 69947c6ae99SBarry Smith dd->neighbors[8] = n8; 70047c6ae99SBarry Smith dd->neighbors[9] = n9; 70147c6ae99SBarry Smith dd->neighbors[10] = n10; 70247c6ae99SBarry Smith dd->neighbors[11] = n11; 70347c6ae99SBarry Smith dd->neighbors[12] = n12; 70447c6ae99SBarry Smith dd->neighbors[13] = rank; 70547c6ae99SBarry Smith dd->neighbors[14] = n14; 70647c6ae99SBarry Smith dd->neighbors[15] = n15; 70747c6ae99SBarry Smith dd->neighbors[16] = n16; 70847c6ae99SBarry Smith dd->neighbors[17] = n17; 70947c6ae99SBarry Smith dd->neighbors[18] = n18; 71047c6ae99SBarry Smith dd->neighbors[19] = n19; 71147c6ae99SBarry Smith dd->neighbors[20] = n20; 71247c6ae99SBarry Smith dd->neighbors[21] = n21; 71347c6ae99SBarry Smith dd->neighbors[22] = n22; 71447c6ae99SBarry Smith dd->neighbors[23] = n23; 71547c6ae99SBarry Smith dd->neighbors[24] = n24; 71647c6ae99SBarry Smith dd->neighbors[25] = n25; 71747c6ae99SBarry Smith dd->neighbors[26] = n26; 71847c6ae99SBarry Smith 71947c6ae99SBarry Smith /* If star stencil then delete the corner neighbors */ 720aa219208SBarry Smith if (stencil_type == DMDA_STENCIL_STAR) { 72147c6ae99SBarry Smith /* save information about corner neighbors */ 72247c6ae99SBarry Smith sn0 = n0; sn1 = n1; sn2 = n2; sn3 = n3; sn5 = n5; sn6 = n6; sn7 = n7; 72347c6ae99SBarry Smith sn8 = n8; sn9 = n9; sn11 = n11; sn15 = n15; sn17 = n17; sn18 = n18; 72447c6ae99SBarry Smith sn19 = n19; sn20 = n20; sn21 = n21; sn23 = n23; sn24 = n24; sn25 = n25; 72547c6ae99SBarry Smith sn26 = n26; 72647c6ae99SBarry Smith n0 = n1 = n2 = n3 = n5 = n6 = n7 = n8 = n9 = n11 = 72747c6ae99SBarry Smith n15 = n17 = n18 = n19 = n20 = n21 = n23 = n24 = n25 = n26 = -1; 72847c6ae99SBarry Smith } 72947c6ae99SBarry Smith 73047c6ae99SBarry Smith 73147c6ae99SBarry Smith ierr = PetscMalloc((Xe-Xs)*(Ye-Ys)*(Ze-Zs)*sizeof(PetscInt),&idx);CHKERRQ(ierr); 73247c6ae99SBarry Smith ierr = PetscLogObjectMemory(da,(Xe-Xs)*(Ye-Ys)*(Ze-Zs)*sizeof(PetscInt));CHKERRQ(ierr); 73347c6ae99SBarry Smith 73447c6ae99SBarry Smith nn = 0; 73547c6ae99SBarry Smith /* Bottom Level */ 73647c6ae99SBarry Smith for (k=0; k<s_z; k++) { 73747c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 73847c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 739*ce00eea3SSatish Balay x_t = lx[n0 % m]; 74047c6ae99SBarry Smith y_t = ly[(n0 % (m*n))/m]; 74147c6ae99SBarry Smith z_t = lz[n0 / (m*n)]; 74247c6ae99SBarry 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; 74347c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 74447c6ae99SBarry Smith } 74547c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 74647c6ae99SBarry Smith x_t = x; 74747c6ae99SBarry Smith y_t = ly[(n1 % (m*n))/m]; 74847c6ae99SBarry Smith z_t = lz[n1 / (m*n)]; 74947c6ae99SBarry 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; 75047c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 75147c6ae99SBarry Smith } 75247c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 753*ce00eea3SSatish Balay x_t = lx[n2 % m]; 75447c6ae99SBarry Smith y_t = ly[(n2 % (m*n))/m]; 75547c6ae99SBarry Smith z_t = lz[n2 / (m*n)]; 75647c6ae99SBarry 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; 75747c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 75847c6ae99SBarry Smith } 75947c6ae99SBarry Smith } 76047c6ae99SBarry Smith 76147c6ae99SBarry Smith for (i=0; i<y; i++) { 76247c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 763*ce00eea3SSatish Balay x_t = lx[n3 % m]; 76447c6ae99SBarry Smith y_t = y; 76547c6ae99SBarry Smith z_t = lz[n3 / (m*n)]; 76647c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 76747c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 76847c6ae99SBarry Smith } 76947c6ae99SBarry Smith 77047c6ae99SBarry Smith if (n4 >= 0) { /* middle */ 77147c6ae99SBarry Smith x_t = x; 77247c6ae99SBarry Smith y_t = y; 77347c6ae99SBarry Smith z_t = lz[n4 / (m*n)]; 77447c6ae99SBarry Smith s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 77547c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 77647c6ae99SBarry Smith } 77747c6ae99SBarry Smith 77847c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 779*ce00eea3SSatish Balay x_t = lx[n5 % m]; 78047c6ae99SBarry Smith y_t = y; 78147c6ae99SBarry Smith z_t = lz[n5 / (m*n)]; 78247c6ae99SBarry Smith s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 78347c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 78447c6ae99SBarry Smith } 78547c6ae99SBarry Smith } 78647c6ae99SBarry Smith 78747c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 78847c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 789*ce00eea3SSatish Balay x_t = lx[n6 % m]; 79047c6ae99SBarry Smith y_t = ly[(n6 % (m*n))/m]; 79147c6ae99SBarry Smith z_t = lz[n6 / (m*n)]; 79247c6ae99SBarry Smith s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 79347c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 79447c6ae99SBarry Smith } 79547c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 79647c6ae99SBarry Smith x_t = x; 79747c6ae99SBarry Smith y_t = ly[(n7 % (m*n))/m]; 79847c6ae99SBarry Smith z_t = lz[n7 / (m*n)]; 79947c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 80047c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 80147c6ae99SBarry Smith } 80247c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 803*ce00eea3SSatish Balay x_t = lx[n8 % m]; 80447c6ae99SBarry Smith y_t = ly[(n8 % (m*n))/m]; 80547c6ae99SBarry Smith z_t = lz[n8 / (m*n)]; 80647c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 80747c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 80847c6ae99SBarry Smith } 80947c6ae99SBarry Smith } 81047c6ae99SBarry Smith } 81147c6ae99SBarry Smith 81247c6ae99SBarry Smith /* Middle Level */ 81347c6ae99SBarry Smith for (k=0; k<z; k++) { 81447c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 81547c6ae99SBarry Smith if (n9 >= 0) { /* left below */ 816*ce00eea3SSatish Balay x_t = lx[n9 % m]; 81747c6ae99SBarry Smith y_t = ly[(n9 % (m*n))/m]; 81847c6ae99SBarry Smith /* z_t = z; */ 81947c6ae99SBarry Smith s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 82047c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 82147c6ae99SBarry Smith } 82247c6ae99SBarry Smith if (n10 >= 0) { /* directly below */ 82347c6ae99SBarry Smith x_t = x; 82447c6ae99SBarry Smith y_t = ly[(n10 % (m*n))/m]; 82547c6ae99SBarry Smith /* z_t = z; */ 82647c6ae99SBarry Smith s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 82747c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 82847c6ae99SBarry Smith } 82947c6ae99SBarry Smith if (n11 >= 0) { /* right below */ 830*ce00eea3SSatish Balay x_t = lx[n11 % m]; 83147c6ae99SBarry Smith y_t = ly[(n11 % (m*n))/m]; 83247c6ae99SBarry Smith /* z_t = z; */ 83347c6ae99SBarry Smith s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 83447c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 83547c6ae99SBarry Smith } 83647c6ae99SBarry Smith } 83747c6ae99SBarry Smith 83847c6ae99SBarry Smith for (i=0; i<y; i++) { 83947c6ae99SBarry Smith if (n12 >= 0) { /* directly left */ 840*ce00eea3SSatish Balay x_t = lx[n12 % m]; 84147c6ae99SBarry Smith y_t = y; 84247c6ae99SBarry Smith /* z_t = z; */ 84347c6ae99SBarry Smith s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t; 84447c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 84547c6ae99SBarry Smith } 84647c6ae99SBarry Smith 84747c6ae99SBarry Smith /* Interior */ 84847c6ae99SBarry Smith s_t = bases[rank] + i*x + k*x*y; 84947c6ae99SBarry Smith for (j=0; j<x; j++) { idx[nn++] = s_t++;} 85047c6ae99SBarry Smith 85147c6ae99SBarry Smith if (n14 >= 0) { /* directly right */ 852*ce00eea3SSatish Balay x_t = lx[n14 % m]; 85347c6ae99SBarry Smith y_t = y; 85447c6ae99SBarry Smith /* z_t = z; */ 85547c6ae99SBarry Smith s_t = bases[n14] + i*x_t + k*x_t*y_t; 85647c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 85747c6ae99SBarry Smith } 85847c6ae99SBarry Smith } 85947c6ae99SBarry Smith 86047c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 86147c6ae99SBarry Smith if (n15 >= 0) { /* left above */ 862*ce00eea3SSatish Balay x_t = lx[n15 % m]; 86347c6ae99SBarry Smith y_t = ly[(n15 % (m*n))/m]; 86447c6ae99SBarry Smith /* z_t = z; */ 86547c6ae99SBarry Smith s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t; 86647c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 86747c6ae99SBarry Smith } 86847c6ae99SBarry Smith if (n16 >= 0) { /* directly above */ 86947c6ae99SBarry Smith x_t = x; 87047c6ae99SBarry Smith y_t = ly[(n16 % (m*n))/m]; 87147c6ae99SBarry Smith /* z_t = z; */ 87247c6ae99SBarry Smith s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t; 87347c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 87447c6ae99SBarry Smith } 87547c6ae99SBarry Smith if (n17 >= 0) { /* right above */ 876*ce00eea3SSatish Balay x_t = lx[n17 % m]; 87747c6ae99SBarry Smith y_t = ly[(n17 % (m*n))/m]; 87847c6ae99SBarry Smith /* z_t = z; */ 87947c6ae99SBarry Smith s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t; 88047c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 88147c6ae99SBarry Smith } 88247c6ae99SBarry Smith } 88347c6ae99SBarry Smith } 88447c6ae99SBarry Smith 88547c6ae99SBarry Smith /* Upper Level */ 88647c6ae99SBarry Smith for (k=0; k<s_z; k++) { 88747c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 88847c6ae99SBarry Smith if (n18 >= 0) { /* left below */ 889*ce00eea3SSatish Balay x_t = lx[n18 % m]; 89047c6ae99SBarry Smith y_t = ly[(n18 % (m*n))/m]; 89147c6ae99SBarry Smith /* z_t = lz[n18 / (m*n)]; */ 89247c6ae99SBarry Smith s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 89347c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 89447c6ae99SBarry Smith } 89547c6ae99SBarry Smith if (n19 >= 0) { /* directly below */ 89647c6ae99SBarry Smith x_t = x; 89747c6ae99SBarry Smith y_t = ly[(n19 % (m*n))/m]; 89847c6ae99SBarry Smith /* z_t = lz[n19 / (m*n)]; */ 89947c6ae99SBarry Smith s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 90047c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 90147c6ae99SBarry Smith } 90247c6ae99SBarry Smith if (n20 >= 0) { /* right below */ 903*ce00eea3SSatish Balay x_t = lx[n20 % m]; 90447c6ae99SBarry Smith y_t = ly[(n20 % (m*n))/m]; 90547c6ae99SBarry Smith /* z_t = lz[n20 / (m*n)]; */ 90647c6ae99SBarry Smith s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 90747c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 90847c6ae99SBarry Smith } 90947c6ae99SBarry Smith } 91047c6ae99SBarry Smith 91147c6ae99SBarry Smith for (i=0; i<y; i++) { 91247c6ae99SBarry Smith if (n21 >= 0) { /* directly left */ 913*ce00eea3SSatish Balay x_t = lx[n21 % m]; 91447c6ae99SBarry Smith y_t = y; 91547c6ae99SBarry Smith /* z_t = lz[n21 / (m*n)]; */ 91647c6ae99SBarry Smith s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t; 91747c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 91847c6ae99SBarry Smith } 91947c6ae99SBarry Smith 92047c6ae99SBarry Smith if (n22 >= 0) { /* middle */ 92147c6ae99SBarry Smith x_t = x; 92247c6ae99SBarry Smith y_t = y; 92347c6ae99SBarry Smith /* z_t = lz[n22 / (m*n)]; */ 92447c6ae99SBarry Smith s_t = bases[n22] + i*x_t + k*x_t*y_t; 92547c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 92647c6ae99SBarry Smith } 92747c6ae99SBarry Smith 92847c6ae99SBarry Smith if (n23 >= 0) { /* directly right */ 929*ce00eea3SSatish Balay x_t = lx[n23 % m]; 93047c6ae99SBarry Smith y_t = y; 93147c6ae99SBarry Smith /* z_t = lz[n23 / (m*n)]; */ 93247c6ae99SBarry Smith s_t = bases[n23] + i*x_t + k*x_t*y_t; 93347c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 93447c6ae99SBarry Smith } 93547c6ae99SBarry Smith } 93647c6ae99SBarry Smith 93747c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 93847c6ae99SBarry Smith if (n24 >= 0) { /* left above */ 939*ce00eea3SSatish Balay x_t = lx[n24 % m]; 94047c6ae99SBarry Smith y_t = ly[(n24 % (m*n))/m]; 94147c6ae99SBarry Smith /* z_t = lz[n24 / (m*n)]; */ 94247c6ae99SBarry Smith s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t; 94347c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 94447c6ae99SBarry Smith } 94547c6ae99SBarry Smith if (n25 >= 0) { /* directly above */ 94647c6ae99SBarry Smith x_t = x; 94747c6ae99SBarry Smith y_t = ly[(n25 % (m*n))/m]; 94847c6ae99SBarry Smith /* z_t = lz[n25 / (m*n)]; */ 94947c6ae99SBarry Smith s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t; 95047c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 95147c6ae99SBarry Smith } 95247c6ae99SBarry Smith if (n26 >= 0) { /* right above */ 953*ce00eea3SSatish Balay x_t = lx[n26 % m]; 95447c6ae99SBarry Smith y_t = ly[(n26 % (m*n))/m]; 95547c6ae99SBarry Smith /* z_t = lz[n26 / (m*n)]; */ 95647c6ae99SBarry Smith s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t; 95747c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 95847c6ae99SBarry Smith } 95947c6ae99SBarry Smith } 96047c6ae99SBarry Smith } 961*ce00eea3SSatish Balay 962*ce00eea3SSatish Balay ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_COPY_VALUES,&from);CHKERRQ(ierr); 96347c6ae99SBarry Smith ierr = VecScatterCreate(global,from,local,to,>ol);CHKERRQ(ierr); 96447c6ae99SBarry Smith ierr = PetscLogObjectParent(da,gtol);CHKERRQ(ierr); 96547c6ae99SBarry Smith ierr = ISDestroy(to);CHKERRQ(ierr); 966*ce00eea3SSatish Balay ierr = ISDestroy(from);CHKERRQ(ierr); 96747c6ae99SBarry Smith 968aa219208SBarry Smith if (stencil_type == DMDA_STENCIL_STAR) { 96947c6ae99SBarry Smith n0 = sn0; n1 = sn1; n2 = sn2; n3 = sn3; n5 = sn5; n6 = sn6; n7 = sn7; 97047c6ae99SBarry Smith n8 = sn8; n9 = sn9; n11 = sn11; n15 = sn15; n17 = sn17; n18 = sn18; 97147c6ae99SBarry Smith n19 = sn19; n20 = sn20; n21 = sn21; n23 = sn23; n24 = sn24; n25 = sn25; 97247c6ae99SBarry Smith n26 = sn26; 973*ce00eea3SSatish Balay } 97447c6ae99SBarry Smith 975*ce00eea3SSatish Balay if ((stencil_type == DMDA_STENCIL_STAR) || 976*ce00eea3SSatish Balay (!DMDAXPeriodic(wrap) && DMDAXGhosted(wrap)) || 977*ce00eea3SSatish Balay (!DMDAYPeriodic(wrap) && DMDAYGhosted(wrap)) || 978*ce00eea3SSatish Balay (!DMDAZPeriodic(wrap) && DMDAZGhosted(wrap))) { 979*ce00eea3SSatish Balay /* 980*ce00eea3SSatish Balay Recompute the local to global mappings, this time keeping the 981*ce00eea3SSatish Balay information about the cross corner processor numbers. 982*ce00eea3SSatish Balay */ 98347c6ae99SBarry Smith nn = 0; 98447c6ae99SBarry Smith /* Bottom Level */ 98547c6ae99SBarry Smith for (k=0; k<s_z; k++) { 98647c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 98747c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 988*ce00eea3SSatish Balay x_t = lx[n0 % m]; 98947c6ae99SBarry Smith y_t = ly[(n0 % (m*n))/m]; 99047c6ae99SBarry Smith z_t = lz[n0 / (m*n)]; 99147c6ae99SBarry 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; 99247c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 993*ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0 && Zs-zs < 0) { 994*ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 99547c6ae99SBarry Smith } 99647c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 99747c6ae99SBarry Smith x_t = x; 99847c6ae99SBarry Smith y_t = ly[(n1 % (m*n))/m]; 99947c6ae99SBarry Smith z_t = lz[n1 / (m*n)]; 100047c6ae99SBarry 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; 100147c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 1002*ce00eea3SSatish Balay } else if (Ys-ys < 0 && Zs-zs < 0) { 1003*ce00eea3SSatish Balay for (j=0; j<x; j++) { idx[nn++] = -1;} 100447c6ae99SBarry Smith } 100547c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 1006*ce00eea3SSatish Balay x_t = lx[n2 % m]; 100747c6ae99SBarry Smith y_t = ly[(n2 % (m*n))/m]; 100847c6ae99SBarry Smith z_t = lz[n2 / (m*n)]; 100947c6ae99SBarry 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; 101047c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1011*ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0 && Zs-zs < 0) { 1012*ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 101347c6ae99SBarry Smith } 101447c6ae99SBarry Smith } 101547c6ae99SBarry Smith 101647c6ae99SBarry Smith for (i=0; i<y; i++) { 101747c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 1018*ce00eea3SSatish Balay x_t = lx[n3 % m]; 101947c6ae99SBarry Smith y_t = y; 102047c6ae99SBarry Smith z_t = lz[n3 / (m*n)]; 102147c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 102247c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1023*ce00eea3SSatish Balay } else if (Xs-xs < 0 && Zs-zs < 0) { 1024*ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 102547c6ae99SBarry Smith } 102647c6ae99SBarry Smith 102747c6ae99SBarry Smith if (n4 >= 0) { /* middle */ 102847c6ae99SBarry Smith x_t = x; 102947c6ae99SBarry Smith y_t = y; 103047c6ae99SBarry Smith z_t = lz[n4 / (m*n)]; 103147c6ae99SBarry Smith s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 103247c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 1033*ce00eea3SSatish Balay } else if (Zs-zs < 0) { 1034*ce00eea3SSatish Balay for (j=0; j<x; j++) { idx[nn++] = -1;} 103547c6ae99SBarry Smith } 103647c6ae99SBarry Smith 103747c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 1038*ce00eea3SSatish Balay x_t = lx[n5 % m]; 103947c6ae99SBarry Smith y_t = y; 104047c6ae99SBarry Smith z_t = lz[n5 / (m*n)]; 104147c6ae99SBarry Smith s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 104247c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1043*ce00eea3SSatish Balay } else if (xe-Xe < 0 && Zs-zs < 0) { 1044*ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 104547c6ae99SBarry Smith } 104647c6ae99SBarry Smith } 104747c6ae99SBarry Smith 104847c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 104947c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 1050*ce00eea3SSatish Balay x_t = lx[n6 % m]; 105147c6ae99SBarry Smith y_t = ly[(n6 % (m*n))/m]; 105247c6ae99SBarry Smith z_t = lz[n6 / (m*n)]; 105347c6ae99SBarry Smith s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 105447c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1055*ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0 && Zs-zs < 0) { 1056*ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 105747c6ae99SBarry Smith } 105847c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 105947c6ae99SBarry Smith x_t = x; 106047c6ae99SBarry Smith y_t = ly[(n7 % (m*n))/m]; 106147c6ae99SBarry Smith z_t = lz[n7 / (m*n)]; 106247c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 106347c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 1064*ce00eea3SSatish Balay } else if (ye-Ye < 0 && Zs-zs < 0) { 1065*ce00eea3SSatish Balay for (j=0; j<x; j++) { idx[nn++] = -1;} 106647c6ae99SBarry Smith } 106747c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 1068*ce00eea3SSatish Balay x_t = lx[n8 % m]; 106947c6ae99SBarry Smith y_t = ly[(n8 % (m*n))/m]; 107047c6ae99SBarry Smith z_t = lz[n8 / (m*n)]; 107147c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 107247c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1073*ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0 && Zs-zs < 0) { 1074*ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 107547c6ae99SBarry Smith } 107647c6ae99SBarry Smith } 107747c6ae99SBarry Smith } 107847c6ae99SBarry Smith 107947c6ae99SBarry Smith /* Middle Level */ 108047c6ae99SBarry Smith for (k=0; k<z; k++) { 108147c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 108247c6ae99SBarry Smith if (n9 >= 0) { /* left below */ 1083*ce00eea3SSatish Balay x_t = lx[n9 % m]; 108447c6ae99SBarry Smith y_t = ly[(n9 % (m*n))/m]; 108547c6ae99SBarry Smith /* z_t = z; */ 108647c6ae99SBarry Smith s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 108747c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1088*ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0) { 1089*ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 109047c6ae99SBarry Smith } 109147c6ae99SBarry Smith if (n10 >= 0) { /* directly below */ 109247c6ae99SBarry Smith x_t = x; 109347c6ae99SBarry Smith y_t = ly[(n10 % (m*n))/m]; 109447c6ae99SBarry Smith /* z_t = z; */ 109547c6ae99SBarry Smith s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 109647c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 1097*ce00eea3SSatish Balay } else if (Ys-ys < 0) { 1098*ce00eea3SSatish Balay for (j=0; j<x; j++) { idx[nn++] = -1;} 109947c6ae99SBarry Smith } 110047c6ae99SBarry Smith if (n11 >= 0) { /* right below */ 1101*ce00eea3SSatish Balay x_t = lx[n11 % m]; 110247c6ae99SBarry Smith y_t = ly[(n11 % (m*n))/m]; 110347c6ae99SBarry Smith /* z_t = z; */ 110447c6ae99SBarry Smith s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 110547c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1106*ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0) { 1107*ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 110847c6ae99SBarry Smith } 110947c6ae99SBarry Smith } 111047c6ae99SBarry Smith 111147c6ae99SBarry Smith for (i=0; i<y; i++) { 111247c6ae99SBarry Smith if (n12 >= 0) { /* directly left */ 1113*ce00eea3SSatish Balay x_t = lx[n12 % m]; 111447c6ae99SBarry Smith y_t = y; 111547c6ae99SBarry Smith /* z_t = z; */ 111647c6ae99SBarry Smith s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t; 111747c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1118*ce00eea3SSatish Balay } else if (Xs-xs < 0) { 1119*ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 112047c6ae99SBarry Smith } 112147c6ae99SBarry Smith 112247c6ae99SBarry Smith /* Interior */ 112347c6ae99SBarry Smith s_t = bases[rank] + i*x + k*x*y; 112447c6ae99SBarry Smith for (j=0; j<x; j++) { idx[nn++] = s_t++;} 112547c6ae99SBarry Smith 112647c6ae99SBarry Smith if (n14 >= 0) { /* directly right */ 1127*ce00eea3SSatish Balay x_t = lx[n14 % m]; 112847c6ae99SBarry Smith y_t = y; 112947c6ae99SBarry Smith /* z_t = z; */ 113047c6ae99SBarry Smith s_t = bases[n14] + i*x_t + k*x_t*y_t; 113147c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1132*ce00eea3SSatish Balay } else if (xe-Xe < 0) { 1133*ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 113447c6ae99SBarry Smith } 113547c6ae99SBarry Smith } 113647c6ae99SBarry Smith 113747c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 113847c6ae99SBarry Smith if (n15 >= 0) { /* left above */ 1139*ce00eea3SSatish Balay x_t = lx[n15 % m]; 114047c6ae99SBarry Smith y_t = ly[(n15 % (m*n))/m]; 114147c6ae99SBarry Smith /* z_t = z; */ 114247c6ae99SBarry Smith s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t; 114347c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1144*ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0) { 1145*ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 114647c6ae99SBarry Smith } 114747c6ae99SBarry Smith if (n16 >= 0) { /* directly above */ 114847c6ae99SBarry Smith x_t = x; 114947c6ae99SBarry Smith y_t = ly[(n16 % (m*n))/m]; 115047c6ae99SBarry Smith /* z_t = z; */ 115147c6ae99SBarry Smith s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t; 115247c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 1153*ce00eea3SSatish Balay } else if (ye-Ye < 0) { 1154*ce00eea3SSatish Balay for (j=0; j<x; j++) { idx[nn++] = -1;} 115547c6ae99SBarry Smith } 115647c6ae99SBarry Smith if (n17 >= 0) { /* right above */ 1157*ce00eea3SSatish Balay x_t = lx[n17 % m]; 115847c6ae99SBarry Smith y_t = ly[(n17 % (m*n))/m]; 115947c6ae99SBarry Smith /* z_t = z; */ 116047c6ae99SBarry Smith s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t; 116147c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1162*ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0) { 1163*ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 116447c6ae99SBarry Smith } 116547c6ae99SBarry Smith } 116647c6ae99SBarry Smith } 116747c6ae99SBarry Smith 116847c6ae99SBarry Smith /* Upper Level */ 116947c6ae99SBarry Smith for (k=0; k<s_z; k++) { 117047c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 117147c6ae99SBarry Smith if (n18 >= 0) { /* left below */ 1172*ce00eea3SSatish Balay x_t = lx[n18 % m]; 117347c6ae99SBarry Smith y_t = ly[(n18 % (m*n))/m]; 117447c6ae99SBarry Smith /* z_t = lz[n18 / (m*n)]; */ 117547c6ae99SBarry Smith s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 117647c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1177*ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0 && ze-Ze < 0) { 1178*ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 117947c6ae99SBarry Smith } 118047c6ae99SBarry Smith if (n19 >= 0) { /* directly below */ 118147c6ae99SBarry Smith x_t = x; 118247c6ae99SBarry Smith y_t = ly[(n19 % (m*n))/m]; 118347c6ae99SBarry Smith /* z_t = lz[n19 / (m*n)]; */ 118447c6ae99SBarry Smith s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 118547c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 1186*ce00eea3SSatish Balay } else if (Ys-ys < 0 && ze-Ze < 0) { 1187*ce00eea3SSatish Balay for (j=0; j<x; j++) { idx[nn++] = -1;} 118847c6ae99SBarry Smith } 118947c6ae99SBarry Smith if (n20 >= 0) { /* right below */ 1190*ce00eea3SSatish Balay x_t = lx[n20 % m]; 119147c6ae99SBarry Smith y_t = ly[(n20 % (m*n))/m]; 119247c6ae99SBarry Smith /* z_t = lz[n20 / (m*n)]; */ 119347c6ae99SBarry Smith s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 119447c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1195*ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0 && ze-Ze < 0) { 1196*ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 119747c6ae99SBarry Smith } 119847c6ae99SBarry Smith } 119947c6ae99SBarry Smith 120047c6ae99SBarry Smith for (i=0; i<y; i++) { 120147c6ae99SBarry Smith if (n21 >= 0) { /* directly left */ 1202*ce00eea3SSatish Balay x_t = lx[n21 % m]; 120347c6ae99SBarry Smith y_t = y; 120447c6ae99SBarry Smith /* z_t = lz[n21 / (m*n)]; */ 120547c6ae99SBarry Smith s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t; 120647c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1207*ce00eea3SSatish Balay } else if (Xs-xs < 0 && ze-Ze < 0) { 1208*ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 120947c6ae99SBarry Smith } 121047c6ae99SBarry Smith 121147c6ae99SBarry Smith if (n22 >= 0) { /* middle */ 121247c6ae99SBarry Smith x_t = x; 121347c6ae99SBarry Smith y_t = y; 121447c6ae99SBarry Smith /* z_t = lz[n22 / (m*n)]; */ 121547c6ae99SBarry Smith s_t = bases[n22] + i*x_t + k*x_t*y_t; 121647c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 1217*ce00eea3SSatish Balay } else if (ze-Ze < 0) { 1218*ce00eea3SSatish Balay for (j=0; j<x; j++) { idx[nn++] = -1;} 121947c6ae99SBarry Smith } 122047c6ae99SBarry Smith 122147c6ae99SBarry Smith if (n23 >= 0) { /* directly right */ 1222*ce00eea3SSatish Balay x_t = lx[n23 % m]; 122347c6ae99SBarry Smith y_t = y; 122447c6ae99SBarry Smith /* z_t = lz[n23 / (m*n)]; */ 122547c6ae99SBarry Smith s_t = bases[n23] + i*x_t + k*x_t*y_t; 122647c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1227*ce00eea3SSatish Balay } else if (xe-Xe < 0 && ze-Ze < 0) { 1228*ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 122947c6ae99SBarry Smith } 123047c6ae99SBarry Smith } 123147c6ae99SBarry Smith 123247c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 123347c6ae99SBarry Smith if (n24 >= 0) { /* left above */ 1234*ce00eea3SSatish Balay x_t = lx[n24 % m]; 123547c6ae99SBarry Smith y_t = ly[(n24 % (m*n))/m]; 123647c6ae99SBarry Smith /* z_t = lz[n24 / (m*n)]; */ 123747c6ae99SBarry Smith s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t; 123847c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1239*ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0 && ze-Ze < 0) { 1240*ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 124147c6ae99SBarry Smith } 124247c6ae99SBarry Smith if (n25 >= 0) { /* directly above */ 124347c6ae99SBarry Smith x_t = x; 124447c6ae99SBarry Smith y_t = ly[(n25 % (m*n))/m]; 124547c6ae99SBarry Smith /* z_t = lz[n25 / (m*n)]; */ 124647c6ae99SBarry Smith s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t; 124747c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 1248*ce00eea3SSatish Balay } else if (ye-Ye < 0 && ze-Ze < 0) { 1249*ce00eea3SSatish Balay for (j=0; j<x; j++) { idx[nn++] = -1;} 125047c6ae99SBarry Smith } 125147c6ae99SBarry Smith if (n26 >= 0) { /* right above */ 1252*ce00eea3SSatish Balay x_t = lx[n26 % m]; 125347c6ae99SBarry Smith y_t = ly[(n26 % (m*n))/m]; 125447c6ae99SBarry Smith /* z_t = lz[n26 / (m*n)]; */ 125547c6ae99SBarry Smith s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t; 125647c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1257*ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0 && ze-Ze < 0) { 1258*ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 125947c6ae99SBarry Smith } 126047c6ae99SBarry Smith } 126147c6ae99SBarry Smith } 126247c6ae99SBarry Smith } 126347c6ae99SBarry Smith /* 126447c6ae99SBarry Smith Set the local to global ordering in the global vector, this allows use 126547c6ae99SBarry Smith of VecSetValuesLocal(). 126647c6ae99SBarry Smith */ 1267*ce00eea3SSatish Balay if (nn != (Xe-Xs)*(Ye-Ys)*(Ze-Zs)) { 1268*ce00eea3SSatish Balay SETERRQ(((PetscObject)da)->comm, PETSC_ERR_SUP,"nn != count_dbg"); 1269*ce00eea3SSatish Balay PetscFunctionReturn(1); 1270*ce00eea3SSatish Balay } 1271*ce00eea3SSatish Balay ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_OWN_POINTER,<ogis);CHKERRQ(ierr); 1272*ce00eea3SSatish Balay ierr = PetscMalloc(nn*dof*sizeof(PetscInt),&idx_cpy);CHKERRQ(ierr); 1273*ce00eea3SSatish Balay /* ierr = PetscLogObjectMemory(da,nn*dof*sizeof(PetscInt));CHKERRQ(ierr); */ 1274*ce00eea3SSatish Balay ierr = ISGetIndices(ltogis, &idx_full); 1275*ce00eea3SSatish Balay ierr = PetscMemcpy(idx_cpy,idx_full,nn*dof*sizeof(PetscInt));CHKERRQ(ierr); 1276*ce00eea3SSatish Balay CHKMEMQ; 1277*ce00eea3SSatish Balay ierr = ISRestoreIndices(ltogis, &idx_full); 1278*ce00eea3SSatish Balay ierr = ISLocalToGlobalMappingCreateIS(ltogis,&da->ltogmap);CHKERRQ(ierr); 1279*ce00eea3SSatish Balay ierr = PetscLogObjectParent(da,da->ltogmap);CHKERRQ(ierr); 1280*ce00eea3SSatish Balay ierr = ISDestroy(ltogis);CHKERRQ(ierr); 12811411c6eeSJed Brown ierr = ISLocalToGlobalMappingBlock(da->ltogmap,dd->w,&da->ltogmapb);CHKERRQ(ierr); 12821411c6eeSJed Brown ierr = PetscLogObjectParent(da,da->ltogmap);CHKERRQ(ierr); 128347c6ae99SBarry Smith 1284*ce00eea3SSatish Balay ierr = PetscFree2(bases,ldims);CHKERRQ(ierr); 1285*ce00eea3SSatish Balay dd->m = m; dd->n = n; dd->p = p; 1286*ce00eea3SSatish Balay /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */ 1287*ce00eea3SSatish Balay dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = zs; dd->ze = ze; 1288*ce00eea3SSatish Balay dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = Zs; dd->Ze = Ze; 1289*ce00eea3SSatish Balay 1290*ce00eea3SSatish Balay ierr = VecDestroy(local);CHKERRQ(ierr); 1291*ce00eea3SSatish Balay ierr = VecDestroy(global);CHKERRQ(ierr); 1292*ce00eea3SSatish Balay 1293*ce00eea3SSatish Balay dd->gtol = gtol; 1294*ce00eea3SSatish Balay dd->ltog = ltog; 1295*ce00eea3SSatish Balay dd->idx = idx_cpy; 1296*ce00eea3SSatish Balay dd->Nl = nn*dof; 1297*ce00eea3SSatish Balay dd->base = base; 1298*ce00eea3SSatish Balay da->ops->view = DMView_DA_3d; 129947c6ae99SBarry Smith dd->ltol = PETSC_NULL; 130047c6ae99SBarry Smith dd->ao = PETSC_NULL; 1301*ce00eea3SSatish Balay 130247c6ae99SBarry Smith PetscFunctionReturn(0); 130347c6ae99SBarry Smith } 1304564755cdSBarry Smith 130547c6ae99SBarry Smith 130647c6ae99SBarry Smith #undef __FUNCT__ 1307aa219208SBarry Smith #define __FUNCT__ "DMDACreate3d" 130847c6ae99SBarry Smith /*@C 1309aa219208SBarry Smith DMDACreate3d - Creates an object that will manage the communication of three-dimensional 131047c6ae99SBarry Smith regular array data that is distributed across some processors. 131147c6ae99SBarry Smith 131247c6ae99SBarry Smith Collective on MPI_Comm 131347c6ae99SBarry Smith 131447c6ae99SBarry Smith Input Parameters: 131547c6ae99SBarry Smith + comm - MPI communicator 131647c6ae99SBarry Smith . wrap - type of periodicity the array should have, if any. Use one 1317aa219208SBarry Smith of DMDA_NONPERIODIC, DMDA_XPERIODIC, DMDA_YPERIODIC, DMDA_ZPERIODIC, DMDA_XYPERIODIC, DMDA_XZPERIODIC, DMDA_YZPERIODIC, DMDA_XYZPERIODIC, or DMDA_XYZGHOSTED. 1318aa219208SBarry Smith . stencil_type - Type of stencil (DMDA_STENCIL_STAR or DMDA_STENCIL_BOX) 131947c6ae99SBarry 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 132047c6ae99SBarry Smith from the command line with -da_grid_x <M> -da_grid_y <N> -da_grid_z <P>) 132147c6ae99SBarry Smith . m,n,p - corresponding number of processors in each dimension 132247c6ae99SBarry Smith (or PETSC_DECIDE to have calculated) 132347c6ae99SBarry Smith . dof - number of degrees of freedom per node 132447c6ae99SBarry Smith . lx, ly, lz - arrays containing the number of nodes in each cell along 132547c6ae99SBarry Smith the x, y, and z coordinates, or PETSC_NULL. If non-null, these 132647c6ae99SBarry Smith must be of length as m,n,p and the corresponding 132747c6ae99SBarry Smith m,n, or p cannot be PETSC_DECIDE. Sum of the lx[] entries must be M, sum of 132847c6ae99SBarry Smith the ly[] must N, sum of the lz[] must be P 132947c6ae99SBarry Smith - s - stencil width 133047c6ae99SBarry Smith 133147c6ae99SBarry Smith Output Parameter: 133247c6ae99SBarry Smith . da - the resulting distributed array object 133347c6ae99SBarry Smith 133447c6ae99SBarry Smith Options Database Key: 1335aa219208SBarry Smith + -da_view - Calls DMView() at the conclusion of DMDACreate3d() 133647c6ae99SBarry Smith . -da_grid_x <nx> - number of grid points in x direction, if M < 0 133747c6ae99SBarry Smith . -da_grid_y <ny> - number of grid points in y direction, if N < 0 133847c6ae99SBarry Smith . -da_grid_z <nz> - number of grid points in z direction, if P < 0 133947c6ae99SBarry Smith . -da_processors_x <MX> - number of processors in x direction 134047c6ae99SBarry Smith . -da_processors_y <MY> - number of processors in y direction 134147c6ae99SBarry Smith . -da_processors_z <MZ> - number of processors in z direction 134247c6ae99SBarry Smith . -da_refine_x - refinement ratio in x direction 134347c6ae99SBarry Smith . -da_refine_y - refinement ratio in y direction 134447c6ae99SBarry Smith - -da_refine_y - refinement ratio in z direction 134547c6ae99SBarry Smith 134647c6ae99SBarry Smith Level: beginner 134747c6ae99SBarry Smith 134847c6ae99SBarry Smith Notes: 1349aa219208SBarry Smith The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the 1350aa219208SBarry Smith standard 7-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes 135147c6ae99SBarry Smith the standard 27-pt stencil. 135247c6ae99SBarry Smith 1353aa219208SBarry Smith The array data itself is NOT stored in the DMDA, it is stored in Vec objects; 1354564755cdSBarry Smith The appropriate vector objects can be obtained with calls to DMCreateGlobalVector() 1355564755cdSBarry Smith and DMCreateLocalVector() and calls to VecDuplicate() if more are needed. 135647c6ae99SBarry Smith 135747c6ae99SBarry Smith .keywords: distributed array, create, three-dimensional 135847c6ae99SBarry Smith 1359aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate2d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(), 1360aa219208SBarry Smith DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMDALocalToLocalBegin(), DMDALocalToLocalEnd(), DMDASetRefinementFactor(), 1361aa219208SBarry Smith DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMDALoad(), DMDAGetOwnershipRanges() 136247c6ae99SBarry Smith 136347c6ae99SBarry Smith @*/ 13647087cfbeSBarry Smith PetscErrorCode DMDACreate3d(MPI_Comm comm,DMDAPeriodicType wrap,DMDAStencilType stencil_type,PetscInt M, 13659a42bb27SBarry 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) 136647c6ae99SBarry Smith { 136747c6ae99SBarry Smith PetscErrorCode ierr; 136847c6ae99SBarry Smith 136947c6ae99SBarry Smith PetscFunctionBegin; 1370aa219208SBarry Smith ierr = DMDACreate(comm, da);CHKERRQ(ierr); 1371aa219208SBarry Smith ierr = DMDASetDim(*da, 3);CHKERRQ(ierr); 1372aa219208SBarry Smith ierr = DMDASetSizes(*da, M, N, P);CHKERRQ(ierr); 1373aa219208SBarry Smith ierr = DMDASetNumProcs(*da, m, n, p);CHKERRQ(ierr); 1374aa219208SBarry Smith ierr = DMDASetPeriodicity(*da, wrap);CHKERRQ(ierr); 1375aa219208SBarry Smith ierr = DMDASetDof(*da, dof);CHKERRQ(ierr); 1376aa219208SBarry Smith ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr); 1377aa219208SBarry Smith ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr); 1378aa219208SBarry Smith ierr = DMDASetOwnershipRanges(*da, lx, ly, lz);CHKERRQ(ierr); 137947c6ae99SBarry Smith /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */ 13809a42bb27SBarry Smith ierr = DMSetFromOptions(*da);CHKERRQ(ierr); 13819a42bb27SBarry Smith ierr = DMSetUp(*da);CHKERRQ(ierr); 13827242296bSJed Brown ierr = DMView_DA_Private(*da);CHKERRQ(ierr); 138347c6ae99SBarry Smith PetscFunctionReturn(0); 138447c6ae99SBarry Smith } 1385