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 */ 13047c6ae99SBarry Smith if (k<0) { plane=dd->P+k; } 13147c6ae99SBarry Smith 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 */ 14047c6ae99SBarry Smith if (y<0) { ycoord = dd->N+y; } 14147c6ae99SBarry Smith 14247c6ae99SBarry Smith if (y>=dd->N) { ycoord = y-dd->N; } 14347c6ae99SBarry Smith xcoord = x; /* Keep x wrap points on drawing */ 14447c6ae99SBarry Smith 14547c6ae99SBarry Smith if (x<xmin) { xcoord = xmax - (xmin-x); } 14647c6ae99SBarry Smith 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 161aa219208SBarry Smith } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for DMDA 1d",((PetscObject)viewer)->type_name); 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; 17688661749SPeter Brune PetscInt o = dd->overlap; 17747c6ae99SBarry Smith const PetscInt dof = dd->w; 17847c6ae99SBarry Smith const PetscInt s = dd->s; 17919fd82e9SBarry Smith DMDABoundaryType bx = dd->bx; 18019fd82e9SBarry Smith DMDABoundaryType by = dd->by; 18119fd82e9SBarry Smith DMDABoundaryType bz = dd->bz; 18219fd82e9SBarry Smith DMDAStencilType stencil_type = dd->stencil_type; 18347c6ae99SBarry Smith PetscInt *lx = dd->lx; 18447c6ae99SBarry Smith PetscInt *ly = dd->ly; 18547c6ae99SBarry Smith PetscInt *lz = dd->lz; 18647c6ae99SBarry Smith MPI_Comm comm; 18747c6ae99SBarry Smith PetscMPIInt rank,size; 188ce00eea3SSatish Balay PetscInt xs = 0,xe,ys = 0,ye,zs = 0,ze,x = 0,y = 0,z = 0; 189ce00eea3SSatish Balay PetscInt Xs,Xe,Ys,Ye,Zs,Ze,IXs,IXe,IYs,IYe,IZs,IZe,start,end,pm; 190ce00eea3SSatish Balay PetscInt left,right,up,down,bottom,top,i,j,k,*idx,*idx_cpy,nn; 191ce00eea3SSatish Balay const PetscInt *idx_full; 19247c6ae99SBarry Smith PetscInt n0,n1,n2,n3,n4,n5,n6,n7,n8,n9,n10,n11,n12,n14; 19347c6ae99SBarry Smith PetscInt n15,n16,n17,n18,n19,n20,n21,n22,n23,n24,n25,n26; 194db87c5ecSEthan Coon PetscInt *bases,*ldims,base,x_t,y_t,z_t,s_t,count,s_x,s_y,s_z; 19547c6ae99SBarry Smith PetscInt sn0 = 0,sn1 = 0,sn2 = 0,sn3 = 0,sn5 = 0,sn6 = 0,sn7 = 0; 19647c6ae99SBarry Smith PetscInt sn8 = 0,sn9 = 0,sn11 = 0,sn15 = 0,sn24 = 0,sn25 = 0,sn26 = 0; 19747c6ae99SBarry Smith PetscInt sn17 = 0,sn18 = 0,sn19 = 0,sn20 = 0,sn21 = 0,sn23 = 0; 19847c6ae99SBarry Smith Vec local,global; 19947c6ae99SBarry Smith VecScatter ltog,gtol; 200ce00eea3SSatish Balay IS to,from,ltogis; 2016f951b95Secoon PetscBool twod; 20247c6ae99SBarry Smith PetscErrorCode ierr; 20347c6ae99SBarry Smith 2046f951b95Secoon 20547c6ae99SBarry Smith PetscFunctionBegin; 20647c6ae99SBarry Smith if (dof < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof); 20747c6ae99SBarry Smith if (s < 0) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s); 20847c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject) da, &comm);CHKERRQ(ierr); 2093855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES) 2103855c12bSBarry 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); 2113855c12bSBarry Smith #endif 2123855c12bSBarry Smith 21347c6ae99SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 21447c6ae99SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 21547c6ae99SBarry Smith 21647c6ae99SBarry Smith dd->dim = 3; 21747c6ae99SBarry Smith ierr = PetscMalloc(dof*sizeof(char*),&dd->fieldname);CHKERRQ(ierr); 21847c6ae99SBarry Smith ierr = PetscMemzero(dd->fieldname,dof*sizeof(char*));CHKERRQ(ierr); 21947c6ae99SBarry Smith 22047c6ae99SBarry Smith if (m != PETSC_DECIDE) { 22147c6ae99SBarry Smith if (m < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m); 22247c6ae99SBarry Smith else if (m > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size); 22347c6ae99SBarry Smith } 22447c6ae99SBarry Smith if (n != PETSC_DECIDE) { 22547c6ae99SBarry Smith if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n); 22647c6ae99SBarry Smith else if (n > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size); 22747c6ae99SBarry Smith } 22847c6ae99SBarry Smith if (p != PETSC_DECIDE) { 22947c6ae99SBarry Smith if (p < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Z direction: %D",p); 23047c6ae99SBarry Smith else if (p > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Z direction: %D %d",p,size); 23147c6ae99SBarry Smith } 2320716a85fSBarry 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); 23347c6ae99SBarry Smith 23447c6ae99SBarry Smith /* Partition the array among the processors */ 23547c6ae99SBarry Smith if (m == PETSC_DECIDE && n != PETSC_DECIDE && p != PETSC_DECIDE) { 23647c6ae99SBarry Smith m = size/(n*p); 23747c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) { 23847c6ae99SBarry Smith n = size/(m*p); 23947c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) { 24047c6ae99SBarry Smith p = size/(m*n); 24147c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) { 24247c6ae99SBarry Smith /* try for squarish distribution */ 2438f1a2a5eSBarry Smith m = (int)(0.5 + sqrt(((double)M)*((double)size)/((double)N*p))); 24447c6ae99SBarry Smith if (!m) m = 1; 24547c6ae99SBarry Smith while (m > 0) { 24647c6ae99SBarry Smith n = size/(m*p); 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 p value: p = %D",p); 25147c6ae99SBarry Smith if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;} 25247c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) { 25347c6ae99SBarry Smith /* try for squarish distribution */ 2548f1a2a5eSBarry Smith m = (int)(0.5 + sqrt(((double)M)*((double)size)/((double)P*n))); 25547c6ae99SBarry Smith if (!m) m = 1; 25647c6ae99SBarry Smith while (m > 0) { 25747c6ae99SBarry Smith p = size/(m*n); 25847c6ae99SBarry Smith if (m*n*p == size) break; 25947c6ae99SBarry Smith m--; 26047c6ae99SBarry Smith } 26147c6ae99SBarry Smith if (!m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad n value: n = %D",n); 26247c6ae99SBarry Smith if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;} 26347c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) { 26447c6ae99SBarry Smith /* try for squarish distribution */ 2658f1a2a5eSBarry Smith n = (int)(0.5 + sqrt(((double)N)*((double)size)/((double)P*m))); 26647c6ae99SBarry Smith if (!n) n = 1; 26747c6ae99SBarry Smith while (n > 0) { 26847c6ae99SBarry Smith p = size/(m*n); 26947c6ae99SBarry Smith if (m*n*p == size) break; 27047c6ae99SBarry Smith n--; 27147c6ae99SBarry Smith } 27247c6ae99SBarry Smith if (!n) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad m value: m = %D",n); 27347c6ae99SBarry Smith if (N > P && n < p) {PetscInt _n = n; n = p; p = _n;} 27447c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) { 27547c6ae99SBarry Smith /* try for squarish distribution */ 2768f1a2a5eSBarry Smith n = (PetscInt)(0.5 + pow(((double)N*N)*((double)size)/((double)P*M),(double)(1./3.))); 27747c6ae99SBarry Smith if (!n) n = 1; 27847c6ae99SBarry Smith while (n > 0) { 27947c6ae99SBarry Smith pm = size/n; 28047c6ae99SBarry Smith if (n*pm == size) break; 28147c6ae99SBarry Smith n--; 28247c6ae99SBarry Smith } 28347c6ae99SBarry Smith if (!n) n = 1; 2848f1a2a5eSBarry Smith m = (PetscInt)(0.5 + sqrt(((double)M)*((double)size)/((double)P*n))); 28547c6ae99SBarry Smith if (!m) m = 1; 28647c6ae99SBarry Smith while (m > 0) { 28747c6ae99SBarry Smith p = size/(m*n); 28847c6ae99SBarry Smith if (m*n*p == size) break; 28947c6ae99SBarry Smith m--; 29047c6ae99SBarry Smith } 29147c6ae99SBarry Smith if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;} 29247c6ae99SBarry Smith } else if (m*n*p != size) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition"); 29347c6ae99SBarry Smith 29447c6ae99SBarry Smith if (m*n*p != size) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_PLIB,"Could not find good partition"); 29547c6ae99SBarry Smith if (M < m) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m); 29647c6ae99SBarry Smith if (N < n) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n); 29747c6ae99SBarry Smith if (P < p) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in z direction is too fine! %D %D",P,p); 29847c6ae99SBarry Smith 29947c6ae99SBarry Smith /* 30047c6ae99SBarry Smith Determine locally owned region 30147c6ae99SBarry Smith [x, y, or z]s is the first local node number, [x, y, z] is the number of local nodes 30247c6ae99SBarry Smith */ 30347c6ae99SBarry Smith 30447c6ae99SBarry Smith if (!lx) { 30547c6ae99SBarry Smith ierr = PetscMalloc(m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr); 30647c6ae99SBarry Smith lx = dd->lx; 30747c6ae99SBarry Smith for (i=0; i<m; i++) { 30847c6ae99SBarry Smith lx[i] = M/m + ((M % m) > (i % m)); 30947c6ae99SBarry Smith } 31047c6ae99SBarry Smith } 31147c6ae99SBarry Smith x = lx[rank % m]; 31247c6ae99SBarry Smith xs = 0; 31347c6ae99SBarry Smith for (i=0; i<(rank%m); i++) { xs += lx[i];} 31488661749SPeter Brune if ((x < s+o) && ((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); 31547c6ae99SBarry Smith 31647c6ae99SBarry Smith if (!ly) { 31747c6ae99SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr); 31847c6ae99SBarry Smith ly = dd->ly; 31947c6ae99SBarry Smith for (i=0; i<n; i++) { 32047c6ae99SBarry Smith ly[i] = N/n + ((N % n) > (i % n)); 32147c6ae99SBarry Smith } 32247c6ae99SBarry Smith } 32347c6ae99SBarry Smith y = ly[(rank % (m*n))/m]; 32488661749SPeter Brune if ((y < s+o) && ((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); 32530729d88SBarry Smith 32647c6ae99SBarry Smith ys = 0; 32747c6ae99SBarry Smith for (i=0; i<(rank % (m*n))/m; i++) { ys += ly[i];} 32847c6ae99SBarry Smith 32947c6ae99SBarry Smith if (!lz) { 33047c6ae99SBarry Smith ierr = PetscMalloc(p*sizeof(PetscInt), &dd->lz);CHKERRQ(ierr); 33147c6ae99SBarry Smith lz = dd->lz; 33247c6ae99SBarry Smith for (i=0; i<p; i++) { 33347c6ae99SBarry Smith lz[i] = P/p + ((P % p) > (i % p)); 33447c6ae99SBarry Smith } 33547c6ae99SBarry Smith } 33647c6ae99SBarry Smith z = lz[rank/(m*n)]; 337bcea557cSEthan Coon 338fdc81ce8SEthan Coon /* note this is different than x- and y-, as we will handle as an important special 339fdc81ce8SEthan Coon case when p=P=1 and DMDA_BOUNDARY_PERIODIC and s > z. This is to deal with 2D problems 340fdc81ce8SEthan Coon in a 3D code. Additional code for this case is noted with "2d case" comments */ 3416f951b95Secoon twod = PETSC_FALSE; 3426f951b95Secoon if (P == 1) { 3436f951b95Secoon twod = PETSC_TRUE; 34488661749SPeter Brune } else if ((z < s+o) && ((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+o); 34547c6ae99SBarry Smith zs = 0; 34647c6ae99SBarry Smith for (i=0; i<(rank/(m*n)); i++) { zs += lz[i];} 34747c6ae99SBarry Smith ye = ys + y; 34847c6ae99SBarry Smith xe = xs + x; 34947c6ae99SBarry Smith ze = zs + z; 35047c6ae99SBarry Smith 35188661749SPeter Brune /* determine ghost region (Xs) and region scattered into (IXs) */ 35288661749SPeter Brune if (xs-s-o > 0) { 35388661749SPeter Brune Xs = xs - s - o; IXs = xs - s - o; 35488661749SPeter Brune } else { 35588661749SPeter Brune if (bx) { 35688661749SPeter Brune Xs = xs - s; 35788661749SPeter Brune } else { 35888661749SPeter Brune Xs = 0; 35988661749SPeter Brune } 36088661749SPeter Brune IXs = 0; 36188661749SPeter Brune } 36288661749SPeter Brune if (xe+s+o <= M) { 36388661749SPeter Brune Xe = xe + s + o; IXe = xe + s + o; 36488661749SPeter Brune } else { 36588661749SPeter Brune if (bx) { 36688661749SPeter Brune Xs = xs - s - o; Xe = xe + s; 36788661749SPeter Brune } else { 36888661749SPeter Brune Xe = M; 36988661749SPeter Brune } 37088661749SPeter Brune IXe = M; 37188661749SPeter Brune } 37247c6ae99SBarry Smith 37388661749SPeter Brune if (bx == DMDA_BOUNDARY_PERIODIC) { 37488661749SPeter Brune IXs = xs - s - o; 37588661749SPeter Brune IXe = xe + s + o; 37688661749SPeter Brune Xs = xs - s - o; 37788661749SPeter Brune Xe = xe + s + o; 37888661749SPeter Brune } 37988661749SPeter Brune 38088661749SPeter Brune if (ys-s-o > 0) { 38188661749SPeter Brune Ys = ys - s - o; IYs = ys - s - o; 38288661749SPeter Brune } else { 38388661749SPeter Brune if (by) { 38488661749SPeter Brune Ys = ys - s; 38588661749SPeter Brune } else { 38688661749SPeter Brune Ys = 0; 38788661749SPeter Brune } 38888661749SPeter Brune IYs = 0; 38988661749SPeter Brune } 39088661749SPeter Brune if (ye+s+o <= N) { 39188661749SPeter Brune Ye = ye + s + o; IYe = ye + s + o; 39288661749SPeter Brune } else { 39388661749SPeter Brune if (by) { 39488661749SPeter Brune Ye = ye + s; 39588661749SPeter Brune } else { 39688661749SPeter Brune Ye = N; 39788661749SPeter Brune } 39888661749SPeter Brune IYe = N; 39988661749SPeter Brune } 40088661749SPeter Brune 40188661749SPeter Brune if (by == DMDA_BOUNDARY_PERIODIC) { 40288661749SPeter Brune IYs = ys - s - o; 40388661749SPeter Brune IYe = ye + s + o; 40488661749SPeter Brune Ys = ys - s - o; 40588661749SPeter Brune Ye = ye + s + o; 40688661749SPeter Brune } 40788661749SPeter Brune 40888661749SPeter Brune if (zs-s-o > 0) { 40988661749SPeter Brune Zs = zs - s - o; IZs = zs - s - o; 41088661749SPeter Brune } else { 41188661749SPeter Brune if (bz) { 41288661749SPeter Brune Zs = zs - s; 41388661749SPeter Brune } else { 41488661749SPeter Brune Zs = 0; 41588661749SPeter Brune } 41688661749SPeter Brune IZs = 0; 41788661749SPeter Brune } 41888661749SPeter Brune if (ze+s+o <= P) { 41988661749SPeter Brune Ze = ze + s + o; IZe = ze + s + o; 42088661749SPeter Brune } else { 42188661749SPeter Brune if (bz) { 42288661749SPeter Brune Ze = ze + s; 42388661749SPeter Brune } else { 42488661749SPeter Brune Ze = P; 42588661749SPeter Brune } 42688661749SPeter Brune IZe = P; 42788661749SPeter Brune } 42888661749SPeter Brune 42988661749SPeter Brune if (bz == DMDA_BOUNDARY_PERIODIC) { 43088661749SPeter Brune IZs = zs - s - o; 43188661749SPeter Brune IZe = ze + s + o; 43288661749SPeter Brune Zs = zs - s - o; 43388661749SPeter Brune Ze = ze + s + o; 43488661749SPeter Brune } 43547c6ae99SBarry Smith 43647c6ae99SBarry Smith /* Resize all X parameters to reflect w */ 43788661749SPeter Brune s_x = s + o; 43888661749SPeter Brune s_y = s + o; 43988661749SPeter Brune s_z = s + o; 44047c6ae99SBarry Smith 44147c6ae99SBarry Smith /* determine starting point of each processor */ 44247c6ae99SBarry Smith nn = x*y*z; 44347c6ae99SBarry Smith ierr = PetscMalloc2(size+1,PetscInt,&bases,size,PetscInt,&ldims);CHKERRQ(ierr); 44447c6ae99SBarry Smith ierr = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr); 44547c6ae99SBarry Smith bases[0] = 0; 44647c6ae99SBarry Smith for (i=1; i<=size; i++) { 44747c6ae99SBarry Smith bases[i] = ldims[i-1]; 44847c6ae99SBarry Smith } 44947c6ae99SBarry Smith for (i=1; i<=size; i++) { 45047c6ae99SBarry Smith bases[i] += bases[i-1]; 45147c6ae99SBarry Smith } 452ce00eea3SSatish Balay base = bases[rank]*dof; 45347c6ae99SBarry Smith 45447c6ae99SBarry Smith /* allocate the base parallel and sequential vectors */ 455ce00eea3SSatish Balay dd->Nlocal = x*y*z*dof; 456778a2246SBarry Smith ierr = VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,0,&global);CHKERRQ(ierr); 457ce00eea3SSatish Balay dd->nlocal = (Xe-Xs)*(Ye-Ys)*(Ze-Zs)*dof; 458778a2246SBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,0,&local);CHKERRQ(ierr); 45947c6ae99SBarry Smith 46047c6ae99SBarry Smith /* generate appropriate vector scatters */ 46147c6ae99SBarry Smith /* local to global inserts non-ghost point region into global */ 46247c6ae99SBarry Smith ierr = VecGetOwnershipRange(global,&start,&end);CHKERRQ(ierr); 463ce00eea3SSatish Balay ierr = ISCreateStride(comm,x*y*z*dof,start,1,&to);CHKERRQ(ierr); 46447c6ae99SBarry Smith 465db87c5ecSEthan Coon count = x*y*z; 466ce00eea3SSatish Balay ierr = PetscMalloc(x*y*z*sizeof(PetscInt),&idx);CHKERRQ(ierr); 467ce00eea3SSatish Balay left = xs - Xs; right = left + x; 46847c6ae99SBarry Smith bottom = ys - Ys; top = bottom + y; 46947c6ae99SBarry Smith down = zs - Zs; up = down + z; 47047c6ae99SBarry Smith count = 0; 47147c6ae99SBarry Smith for (i=down; i<up; i++) { 47247c6ae99SBarry Smith for (j=bottom; j<top; j++) { 473ce00eea3SSatish Balay for (k=left; k<right; k++) { 474ce00eea3SSatish Balay idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 47547c6ae99SBarry Smith } 47647c6ae99SBarry Smith } 47747c6ae99SBarry Smith } 47847c6ae99SBarry Smith 47947c6ae99SBarry Smith ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&from);CHKERRQ(ierr); 48047c6ae99SBarry Smith ierr = VecScatterCreate(local,from,global,to,<og);CHKERRQ(ierr); 48147c6ae99SBarry Smith ierr = PetscLogObjectParent(da,ltog);CHKERRQ(ierr); 482fcfd50ebSBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 483fcfd50ebSBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 48447c6ae99SBarry Smith 485ce00eea3SSatish Balay /* global to local must include ghost points within the domain, 486ce00eea3SSatish Balay but not ghost points outside the domain that aren't periodic */ 48788661749SPeter Brune if (stencil_type == DMDA_STENCIL_BOX || o > 0) { 488db87c5ecSEthan Coon count = (IXe-IXs)*(IYe-IYs)*(IZe-IZs); 489db87c5ecSEthan Coon ierr = PetscMalloc(count*sizeof(PetscInt),&idx);CHKERRQ(ierr); 490ce00eea3SSatish Balay 491ce00eea3SSatish Balay left = IXs - Xs; right = left + (IXe-IXs); 492ce00eea3SSatish Balay bottom = IYs - Ys; top = bottom + (IYe-IYs); 493ce00eea3SSatish Balay down = IZs - Zs; up = down + (IZe-IZs); 494ce00eea3SSatish Balay count = 0; 495ce00eea3SSatish Balay for (i=down; i<up; i++) { 496ce00eea3SSatish Balay for (j=bottom; j<top; j++) { 497ce00eea3SSatish Balay for (k=left; k<right; k++) { 498ce00eea3SSatish Balay idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 499ce00eea3SSatish Balay } 500ce00eea3SSatish Balay } 501ce00eea3SSatish Balay } 502ce00eea3SSatish Balay ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr); 503ce00eea3SSatish Balay 50447c6ae99SBarry Smith } else { 50547c6ae99SBarry Smith /* This is way ugly! We need to list the funny cross type region */ 506db87c5ecSEthan Coon count = ((ys-IYs) + (IYe-ye))*x*z + ((xs-IXs) + (IXe-xe))*y*z + ((zs-IZs) + (IZe-ze))*x*y + x*y*z; 507db87c5ecSEthan Coon ierr = PetscMalloc(count*sizeof(PetscInt),&idx);CHKERRQ(ierr); 508ce00eea3SSatish Balay 509ce00eea3SSatish Balay left = xs - Xs; right = left + x; 51047c6ae99SBarry Smith bottom = ys - Ys; top = bottom + y; 51147c6ae99SBarry Smith down = zs - Zs; up = down + z; 51247c6ae99SBarry Smith count = 0; 513ce00eea3SSatish Balay /* the bottom chunck */ 514ce00eea3SSatish Balay for (i=(IZs-Zs); i<down; i++) { 51547c6ae99SBarry Smith for (j=bottom; j<top; j++) { 516ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 51747c6ae99SBarry Smith } 51847c6ae99SBarry Smith } 51947c6ae99SBarry Smith /* the middle piece */ 52047c6ae99SBarry Smith for (i=down; i<up; i++) { 52147c6ae99SBarry Smith /* front */ 522ce00eea3SSatish Balay for (j=(IYs-Ys); j<bottom; j++) { 523ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 52447c6ae99SBarry Smith } 52547c6ae99SBarry Smith /* middle */ 52647c6ae99SBarry Smith for (j=bottom; j<top; j++) { 527ce00eea3SSatish Balay for (k=IXs-Xs; k<IXe-Xs; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 52847c6ae99SBarry Smith } 52947c6ae99SBarry Smith /* back */ 530ce00eea3SSatish Balay for (j=top; j<top+IYe-ye; j++) { 531ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 53247c6ae99SBarry Smith } 53347c6ae99SBarry Smith } 53447c6ae99SBarry Smith /* the top piece */ 535ce00eea3SSatish Balay for (i=up; i<up+IZe-ze; i++) { 53647c6ae99SBarry Smith for (j=bottom; j<top; j++) { 537ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 53847c6ae99SBarry Smith } 53947c6ae99SBarry Smith } 54047c6ae99SBarry Smith ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr); 54147c6ae99SBarry Smith } 54247c6ae99SBarry Smith 54347c6ae99SBarry Smith /* determine who lies on each side of use stored in n24 n25 n26 54447c6ae99SBarry Smith n21 n22 n23 54547c6ae99SBarry Smith n18 n19 n20 54647c6ae99SBarry Smith 54747c6ae99SBarry Smith n15 n16 n17 54847c6ae99SBarry Smith n12 n14 54947c6ae99SBarry Smith n9 n10 n11 55047c6ae99SBarry Smith 55147c6ae99SBarry Smith n6 n7 n8 55247c6ae99SBarry Smith n3 n4 n5 55347c6ae99SBarry Smith n0 n1 n2 55447c6ae99SBarry Smith */ 55547c6ae99SBarry Smith 55647c6ae99SBarry Smith /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */ 55747c6ae99SBarry Smith /* Assume Nodes are Internal to the Cube */ 55847c6ae99SBarry Smith n0 = rank - m*n - m - 1; 55947c6ae99SBarry Smith n1 = rank - m*n - m; 56047c6ae99SBarry Smith n2 = rank - m*n - m + 1; 56147c6ae99SBarry Smith n3 = rank - m*n -1; 56247c6ae99SBarry Smith n4 = rank - m*n; 56347c6ae99SBarry Smith n5 = rank - m*n + 1; 56447c6ae99SBarry Smith n6 = rank - m*n + m - 1; 56547c6ae99SBarry Smith n7 = rank - m*n + m; 56647c6ae99SBarry Smith n8 = rank - m*n + m + 1; 56747c6ae99SBarry Smith 56847c6ae99SBarry Smith n9 = rank - m - 1; 56947c6ae99SBarry Smith n10 = rank - m; 57047c6ae99SBarry Smith n11 = rank - m + 1; 57147c6ae99SBarry Smith n12 = rank - 1; 57247c6ae99SBarry Smith n14 = rank + 1; 57347c6ae99SBarry Smith n15 = rank + m - 1; 57447c6ae99SBarry Smith n16 = rank + m; 57547c6ae99SBarry Smith n17 = rank + m + 1; 57647c6ae99SBarry Smith 57747c6ae99SBarry Smith n18 = rank + m*n - m - 1; 57847c6ae99SBarry Smith n19 = rank + m*n - m; 57947c6ae99SBarry Smith n20 = rank + m*n - m + 1; 58047c6ae99SBarry Smith n21 = rank + m*n - 1; 58147c6ae99SBarry Smith n22 = rank + m*n; 58247c6ae99SBarry Smith n23 = rank + m*n + 1; 58347c6ae99SBarry Smith n24 = rank + m*n + m - 1; 58447c6ae99SBarry Smith n25 = rank + m*n + m; 58547c6ae99SBarry Smith n26 = rank + m*n + m + 1; 58647c6ae99SBarry Smith 58747c6ae99SBarry Smith /* Assume Pieces are on Faces of Cube */ 58847c6ae99SBarry Smith 58947c6ae99SBarry Smith if (xs == 0) { /* First assume not corner or edge */ 59047c6ae99SBarry Smith n0 = rank -1 - (m*n); 59147c6ae99SBarry Smith n3 = rank + m -1 - (m*n); 59247c6ae99SBarry Smith n6 = rank + 2*m -1 - (m*n); 59347c6ae99SBarry Smith n9 = rank -1; 59447c6ae99SBarry Smith n12 = rank + m -1; 59547c6ae99SBarry Smith n15 = rank + 2*m -1; 59647c6ae99SBarry Smith n18 = rank -1 + (m*n); 59747c6ae99SBarry Smith n21 = rank + m -1 + (m*n); 59847c6ae99SBarry Smith n24 = rank + 2*m -1 + (m*n); 59947c6ae99SBarry Smith } 60047c6ae99SBarry Smith 601ce00eea3SSatish Balay if (xe == M) { /* First assume not corner or edge */ 60247c6ae99SBarry Smith n2 = rank -2*m +1 - (m*n); 60347c6ae99SBarry Smith n5 = rank - m +1 - (m*n); 60447c6ae99SBarry Smith n8 = rank +1 - (m*n); 60547c6ae99SBarry Smith n11 = rank -2*m +1; 60647c6ae99SBarry Smith n14 = rank - m +1; 60747c6ae99SBarry Smith n17 = rank +1; 60847c6ae99SBarry Smith n20 = rank -2*m +1 + (m*n); 60947c6ae99SBarry Smith n23 = rank - m +1 + (m*n); 61047c6ae99SBarry Smith n26 = rank +1 + (m*n); 61147c6ae99SBarry Smith } 61247c6ae99SBarry Smith 61347c6ae99SBarry Smith if (ys==0) { /* First assume not corner or edge */ 61447c6ae99SBarry Smith n0 = rank + m * (n-1) -1 - (m*n); 61547c6ae99SBarry Smith n1 = rank + m * (n-1) - (m*n); 61647c6ae99SBarry Smith n2 = rank + m * (n-1) +1 - (m*n); 61747c6ae99SBarry Smith n9 = rank + m * (n-1) -1; 61847c6ae99SBarry Smith n10 = rank + m * (n-1); 61947c6ae99SBarry Smith n11 = rank + m * (n-1) +1; 62047c6ae99SBarry Smith n18 = rank + m * (n-1) -1 + (m*n); 62147c6ae99SBarry Smith n19 = rank + m * (n-1) + (m*n); 62247c6ae99SBarry Smith n20 = rank + m * (n-1) +1 + (m*n); 62347c6ae99SBarry Smith } 62447c6ae99SBarry Smith 62547c6ae99SBarry Smith if (ye == N) { /* First assume not corner or edge */ 62647c6ae99SBarry Smith n6 = rank - m * (n-1) -1 - (m*n); 62747c6ae99SBarry Smith n7 = rank - m * (n-1) - (m*n); 62847c6ae99SBarry Smith n8 = rank - m * (n-1) +1 - (m*n); 62947c6ae99SBarry Smith n15 = rank - m * (n-1) -1; 63047c6ae99SBarry Smith n16 = rank - m * (n-1); 63147c6ae99SBarry Smith n17 = rank - m * (n-1) +1; 63247c6ae99SBarry Smith n24 = rank - m * (n-1) -1 + (m*n); 63347c6ae99SBarry Smith n25 = rank - m * (n-1) + (m*n); 63447c6ae99SBarry Smith n26 = rank - m * (n-1) +1 + (m*n); 63547c6ae99SBarry Smith } 63647c6ae99SBarry Smith 63747c6ae99SBarry Smith if (zs == 0) { /* First assume not corner or edge */ 63847c6ae99SBarry Smith n0 = size - (m*n) + rank - m - 1; 63947c6ae99SBarry Smith n1 = size - (m*n) + rank - m; 64047c6ae99SBarry Smith n2 = size - (m*n) + rank - m + 1; 64147c6ae99SBarry Smith n3 = size - (m*n) + rank - 1; 64247c6ae99SBarry Smith n4 = size - (m*n) + rank; 64347c6ae99SBarry Smith n5 = size - (m*n) + rank + 1; 64447c6ae99SBarry Smith n6 = size - (m*n) + rank + m - 1; 64547c6ae99SBarry Smith n7 = size - (m*n) + rank + m ; 64647c6ae99SBarry Smith n8 = size - (m*n) + rank + m + 1; 64747c6ae99SBarry Smith } 64847c6ae99SBarry Smith 64947c6ae99SBarry Smith if (ze == P) { /* First assume not corner or edge */ 65047c6ae99SBarry Smith n18 = (m*n) - (size-rank) - m - 1; 65147c6ae99SBarry Smith n19 = (m*n) - (size-rank) - m; 65247c6ae99SBarry Smith n20 = (m*n) - (size-rank) - m + 1; 65347c6ae99SBarry Smith n21 = (m*n) - (size-rank) - 1; 65447c6ae99SBarry Smith n22 = (m*n) - (size-rank); 65547c6ae99SBarry Smith n23 = (m*n) - (size-rank) + 1; 65647c6ae99SBarry Smith n24 = (m*n) - (size-rank) + m - 1; 65747c6ae99SBarry Smith n25 = (m*n) - (size-rank) + m; 65847c6ae99SBarry Smith n26 = (m*n) - (size-rank) + m + 1; 65947c6ae99SBarry Smith } 66047c6ae99SBarry Smith 66147c6ae99SBarry Smith if ((xs==0) && (zs==0)) { /* Assume an edge, not corner */ 66247c6ae99SBarry Smith n0 = size - m*n + rank + m-1 - m; 66347c6ae99SBarry Smith n3 = size - m*n + rank + m-1; 66447c6ae99SBarry Smith n6 = size - m*n + rank + m-1 + m; 66547c6ae99SBarry Smith } 66647c6ae99SBarry Smith 66747c6ae99SBarry Smith if ((xs==0) && (ze==P)) { /* Assume an edge, not corner */ 66847c6ae99SBarry Smith n18 = m*n - (size - rank) + m-1 - m; 66947c6ae99SBarry Smith n21 = m*n - (size - rank) + m-1; 67047c6ae99SBarry Smith n24 = m*n - (size - rank) + m-1 + m; 67147c6ae99SBarry Smith } 67247c6ae99SBarry Smith 67347c6ae99SBarry Smith if ((xs==0) && (ys==0)) { /* Assume an edge, not corner */ 67447c6ae99SBarry Smith n0 = rank + m*n -1 - m*n; 67547c6ae99SBarry Smith n9 = rank + m*n -1; 67647c6ae99SBarry Smith n18 = rank + m*n -1 + m*n; 67747c6ae99SBarry Smith } 67847c6ae99SBarry Smith 67947c6ae99SBarry Smith if ((xs==0) && (ye==N)) { /* Assume an edge, not corner */ 68047c6ae99SBarry Smith n6 = rank - m*(n-1) + m-1 - m*n; 68147c6ae99SBarry Smith n15 = rank - m*(n-1) + m-1; 68247c6ae99SBarry Smith n24 = rank - m*(n-1) + m-1 + m*n; 68347c6ae99SBarry Smith } 68447c6ae99SBarry Smith 685ce00eea3SSatish Balay if ((xe==M) && (zs==0)) { /* Assume an edge, not corner */ 68647c6ae99SBarry Smith n2 = size - (m*n-rank) - (m-1) - m; 68747c6ae99SBarry Smith n5 = size - (m*n-rank) - (m-1); 68847c6ae99SBarry Smith n8 = size - (m*n-rank) - (m-1) + m; 68947c6ae99SBarry Smith } 69047c6ae99SBarry Smith 691ce00eea3SSatish Balay if ((xe==M) && (ze==P)) { /* Assume an edge, not corner */ 69247c6ae99SBarry Smith n20 = m*n - (size - rank) - (m-1) - m; 69347c6ae99SBarry Smith n23 = m*n - (size - rank) - (m-1); 69447c6ae99SBarry Smith n26 = m*n - (size - rank) - (m-1) + m; 69547c6ae99SBarry Smith } 69647c6ae99SBarry Smith 697ce00eea3SSatish Balay if ((xe==M) && (ys==0)) { /* Assume an edge, not corner */ 69847c6ae99SBarry Smith n2 = rank + m*(n-1) - (m-1) - m*n; 69947c6ae99SBarry Smith n11 = rank + m*(n-1) - (m-1); 70047c6ae99SBarry Smith n20 = rank + m*(n-1) - (m-1) + m*n; 70147c6ae99SBarry Smith } 70247c6ae99SBarry Smith 703ce00eea3SSatish Balay if ((xe==M) && (ye==N)) { /* Assume an edge, not corner */ 70447c6ae99SBarry Smith n8 = rank - m*n +1 - m*n; 70547c6ae99SBarry Smith n17 = rank - m*n +1; 70647c6ae99SBarry Smith n26 = rank - m*n +1 + m*n; 70747c6ae99SBarry Smith } 70847c6ae99SBarry Smith 70947c6ae99SBarry Smith if ((ys==0) && (zs==0)) { /* Assume an edge, not corner */ 71047c6ae99SBarry Smith n0 = size - m + rank -1; 71147c6ae99SBarry Smith n1 = size - m + rank; 71247c6ae99SBarry Smith n2 = size - m + rank +1; 71347c6ae99SBarry Smith } 71447c6ae99SBarry Smith 71547c6ae99SBarry Smith if ((ys==0) && (ze==P)) { /* Assume an edge, not corner */ 71647c6ae99SBarry Smith n18 = m*n - (size - rank) + m*(n-1) -1; 71747c6ae99SBarry Smith n19 = m*n - (size - rank) + m*(n-1); 71847c6ae99SBarry Smith n20 = m*n - (size - rank) + m*(n-1) +1; 71947c6ae99SBarry Smith } 72047c6ae99SBarry Smith 72147c6ae99SBarry Smith if ((ye==N) && (zs==0)) { /* Assume an edge, not corner */ 72247c6ae99SBarry Smith n6 = size - (m*n-rank) - m * (n-1) -1; 72347c6ae99SBarry Smith n7 = size - (m*n-rank) - m * (n-1); 72447c6ae99SBarry Smith n8 = size - (m*n-rank) - m * (n-1) +1; 72547c6ae99SBarry Smith } 72647c6ae99SBarry Smith 72747c6ae99SBarry Smith if ((ye==N) && (ze==P)) { /* Assume an edge, not corner */ 72847c6ae99SBarry Smith n24 = rank - (size-m) -1; 72947c6ae99SBarry Smith n25 = rank - (size-m); 73047c6ae99SBarry Smith n26 = rank - (size-m) +1; 73147c6ae99SBarry Smith } 73247c6ae99SBarry Smith 73347c6ae99SBarry Smith /* Check for Corners */ 73447c6ae99SBarry Smith if ((xs==0) && (ys==0) && (zs==0)) { n0 = size -1;} 73547c6ae99SBarry Smith if ((xs==0) && (ys==0) && (ze==P)) { n18 = m*n-1;} 73647c6ae99SBarry Smith if ((xs==0) && (ye==N) && (zs==0)) { n6 = (size-1)-m*(n-1);} 73747c6ae99SBarry Smith if ((xs==0) && (ye==N) && (ze==P)) { n24 = m-1;} 738ce00eea3SSatish Balay if ((xe==M) && (ys==0) && (zs==0)) { n2 = size-m;} 739ce00eea3SSatish Balay if ((xe==M) && (ys==0) && (ze==P)) { n20 = m*n-m;} 740ce00eea3SSatish Balay if ((xe==M) && (ye==N) && (zs==0)) { n8 = size-m*n;} 741ce00eea3SSatish Balay if ((xe==M) && (ye==N) && (ze==P)) { n26 = 0;} 74247c6ae99SBarry Smith 74347c6ae99SBarry Smith /* Check for when not X,Y, and Z Periodic */ 74447c6ae99SBarry Smith 74547c6ae99SBarry Smith /* If not X periodic */ 7461321219cSEthan Coon if (bx != DMDA_BOUNDARY_PERIODIC) { 74747c6ae99SBarry Smith if (xs==0) {n0 = n3 = n6 = n9 = n12 = n15 = n18 = n21 = n24 = -2;} 748ce00eea3SSatish Balay if (xe==M) {n2 = n5 = n8 = n11 = n14 = n17 = n20 = n23 = n26 = -2;} 74947c6ae99SBarry Smith } 75047c6ae99SBarry Smith 75147c6ae99SBarry Smith /* If not Y periodic */ 7521321219cSEthan Coon if (by != DMDA_BOUNDARY_PERIODIC) { 75347c6ae99SBarry Smith if (ys==0) {n0 = n1 = n2 = n9 = n10 = n11 = n18 = n19 = n20 = -2;} 75447c6ae99SBarry Smith if (ye==N) {n6 = n7 = n8 = n15 = n16 = n17 = n24 = n25 = n26 = -2;} 75547c6ae99SBarry Smith } 75647c6ae99SBarry Smith 75747c6ae99SBarry Smith /* If not Z periodic */ 7581321219cSEthan Coon if (bz != DMDA_BOUNDARY_PERIODIC) { 75947c6ae99SBarry Smith if (zs==0) {n0 = n1 = n2 = n3 = n4 = n5 = n6 = n7 = n8 = -2;} 76047c6ae99SBarry Smith if (ze==P) {n18 = n19 = n20 = n21 = n22 = n23 = n24 = n25 = n26 = -2;} 76147c6ae99SBarry Smith } 76247c6ae99SBarry Smith 76347c6ae99SBarry Smith ierr = PetscMalloc(27*sizeof(PetscInt),&dd->neighbors);CHKERRQ(ierr); 76447c6ae99SBarry Smith dd->neighbors[0] = n0; 76547c6ae99SBarry Smith dd->neighbors[1] = n1; 76647c6ae99SBarry Smith dd->neighbors[2] = n2; 76747c6ae99SBarry Smith dd->neighbors[3] = n3; 76847c6ae99SBarry Smith dd->neighbors[4] = n4; 76947c6ae99SBarry Smith dd->neighbors[5] = n5; 77047c6ae99SBarry Smith dd->neighbors[6] = n6; 77147c6ae99SBarry Smith dd->neighbors[7] = n7; 77247c6ae99SBarry Smith dd->neighbors[8] = n8; 77347c6ae99SBarry Smith dd->neighbors[9] = n9; 77447c6ae99SBarry Smith dd->neighbors[10] = n10; 77547c6ae99SBarry Smith dd->neighbors[11] = n11; 77647c6ae99SBarry Smith dd->neighbors[12] = n12; 77747c6ae99SBarry Smith dd->neighbors[13] = rank; 77847c6ae99SBarry Smith dd->neighbors[14] = n14; 77947c6ae99SBarry Smith dd->neighbors[15] = n15; 78047c6ae99SBarry Smith dd->neighbors[16] = n16; 78147c6ae99SBarry Smith dd->neighbors[17] = n17; 78247c6ae99SBarry Smith dd->neighbors[18] = n18; 78347c6ae99SBarry Smith dd->neighbors[19] = n19; 78447c6ae99SBarry Smith dd->neighbors[20] = n20; 78547c6ae99SBarry Smith dd->neighbors[21] = n21; 78647c6ae99SBarry Smith dd->neighbors[22] = n22; 78747c6ae99SBarry Smith dd->neighbors[23] = n23; 78847c6ae99SBarry Smith dd->neighbors[24] = n24; 78947c6ae99SBarry Smith dd->neighbors[25] = n25; 79047c6ae99SBarry Smith dd->neighbors[26] = n26; 79147c6ae99SBarry Smith 79247c6ae99SBarry Smith /* If star stencil then delete the corner neighbors */ 79388661749SPeter Brune if (stencil_type == DMDA_STENCIL_STAR && o == 0) { 79447c6ae99SBarry Smith /* save information about corner neighbors */ 79547c6ae99SBarry Smith sn0 = n0; sn1 = n1; sn2 = n2; sn3 = n3; sn5 = n5; sn6 = n6; sn7 = n7; 79647c6ae99SBarry Smith sn8 = n8; sn9 = n9; sn11 = n11; sn15 = n15; sn17 = n17; sn18 = n18; 79747c6ae99SBarry Smith sn19 = n19; sn20 = n20; sn21 = n21; sn23 = n23; sn24 = n24; sn25 = n25; 79847c6ae99SBarry Smith sn26 = n26; 79947c6ae99SBarry Smith n0 = n1 = n2 = n3 = n5 = n6 = n7 = n8 = n9 = n11 = 80047c6ae99SBarry Smith n15 = n17 = n18 = n19 = n20 = n21 = n23 = n24 = n25 = n26 = -1; 80147c6ae99SBarry Smith } 80247c6ae99SBarry Smith 80347c6ae99SBarry Smith 80447c6ae99SBarry Smith ierr = PetscMalloc((Xe-Xs)*(Ye-Ys)*(Ze-Zs)*sizeof(PetscInt),&idx);CHKERRQ(ierr); 80547c6ae99SBarry Smith ierr = PetscLogObjectMemory(da,(Xe-Xs)*(Ye-Ys)*(Ze-Zs)*sizeof(PetscInt));CHKERRQ(ierr); 80647c6ae99SBarry Smith 80747c6ae99SBarry Smith nn = 0; 80847c6ae99SBarry Smith /* Bottom Level */ 80947c6ae99SBarry Smith for (k=0; k<s_z; k++) { 81047c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 81147c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 812ce00eea3SSatish Balay x_t = lx[n0 % m]; 81347c6ae99SBarry Smith y_t = ly[(n0 % (m*n))/m]; 81447c6ae99SBarry Smith z_t = lz[n0 / (m*n)]; 81547c6ae99SBarry 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; 8166f951b95Secoon if (twod && (s_t < 0)) {s_t = bases[n0] + x_t*y_t*z_t - (s_y-i)*x_t - s_x;} /* 2D case */ 81747c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 81847c6ae99SBarry Smith } 81947c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 82047c6ae99SBarry Smith x_t = x; 82147c6ae99SBarry Smith y_t = ly[(n1 % (m*n))/m]; 82247c6ae99SBarry Smith z_t = lz[n1 / (m*n)]; 82347c6ae99SBarry 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; 8246f951b95Secoon if (twod && (s_t < 0)) {s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t;} /* 2D case */ 82547c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 82647c6ae99SBarry Smith } 82747c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 828ce00eea3SSatish Balay x_t = lx[n2 % m]; 82947c6ae99SBarry Smith y_t = ly[(n2 % (m*n))/m]; 83047c6ae99SBarry Smith z_t = lz[n2 / (m*n)]; 83147c6ae99SBarry 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; 8326f951b95Secoon if (twod && (s_t < 0)) {s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t;} /* 2D case */ 83347c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 83447c6ae99SBarry Smith } 83547c6ae99SBarry Smith } 83647c6ae99SBarry Smith 83747c6ae99SBarry Smith for (i=0; i<y; i++) { 83847c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 839ce00eea3SSatish Balay x_t = lx[n3 % m]; 84047c6ae99SBarry Smith y_t = y; 84147c6ae99SBarry Smith z_t = lz[n3 / (m*n)]; 84247c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8436f951b95Secoon 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 */ 84447c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 84547c6ae99SBarry Smith } 84647c6ae99SBarry Smith 84747c6ae99SBarry Smith if (n4 >= 0) { /* middle */ 84847c6ae99SBarry Smith x_t = x; 84947c6ae99SBarry Smith y_t = y; 85047c6ae99SBarry Smith z_t = lz[n4 / (m*n)]; 85147c6ae99SBarry Smith s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8526f951b95Secoon if (twod && (s_t < 0)) {s_t = bases[n4] + i*x_t + x_t*y_t*z_t - x_t*y_t;} /* 2D case */ 85347c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 85447c6ae99SBarry Smith } 85547c6ae99SBarry Smith 85647c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 857ce00eea3SSatish Balay x_t = lx[n5 % m]; 85847c6ae99SBarry Smith y_t = y; 85947c6ae99SBarry Smith z_t = lz[n5 / (m*n)]; 86047c6ae99SBarry Smith s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8616f951b95Secoon if (twod && (s_t < 0)) {s_t = bases[n5] + i*x_t + x_t*y_t*z_t - x_t*y_t;} /* 2D case */ 86247c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 86347c6ae99SBarry Smith } 86447c6ae99SBarry Smith } 86547c6ae99SBarry Smith 86647c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 86747c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 868ce00eea3SSatish Balay x_t = lx[n6 % m]; 86947c6ae99SBarry Smith y_t = ly[(n6 % (m*n))/m]; 87047c6ae99SBarry Smith z_t = lz[n6 / (m*n)]; 87147c6ae99SBarry Smith s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8726f951b95Secoon 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 */ 87347c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 87447c6ae99SBarry Smith } 87547c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 87647c6ae99SBarry Smith x_t = x; 87747c6ae99SBarry Smith y_t = ly[(n7 % (m*n))/m]; 87847c6ae99SBarry Smith z_t = lz[n7 / (m*n)]; 87947c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8806f951b95Secoon 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 */ 88147c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 88247c6ae99SBarry Smith } 88347c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 884ce00eea3SSatish Balay x_t = lx[n8 % m]; 88547c6ae99SBarry Smith y_t = ly[(n8 % (m*n))/m]; 88647c6ae99SBarry Smith z_t = lz[n8 / (m*n)]; 88747c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8886f951b95Secoon 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 */ 88947c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 89047c6ae99SBarry Smith } 89147c6ae99SBarry Smith } 89247c6ae99SBarry Smith } 89347c6ae99SBarry Smith 89447c6ae99SBarry Smith /* Middle Level */ 89547c6ae99SBarry Smith for (k=0; k<z; k++) { 89647c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 89747c6ae99SBarry Smith if (n9 >= 0) { /* left below */ 898ce00eea3SSatish Balay x_t = lx[n9 % m]; 89947c6ae99SBarry Smith y_t = ly[(n9 % (m*n))/m]; 90047c6ae99SBarry Smith /* z_t = z; */ 90147c6ae99SBarry Smith s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 90247c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 90347c6ae99SBarry Smith } 90447c6ae99SBarry Smith if (n10 >= 0) { /* directly below */ 90547c6ae99SBarry Smith x_t = x; 90647c6ae99SBarry Smith y_t = ly[(n10 % (m*n))/m]; 90747c6ae99SBarry Smith /* z_t = z; */ 90847c6ae99SBarry Smith s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 90947c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 91047c6ae99SBarry Smith } 91147c6ae99SBarry Smith if (n11 >= 0) { /* right below */ 912ce00eea3SSatish Balay x_t = lx[n11 % m]; 91347c6ae99SBarry Smith y_t = ly[(n11 % (m*n))/m]; 91447c6ae99SBarry Smith /* z_t = z; */ 91547c6ae99SBarry Smith s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 91647c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 91747c6ae99SBarry Smith } 91847c6ae99SBarry Smith } 91947c6ae99SBarry Smith 92047c6ae99SBarry Smith for (i=0; i<y; i++) { 92147c6ae99SBarry Smith if (n12 >= 0) { /* directly left */ 922ce00eea3SSatish Balay x_t = lx[n12 % m]; 92347c6ae99SBarry Smith y_t = y; 92447c6ae99SBarry Smith /* z_t = z; */ 92547c6ae99SBarry Smith s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t; 92647c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 92747c6ae99SBarry Smith } 92847c6ae99SBarry Smith 92947c6ae99SBarry Smith /* Interior */ 93047c6ae99SBarry Smith s_t = bases[rank] + i*x + k*x*y; 93147c6ae99SBarry Smith for (j=0; j<x; j++) { idx[nn++] = s_t++;} 93247c6ae99SBarry Smith 93347c6ae99SBarry Smith if (n14 >= 0) { /* directly right */ 934ce00eea3SSatish Balay x_t = lx[n14 % m]; 93547c6ae99SBarry Smith y_t = y; 93647c6ae99SBarry Smith /* z_t = z; */ 93747c6ae99SBarry Smith s_t = bases[n14] + i*x_t + k*x_t*y_t; 93847c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 93947c6ae99SBarry Smith } 94047c6ae99SBarry Smith } 94147c6ae99SBarry Smith 94247c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 94347c6ae99SBarry Smith if (n15 >= 0) { /* left above */ 944ce00eea3SSatish Balay x_t = lx[n15 % m]; 94547c6ae99SBarry Smith y_t = ly[(n15 % (m*n))/m]; 94647c6ae99SBarry Smith /* z_t = z; */ 94747c6ae99SBarry Smith s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t; 94847c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 94947c6ae99SBarry Smith } 95047c6ae99SBarry Smith if (n16 >= 0) { /* directly above */ 95147c6ae99SBarry Smith x_t = x; 95247c6ae99SBarry Smith y_t = ly[(n16 % (m*n))/m]; 95347c6ae99SBarry Smith /* z_t = z; */ 95447c6ae99SBarry Smith s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t; 95547c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 95647c6ae99SBarry Smith } 95747c6ae99SBarry Smith if (n17 >= 0) { /* right above */ 958ce00eea3SSatish Balay x_t = lx[n17 % m]; 95947c6ae99SBarry Smith y_t = ly[(n17 % (m*n))/m]; 96047c6ae99SBarry Smith /* z_t = z; */ 96147c6ae99SBarry Smith s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t; 96247c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 96347c6ae99SBarry Smith } 96447c6ae99SBarry Smith } 96547c6ae99SBarry Smith } 96647c6ae99SBarry Smith 96747c6ae99SBarry Smith /* Upper Level */ 96847c6ae99SBarry Smith for (k=0; k<s_z; k++) { 96947c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 97047c6ae99SBarry Smith if (n18 >= 0) { /* left below */ 971ce00eea3SSatish Balay x_t = lx[n18 % m]; 97247c6ae99SBarry Smith y_t = ly[(n18 % (m*n))/m]; 97347c6ae99SBarry Smith /* z_t = lz[n18 / (m*n)]; */ 97447c6ae99SBarry Smith s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 9756f951b95Secoon if (twod && (s_t >= M*N*P)) {s_t = bases[n18] - (s_y-i)*x_t -s_x + x_t*y_t;} /* 2d case */ 97647c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 97747c6ae99SBarry Smith } 97847c6ae99SBarry Smith if (n19 >= 0) { /* directly below */ 97947c6ae99SBarry Smith x_t = x; 98047c6ae99SBarry Smith y_t = ly[(n19 % (m*n))/m]; 98147c6ae99SBarry Smith /* z_t = lz[n19 / (m*n)]; */ 98247c6ae99SBarry Smith s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 9836f951b95Secoon if (twod && (s_t >= M*N*P)) {s_t = bases[n19] - (s_y+1-i)*x_t + x_t*y_t;} /* 2d case */ 98447c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 98547c6ae99SBarry Smith } 98647c6ae99SBarry Smith if (n20 >= 0) { /* right below */ 987ce00eea3SSatish Balay x_t = lx[n20 % m]; 98847c6ae99SBarry Smith y_t = ly[(n20 % (m*n))/m]; 98947c6ae99SBarry Smith /* z_t = lz[n20 / (m*n)]; */ 99047c6ae99SBarry Smith s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 9916f951b95Secoon if (twod && (s_t >= M*N*P)) {s_t = bases[n20] - (s_y+1-i)*x_t + x_t*y_t;} /* 2d case */ 99247c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 99347c6ae99SBarry Smith } 99447c6ae99SBarry Smith } 99547c6ae99SBarry Smith 99647c6ae99SBarry Smith for (i=0; i<y; i++) { 99747c6ae99SBarry Smith if (n21 >= 0) { /* directly left */ 998ce00eea3SSatish Balay x_t = lx[n21 % m]; 99947c6ae99SBarry Smith y_t = y; 100047c6ae99SBarry Smith /* z_t = lz[n21 / (m*n)]; */ 100147c6ae99SBarry Smith s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t; 10026f951b95Secoon if (twod && (s_t >= M*N*P)) {s_t = bases[n21] + (i+1)*x_t - s_x;} /* 2d case */ 100347c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 100447c6ae99SBarry Smith } 100547c6ae99SBarry Smith 100647c6ae99SBarry Smith if (n22 >= 0) { /* middle */ 100747c6ae99SBarry Smith x_t = x; 100847c6ae99SBarry Smith y_t = y; 100947c6ae99SBarry Smith /* z_t = lz[n22 / (m*n)]; */ 101047c6ae99SBarry Smith s_t = bases[n22] + i*x_t + k*x_t*y_t; 10116f951b95Secoon if (twod && (s_t >= M*N*P)) {s_t = bases[n22] + i*x_t;} /* 2d case */ 101247c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 101347c6ae99SBarry Smith } 101447c6ae99SBarry Smith 101547c6ae99SBarry Smith if (n23 >= 0) { /* directly right */ 1016ce00eea3SSatish Balay x_t = lx[n23 % m]; 101747c6ae99SBarry Smith y_t = y; 101847c6ae99SBarry Smith /* z_t = lz[n23 / (m*n)]; */ 101947c6ae99SBarry Smith s_t = bases[n23] + i*x_t + k*x_t*y_t; 10206f951b95Secoon if (twod && (s_t >= M*N*P)) {s_t = bases[n23] + i*x_t;} /* 2d case */ 102147c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 102247c6ae99SBarry Smith } 102347c6ae99SBarry Smith } 102447c6ae99SBarry Smith 102547c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 102647c6ae99SBarry Smith if (n24 >= 0) { /* left above */ 1027ce00eea3SSatish Balay x_t = lx[n24 % m]; 102847c6ae99SBarry Smith y_t = ly[(n24 % (m*n))/m]; 102947c6ae99SBarry Smith /* z_t = lz[n24 / (m*n)]; */ 103047c6ae99SBarry Smith s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t; 10316f951b95Secoon if (twod && (s_t >= M*N*P)) {s_t = bases[n24] + i*x_t - s_x;} /* 2d case */ 103247c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 103347c6ae99SBarry Smith } 103447c6ae99SBarry Smith if (n25 >= 0) { /* directly above */ 103547c6ae99SBarry Smith x_t = x; 103647c6ae99SBarry Smith y_t = ly[(n25 % (m*n))/m]; 103747c6ae99SBarry Smith /* z_t = lz[n25 / (m*n)]; */ 103847c6ae99SBarry Smith s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t; 10396f951b95Secoon if (twod && (s_t >= M*N*P)) {s_t = bases[n25] + (i-1)*x_t;} /* 2d case */ 104047c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 104147c6ae99SBarry Smith } 104247c6ae99SBarry Smith if (n26 >= 0) { /* right above */ 1043ce00eea3SSatish Balay x_t = lx[n26 % m]; 104447c6ae99SBarry Smith y_t = ly[(n26 % (m*n))/m]; 104547c6ae99SBarry Smith /* z_t = lz[n26 / (m*n)]; */ 104647c6ae99SBarry Smith s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t; 10476f951b95Secoon if (twod && (s_t >= M*N*P)) {s_t = bases[n26] + (i-1)*x_t;} /* 2d case */ 104847c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 104947c6ae99SBarry Smith } 105047c6ae99SBarry Smith } 105147c6ae99SBarry Smith } 1052ce00eea3SSatish Balay 1053ce00eea3SSatish Balay ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_COPY_VALUES,&from);CHKERRQ(ierr); 105447c6ae99SBarry Smith ierr = VecScatterCreate(global,from,local,to,>ol);CHKERRQ(ierr); 105547c6ae99SBarry Smith ierr = PetscLogObjectParent(da,gtol);CHKERRQ(ierr); 1056fcfd50ebSBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 1057fcfd50ebSBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 105847c6ae99SBarry Smith 105988661749SPeter Brune if (stencil_type == DMDA_STENCIL_STAR && o == 0) { 106047c6ae99SBarry Smith n0 = sn0; n1 = sn1; n2 = sn2; n3 = sn3; n5 = sn5; n6 = sn6; n7 = sn7; 106147c6ae99SBarry Smith n8 = sn8; n9 = sn9; n11 = sn11; n15 = sn15; n17 = sn17; n18 = sn18; 106247c6ae99SBarry Smith n19 = sn19; n20 = sn20; n21 = sn21; n23 = sn23; n24 = sn24; n25 = sn25; 106347c6ae99SBarry Smith n26 = sn26; 1064ce00eea3SSatish Balay } 106547c6ae99SBarry Smith 106688661749SPeter Brune if (((stencil_type == DMDA_STENCIL_STAR) || 10671321219cSEthan Coon (bx != DMDA_BOUNDARY_PERIODIC && bx) || 10681321219cSEthan Coon (by != DMDA_BOUNDARY_PERIODIC && by) || 106988661749SPeter Brune (bz != DMDA_BOUNDARY_PERIODIC && bz)) && o == 0) { 1070ce00eea3SSatish Balay /* 1071ce00eea3SSatish Balay Recompute the local to global mappings, this time keeping the 1072ce00eea3SSatish Balay information about the cross corner processor numbers. 1073ce00eea3SSatish Balay */ 107447c6ae99SBarry Smith nn = 0; 107547c6ae99SBarry Smith /* Bottom Level */ 107647c6ae99SBarry Smith for (k=0; k<s_z; k++) { 107747c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 107847c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 1079ce00eea3SSatish Balay x_t = lx[n0 % m]; 108047c6ae99SBarry Smith y_t = ly[(n0 % (m*n))/m]; 108147c6ae99SBarry Smith z_t = lz[n0 / (m*n)]; 108247c6ae99SBarry 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; 108347c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1084ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0 && Zs-zs < 0) { 1085ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 108647c6ae99SBarry Smith } 108747c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 108847c6ae99SBarry Smith x_t = x; 108947c6ae99SBarry Smith y_t = ly[(n1 % (m*n))/m]; 109047c6ae99SBarry Smith z_t = lz[n1 / (m*n)]; 109147c6ae99SBarry 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; 109247c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 1093ce00eea3SSatish Balay } else if (Ys-ys < 0 && Zs-zs < 0) { 1094ce00eea3SSatish Balay for (j=0; j<x; j++) { idx[nn++] = -1;} 109547c6ae99SBarry Smith } 109647c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 1097ce00eea3SSatish Balay x_t = lx[n2 % m]; 109847c6ae99SBarry Smith y_t = ly[(n2 % (m*n))/m]; 109947c6ae99SBarry Smith z_t = lz[n2 / (m*n)]; 110047c6ae99SBarry 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; 110147c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1102ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0 && Zs-zs < 0) { 1103ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 110447c6ae99SBarry Smith } 110547c6ae99SBarry Smith } 110647c6ae99SBarry Smith 110747c6ae99SBarry Smith for (i=0; i<y; i++) { 110847c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 1109ce00eea3SSatish Balay x_t = lx[n3 % m]; 111047c6ae99SBarry Smith y_t = y; 111147c6ae99SBarry Smith z_t = lz[n3 / (m*n)]; 111247c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 111347c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1114ce00eea3SSatish Balay } else if (Xs-xs < 0 && Zs-zs < 0) { 1115ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 111647c6ae99SBarry Smith } 111747c6ae99SBarry Smith 111847c6ae99SBarry Smith if (n4 >= 0) { /* middle */ 111947c6ae99SBarry Smith x_t = x; 112047c6ae99SBarry Smith y_t = y; 112147c6ae99SBarry Smith z_t = lz[n4 / (m*n)]; 112247c6ae99SBarry Smith s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 112347c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 1124ce00eea3SSatish Balay } else if (Zs-zs < 0) { 1125ce00eea3SSatish Balay for (j=0; j<x; j++) { idx[nn++] = -1;} 112647c6ae99SBarry Smith } 112747c6ae99SBarry Smith 112847c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 1129ce00eea3SSatish Balay x_t = lx[n5 % m]; 113047c6ae99SBarry Smith y_t = y; 113147c6ae99SBarry Smith z_t = lz[n5 / (m*n)]; 113247c6ae99SBarry Smith s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 113347c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1134ce00eea3SSatish Balay } else if (xe-Xe < 0 && Zs-zs < 0) { 1135ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 113647c6ae99SBarry Smith } 113747c6ae99SBarry Smith } 113847c6ae99SBarry Smith 113947c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 114047c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 1141ce00eea3SSatish Balay x_t = lx[n6 % m]; 114247c6ae99SBarry Smith y_t = ly[(n6 % (m*n))/m]; 114347c6ae99SBarry Smith z_t = lz[n6 / (m*n)]; 114447c6ae99SBarry Smith s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 114547c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1146ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0 && Zs-zs < 0) { 1147ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 114847c6ae99SBarry Smith } 114947c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 115047c6ae99SBarry Smith x_t = x; 115147c6ae99SBarry Smith y_t = ly[(n7 % (m*n))/m]; 115247c6ae99SBarry Smith z_t = lz[n7 / (m*n)]; 115347c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 115447c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 1155ce00eea3SSatish Balay } else if (ye-Ye < 0 && Zs-zs < 0) { 1156ce00eea3SSatish Balay for (j=0; j<x; j++) { idx[nn++] = -1;} 115747c6ae99SBarry Smith } 115847c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 1159ce00eea3SSatish Balay x_t = lx[n8 % m]; 116047c6ae99SBarry Smith y_t = ly[(n8 % (m*n))/m]; 116147c6ae99SBarry Smith z_t = lz[n8 / (m*n)]; 116247c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 116347c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1164ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0 && Zs-zs < 0) { 1165ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 116647c6ae99SBarry Smith } 116747c6ae99SBarry Smith } 116847c6ae99SBarry Smith } 116947c6ae99SBarry Smith 117047c6ae99SBarry Smith /* Middle Level */ 117147c6ae99SBarry Smith for (k=0; k<z; k++) { 117247c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 117347c6ae99SBarry Smith if (n9 >= 0) { /* left below */ 1174ce00eea3SSatish Balay x_t = lx[n9 % m]; 117547c6ae99SBarry Smith y_t = ly[(n9 % (m*n))/m]; 117647c6ae99SBarry Smith /* z_t = z; */ 117747c6ae99SBarry Smith s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 117847c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1179ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0) { 1180ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 118147c6ae99SBarry Smith } 118247c6ae99SBarry Smith if (n10 >= 0) { /* directly below */ 118347c6ae99SBarry Smith x_t = x; 118447c6ae99SBarry Smith y_t = ly[(n10 % (m*n))/m]; 118547c6ae99SBarry Smith /* z_t = z; */ 118647c6ae99SBarry Smith s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 118747c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 1188ce00eea3SSatish Balay } else if (Ys-ys < 0) { 1189ce00eea3SSatish Balay for (j=0; j<x; j++) { idx[nn++] = -1;} 119047c6ae99SBarry Smith } 119147c6ae99SBarry Smith if (n11 >= 0) { /* right below */ 1192ce00eea3SSatish Balay x_t = lx[n11 % m]; 119347c6ae99SBarry Smith y_t = ly[(n11 % (m*n))/m]; 119447c6ae99SBarry Smith /* z_t = z; */ 119547c6ae99SBarry Smith s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 119647c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1197ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0) { 1198ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 119947c6ae99SBarry Smith } 120047c6ae99SBarry Smith } 120147c6ae99SBarry Smith 120247c6ae99SBarry Smith for (i=0; i<y; i++) { 120347c6ae99SBarry Smith if (n12 >= 0) { /* directly left */ 1204ce00eea3SSatish Balay x_t = lx[n12 % m]; 120547c6ae99SBarry Smith y_t = y; 120647c6ae99SBarry Smith /* z_t = z; */ 120747c6ae99SBarry Smith s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t; 120847c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1209ce00eea3SSatish Balay } else if (Xs-xs < 0) { 1210ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 121147c6ae99SBarry Smith } 121247c6ae99SBarry Smith 121347c6ae99SBarry Smith /* Interior */ 121447c6ae99SBarry Smith s_t = bases[rank] + i*x + k*x*y; 121547c6ae99SBarry Smith for (j=0; j<x; j++) { idx[nn++] = s_t++;} 121647c6ae99SBarry Smith 121747c6ae99SBarry Smith if (n14 >= 0) { /* directly right */ 1218ce00eea3SSatish Balay x_t = lx[n14 % m]; 121947c6ae99SBarry Smith y_t = y; 122047c6ae99SBarry Smith /* z_t = z; */ 122147c6ae99SBarry Smith s_t = bases[n14] + i*x_t + k*x_t*y_t; 122247c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1223ce00eea3SSatish Balay } else if (xe-Xe < 0) { 1224ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 122547c6ae99SBarry Smith } 122647c6ae99SBarry Smith } 122747c6ae99SBarry Smith 122847c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 122947c6ae99SBarry Smith if (n15 >= 0) { /* left above */ 1230ce00eea3SSatish Balay x_t = lx[n15 % m]; 123147c6ae99SBarry Smith y_t = ly[(n15 % (m*n))/m]; 123247c6ae99SBarry Smith /* z_t = z; */ 123347c6ae99SBarry Smith s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t; 123447c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1235ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0) { 1236ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 123747c6ae99SBarry Smith } 123847c6ae99SBarry Smith if (n16 >= 0) { /* directly above */ 123947c6ae99SBarry Smith x_t = x; 124047c6ae99SBarry Smith y_t = ly[(n16 % (m*n))/m]; 124147c6ae99SBarry Smith /* z_t = z; */ 124247c6ae99SBarry Smith s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t; 124347c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 1244ce00eea3SSatish Balay } else if (ye-Ye < 0) { 1245ce00eea3SSatish Balay for (j=0; j<x; j++) { idx[nn++] = -1;} 124647c6ae99SBarry Smith } 124747c6ae99SBarry Smith if (n17 >= 0) { /* right above */ 1248ce00eea3SSatish Balay x_t = lx[n17 % m]; 124947c6ae99SBarry Smith y_t = ly[(n17 % (m*n))/m]; 125047c6ae99SBarry Smith /* z_t = z; */ 125147c6ae99SBarry Smith s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t; 125247c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1253ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0) { 1254ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 125547c6ae99SBarry Smith } 125647c6ae99SBarry Smith } 125747c6ae99SBarry Smith } 125847c6ae99SBarry Smith 125947c6ae99SBarry Smith /* Upper Level */ 126047c6ae99SBarry Smith for (k=0; k<s_z; k++) { 126147c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 126247c6ae99SBarry Smith if (n18 >= 0) { /* left below */ 1263ce00eea3SSatish Balay x_t = lx[n18 % m]; 126447c6ae99SBarry Smith y_t = ly[(n18 % (m*n))/m]; 126547c6ae99SBarry Smith /* z_t = lz[n18 / (m*n)]; */ 126647c6ae99SBarry Smith s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 126747c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1268ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0 && ze-Ze < 0) { 1269ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 127047c6ae99SBarry Smith } 127147c6ae99SBarry Smith if (n19 >= 0) { /* directly below */ 127247c6ae99SBarry Smith x_t = x; 127347c6ae99SBarry Smith y_t = ly[(n19 % (m*n))/m]; 127447c6ae99SBarry Smith /* z_t = lz[n19 / (m*n)]; */ 127547c6ae99SBarry Smith s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 127647c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 1277ce00eea3SSatish Balay } else if (Ys-ys < 0 && ze-Ze < 0) { 1278ce00eea3SSatish Balay for (j=0; j<x; j++) { idx[nn++] = -1;} 127947c6ae99SBarry Smith } 128047c6ae99SBarry Smith if (n20 >= 0) { /* right below */ 1281ce00eea3SSatish Balay x_t = lx[n20 % m]; 128247c6ae99SBarry Smith y_t = ly[(n20 % (m*n))/m]; 128347c6ae99SBarry Smith /* z_t = lz[n20 / (m*n)]; */ 128447c6ae99SBarry Smith s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 128547c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1286ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0 && ze-Ze < 0) { 1287ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 128847c6ae99SBarry Smith } 128947c6ae99SBarry Smith } 129047c6ae99SBarry Smith 129147c6ae99SBarry Smith for (i=0; i<y; i++) { 129247c6ae99SBarry Smith if (n21 >= 0) { /* directly left */ 1293ce00eea3SSatish Balay x_t = lx[n21 % m]; 129447c6ae99SBarry Smith y_t = y; 129547c6ae99SBarry Smith /* z_t = lz[n21 / (m*n)]; */ 129647c6ae99SBarry Smith s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t; 129747c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1298ce00eea3SSatish Balay } else if (Xs-xs < 0 && ze-Ze < 0) { 1299ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 130047c6ae99SBarry Smith } 130147c6ae99SBarry Smith 130247c6ae99SBarry Smith if (n22 >= 0) { /* middle */ 130347c6ae99SBarry Smith x_t = x; 130447c6ae99SBarry Smith y_t = y; 130547c6ae99SBarry Smith /* z_t = lz[n22 / (m*n)]; */ 130647c6ae99SBarry Smith s_t = bases[n22] + i*x_t + k*x_t*y_t; 130747c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 1308ce00eea3SSatish Balay } else if (ze-Ze < 0) { 1309ce00eea3SSatish Balay for (j=0; j<x; j++) { idx[nn++] = -1;} 131047c6ae99SBarry 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; 131747c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1318ce00eea3SSatish Balay } else if (xe-Xe < 0 && ze-Ze < 0) { 1319ce00eea3SSatish Balay 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; 132947c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1330ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0 && ze-Ze < 0) { 1331ce00eea3SSatish Balay 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; 133847c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 1339ce00eea3SSatish Balay } else if (ye-Ye < 0 && ze-Ze < 0) { 1340ce00eea3SSatish Balay 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; 134747c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1348ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0 && ze-Ze < 0) { 1349ce00eea3SSatish Balay 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); 1361ce00eea3SSatish Balay ierr = ISGetIndices(ltogis, &idx_full); 1362ce00eea3SSatish Balay ierr = PetscMemcpy(idx_cpy,idx_full,nn*dof*sizeof(PetscInt));CHKERRQ(ierr); 1363ce00eea3SSatish Balay ierr = ISRestoreIndices(ltogis, &idx_full); 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; 1387ce00eea3SSatish Balay 138847c6ae99SBarry Smith PetscFunctionReturn(0); 138947c6ae99SBarry Smith } 1390564755cdSBarry Smith 139147c6ae99SBarry Smith 139247c6ae99SBarry Smith #undef __FUNCT__ 1393aa219208SBarry Smith #define __FUNCT__ "DMDACreate3d" 139447c6ae99SBarry Smith /*@C 1395aa219208SBarry Smith DMDACreate3d - Creates an object that will manage the communication of three-dimensional 139647c6ae99SBarry Smith regular array data that is distributed across some processors. 139747c6ae99SBarry Smith 139847c6ae99SBarry Smith Collective on MPI_Comm 139947c6ae99SBarry Smith 140047c6ae99SBarry Smith Input Parameters: 140147c6ae99SBarry Smith + comm - MPI communicator 14021321219cSEthan Coon . bx,by,bz - type of ghost nodes the array have. 14031321219cSEthan Coon Use one of DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_PERIODIC. 1404aa219208SBarry Smith . stencil_type - Type of stencil (DMDA_STENCIL_STAR or DMDA_STENCIL_BOX) 140547c6ae99SBarry Smith . M,N,P - global dimension in each direction of the array (use -M, -N, and or -P to indicate that it may be set to a different value 140647c6ae99SBarry Smith from the command line with -da_grid_x <M> -da_grid_y <N> -da_grid_z <P>) 140747c6ae99SBarry Smith . m,n,p - corresponding number of processors in each dimension 140847c6ae99SBarry Smith (or PETSC_DECIDE to have calculated) 140947c6ae99SBarry Smith . dof - number of degrees of freedom per node 141010d7c030SMatthew G Knepley . s - stencil width 141110d7c030SMatthew G Knepley - lx, ly, lz - arrays containing the number of nodes in each cell along 141247c6ae99SBarry Smith the x, y, and z coordinates, or PETSC_NULL. If non-null, these 141347c6ae99SBarry Smith must be of length as m,n,p and the corresponding 141447c6ae99SBarry Smith m,n, or p cannot be PETSC_DECIDE. Sum of the lx[] entries must be M, sum of 141547c6ae99SBarry Smith the ly[] must N, sum of the lz[] must be P 141647c6ae99SBarry Smith 141747c6ae99SBarry Smith Output Parameter: 141847c6ae99SBarry Smith . da - the resulting distributed array object 141947c6ae99SBarry Smith 142047c6ae99SBarry Smith Options Database Key: 1421aa219208SBarry Smith + -da_view - Calls DMView() at the conclusion of DMDACreate3d() 142247c6ae99SBarry Smith . -da_grid_x <nx> - number of grid points in x direction, if M < 0 142347c6ae99SBarry Smith . -da_grid_y <ny> - number of grid points in y direction, if N < 0 142447c6ae99SBarry Smith . -da_grid_z <nz> - number of grid points in z direction, if P < 0 142547c6ae99SBarry Smith . -da_processors_x <MX> - number of processors in x direction 142647c6ae99SBarry Smith . -da_processors_y <MY> - number of processors in y direction 142747c6ae99SBarry Smith . -da_processors_z <MZ> - number of processors in z direction 1428e0f5d30fSBarry Smith . -da_refine_x <rx> - refinement ratio in x direction 1429e0f5d30fSBarry Smith . -da_refine_y <ry> - refinement ratio in y direction 1430e0f5d30fSBarry Smith . -da_refine_z <rz>- refinement ratio in z directio 1431e0f5d30fSBarry Smith - -da_refine <n> - refine the DMDA n times before creating it, , if M, N, or P < 0 143247c6ae99SBarry Smith 143347c6ae99SBarry Smith Level: beginner 143447c6ae99SBarry Smith 143547c6ae99SBarry Smith Notes: 1436aa219208SBarry Smith The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the 1437aa219208SBarry Smith standard 7-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes 143847c6ae99SBarry Smith the standard 27-pt stencil. 143947c6ae99SBarry Smith 1440aa219208SBarry Smith The array data itself is NOT stored in the DMDA, it is stored in Vec objects; 1441564755cdSBarry Smith The appropriate vector objects can be obtained with calls to DMCreateGlobalVector() 1442564755cdSBarry Smith and DMCreateLocalVector() and calls to VecDuplicate() if more are needed. 144347c6ae99SBarry Smith 144447c6ae99SBarry Smith .keywords: distributed array, create, three-dimensional 144547c6ae99SBarry Smith 1446aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate2d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(), 1447aa219208SBarry Smith DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMDALocalToLocalBegin(), DMDALocalToLocalEnd(), DMDASetRefinementFactor(), 1448d461ba97SBarry Smith DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges() 144947c6ae99SBarry Smith 145047c6ae99SBarry Smith @*/ 14511321219cSEthan Coon PetscErrorCode DMDACreate3d(MPI_Comm comm,DMDABoundaryType bx,DMDABoundaryType by,DMDABoundaryType bz,DMDAStencilType stencil_type,PetscInt M, 14529a42bb27SBarry Smith PetscInt N,PetscInt P,PetscInt m,PetscInt n,PetscInt p,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],const PetscInt lz[],DM *da) 145347c6ae99SBarry Smith { 145447c6ae99SBarry Smith PetscErrorCode ierr; 145547c6ae99SBarry Smith 145647c6ae99SBarry Smith PetscFunctionBegin; 1457aa219208SBarry Smith ierr = DMDACreate(comm, da);CHKERRQ(ierr); 1458aa219208SBarry Smith ierr = DMDASetDim(*da, 3);CHKERRQ(ierr); 1459aa219208SBarry Smith ierr = DMDASetSizes(*da, M, N, P);CHKERRQ(ierr); 1460aa219208SBarry Smith ierr = DMDASetNumProcs(*da, m, n, p);CHKERRQ(ierr); 14611321219cSEthan Coon ierr = DMDASetBoundaryType(*da, bx, by, bz);CHKERRQ(ierr); 1462aa219208SBarry Smith ierr = DMDASetDof(*da, dof);CHKERRQ(ierr); 1463aa219208SBarry Smith ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr); 1464aa219208SBarry Smith ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr); 1465aa219208SBarry Smith ierr = DMDASetOwnershipRanges(*da, lx, ly, lz);CHKERRQ(ierr); 146647c6ae99SBarry Smith /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */ 14679a42bb27SBarry Smith ierr = DMSetFromOptions(*da);CHKERRQ(ierr); 14689a42bb27SBarry Smith ierr = DMSetUp(*da);CHKERRQ(ierr); 1469*ca266f36SBarry Smith ierr = DMViewFromOptions(*da,"-dm_view");CHKERRQ(ierr); 147047c6ae99SBarry Smith PetscFunctionReturn(0); 147147c6ae99SBarry Smith } 1472