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; 17647c6ae99SBarry Smith const PetscInt dof = dd->w; 17747c6ae99SBarry Smith const PetscInt s = dd->s; 17819fd82e9SBarry Smith DMDABoundaryType bx = dd->bx; 17919fd82e9SBarry Smith DMDABoundaryType by = dd->by; 18019fd82e9SBarry Smith DMDABoundaryType bz = dd->bz; 18119fd82e9SBarry Smith DMDAStencilType stencil_type = dd->stencil_type; 18247c6ae99SBarry Smith PetscInt *lx = dd->lx; 18347c6ae99SBarry Smith PetscInt *ly = dd->ly; 18447c6ae99SBarry Smith PetscInt *lz = dd->lz; 18547c6ae99SBarry Smith MPI_Comm comm; 18647c6ae99SBarry Smith PetscMPIInt rank,size; 187ce00eea3SSatish Balay PetscInt xs = 0,xe,ys = 0,ye,zs = 0,ze,x = 0,y = 0,z = 0; 188ce00eea3SSatish Balay PetscInt Xs,Xe,Ys,Ye,Zs,Ze,IXs,IXe,IYs,IYe,IZs,IZe,start,end,pm; 189ce00eea3SSatish Balay PetscInt left,right,up,down,bottom,top,i,j,k,*idx,*idx_cpy,nn; 190ce00eea3SSatish Balay const PetscInt *idx_full; 19147c6ae99SBarry Smith PetscInt n0,n1,n2,n3,n4,n5,n6,n7,n8,n9,n10,n11,n12,n14; 19247c6ae99SBarry Smith PetscInt n15,n16,n17,n18,n19,n20,n21,n22,n23,n24,n25,n26; 193db87c5ecSEthan Coon PetscInt *bases,*ldims,base,x_t,y_t,z_t,s_t,count,s_x,s_y,s_z; 19447c6ae99SBarry Smith PetscInt sn0 = 0,sn1 = 0,sn2 = 0,sn3 = 0,sn5 = 0,sn6 = 0,sn7 = 0; 19547c6ae99SBarry Smith PetscInt sn8 = 0,sn9 = 0,sn11 = 0,sn15 = 0,sn24 = 0,sn25 = 0,sn26 = 0; 19647c6ae99SBarry Smith PetscInt sn17 = 0,sn18 = 0,sn19 = 0,sn20 = 0,sn21 = 0,sn23 = 0; 19747c6ae99SBarry Smith Vec local,global; 19847c6ae99SBarry Smith VecScatter ltog,gtol; 199ce00eea3SSatish Balay IS to,from,ltogis; 2006f951b95Secoon PetscBool twod; 20147c6ae99SBarry Smith PetscErrorCode ierr; 20247c6ae99SBarry Smith 2036f951b95Secoon 20447c6ae99SBarry Smith PetscFunctionBegin; 205c2859e5eSBarry Smith if (stencil_type == DMDA_STENCIL_BOX && (bx == DMDA_BOUNDARY_MIRROR || by == DMDA_BOUNDARY_MIRROR || bz == DMDA_BOUNDARY_MIRROR)) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Mirror boundary and box stencil"); 20647c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject) da, &comm);CHKERRQ(ierr); 2073855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES) 2083855c12bSBarry Smith if (((Petsc64bitInt) M)*((Petsc64bitInt) N)*((Petsc64bitInt) P)*((Petsc64bitInt) dof) > (Petsc64bitInt) PETSC_MPI_INT_MAX) SETERRQ3(comm,PETSC_ERR_INT_OVERFLOW,"Mesh of %D by %D by %D (dof) is too large for 32 bit indices",M,N,dof); 2093855c12bSBarry Smith #endif 2103855c12bSBarry Smith 21147c6ae99SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 21247c6ae99SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 21347c6ae99SBarry Smith 21447c6ae99SBarry Smith if (m != PETSC_DECIDE) { 21547c6ae99SBarry Smith if (m < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m); 21647c6ae99SBarry Smith else if (m > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size); 21747c6ae99SBarry Smith } 21847c6ae99SBarry Smith if (n != PETSC_DECIDE) { 21947c6ae99SBarry Smith if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n); 22047c6ae99SBarry Smith else if (n > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size); 22147c6ae99SBarry Smith } 22247c6ae99SBarry Smith if (p != PETSC_DECIDE) { 22347c6ae99SBarry Smith if (p < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Z direction: %D",p); 22447c6ae99SBarry Smith else if (p > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Z direction: %D %d",p,size); 22547c6ae99SBarry Smith } 2260716a85fSBarry Smith if ((m > 0) && (n > 0) && (p > 0) && (m*n*p != size)) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"m %D * n %D * p %D != size %d",m,n,p,size); 22747c6ae99SBarry Smith 22847c6ae99SBarry Smith /* Partition the array among the processors */ 22947c6ae99SBarry Smith if (m == PETSC_DECIDE && n != PETSC_DECIDE && p != PETSC_DECIDE) { 23047c6ae99SBarry Smith m = size/(n*p); 23147c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) { 23247c6ae99SBarry Smith n = size/(m*p); 23347c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) { 23447c6ae99SBarry Smith p = size/(m*n); 23547c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) { 23647c6ae99SBarry Smith /* try for squarish distribution */ 2378f1a2a5eSBarry Smith m = (int)(0.5 + sqrt(((double)M)*((double)size)/((double)N*p))); 23847c6ae99SBarry Smith if (!m) m = 1; 23947c6ae99SBarry Smith while (m > 0) { 24047c6ae99SBarry Smith n = size/(m*p); 24147c6ae99SBarry Smith if (m*n*p == size) break; 24247c6ae99SBarry Smith m--; 24347c6ae99SBarry Smith } 24447c6ae99SBarry Smith if (!m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad p value: p = %D",p); 24547c6ae99SBarry Smith if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;} 24647c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) { 24747c6ae99SBarry Smith /* try for squarish distribution */ 2488f1a2a5eSBarry Smith m = (int)(0.5 + sqrt(((double)M)*((double)size)/((double)P*n))); 24947c6ae99SBarry Smith if (!m) m = 1; 25047c6ae99SBarry Smith while (m > 0) { 25147c6ae99SBarry Smith p = size/(m*n); 25247c6ae99SBarry Smith if (m*n*p == size) break; 25347c6ae99SBarry Smith m--; 25447c6ae99SBarry Smith } 25547c6ae99SBarry Smith if (!m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad n value: n = %D",n); 25647c6ae99SBarry Smith if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;} 25747c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) { 25847c6ae99SBarry Smith /* try for squarish distribution */ 2598f1a2a5eSBarry Smith n = (int)(0.5 + sqrt(((double)N)*((double)size)/((double)P*m))); 26047c6ae99SBarry Smith if (!n) n = 1; 26147c6ae99SBarry Smith while (n > 0) { 26247c6ae99SBarry Smith p = size/(m*n); 26347c6ae99SBarry Smith if (m*n*p == size) break; 26447c6ae99SBarry Smith n--; 26547c6ae99SBarry Smith } 26647c6ae99SBarry Smith if (!n) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad m value: m = %D",n); 26747c6ae99SBarry Smith if (N > P && n < p) {PetscInt _n = n; n = p; p = _n;} 26847c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) { 26947c6ae99SBarry Smith /* try for squarish distribution */ 2708f1a2a5eSBarry Smith n = (PetscInt)(0.5 + pow(((double)N*N)*((double)size)/((double)P*M),(double)(1./3.))); 27147c6ae99SBarry Smith if (!n) n = 1; 27247c6ae99SBarry Smith while (n > 0) { 27347c6ae99SBarry Smith pm = size/n; 27447c6ae99SBarry Smith if (n*pm == size) break; 27547c6ae99SBarry Smith n--; 27647c6ae99SBarry Smith } 27747c6ae99SBarry Smith if (!n) n = 1; 2788f1a2a5eSBarry Smith m = (PetscInt)(0.5 + sqrt(((double)M)*((double)size)/((double)P*n))); 27947c6ae99SBarry Smith if (!m) m = 1; 28047c6ae99SBarry Smith while (m > 0) { 28147c6ae99SBarry Smith p = size/(m*n); 28247c6ae99SBarry Smith if (m*n*p == size) break; 28347c6ae99SBarry Smith m--; 28447c6ae99SBarry Smith } 28547c6ae99SBarry Smith if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;} 28647c6ae99SBarry Smith } else if (m*n*p != size) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition"); 28747c6ae99SBarry Smith 28847c6ae99SBarry Smith if (m*n*p != size) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_PLIB,"Could not find good partition"); 28947c6ae99SBarry Smith if (M < m) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m); 29047c6ae99SBarry Smith if (N < n) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n); 29147c6ae99SBarry Smith if (P < p) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in z direction is too fine! %D %D",P,p); 29247c6ae99SBarry Smith 29347c6ae99SBarry Smith /* 29447c6ae99SBarry Smith Determine locally owned region 29547c6ae99SBarry Smith [x, y, or z]s is the first local node number, [x, y, z] is the number of local nodes 29647c6ae99SBarry Smith */ 29747c6ae99SBarry Smith 29847c6ae99SBarry Smith if (!lx) { 29947c6ae99SBarry Smith ierr = PetscMalloc(m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr); 30047c6ae99SBarry Smith lx = dd->lx; 30147c6ae99SBarry Smith for (i=0; i<m; i++) { 30247c6ae99SBarry Smith lx[i] = M/m + ((M % m) > (i % m)); 30347c6ae99SBarry Smith } 30447c6ae99SBarry Smith } 30547c6ae99SBarry Smith x = lx[rank % m]; 30647c6ae99SBarry Smith xs = 0; 30747c6ae99SBarry Smith for (i=0; i<(rank%m); i++) { xs += lx[i];} 308*d9c9ebe5SPeter Brune if ((x < s) && ((m > 1) || (bx == DMDA_BOUNDARY_PERIODIC))) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local x-width of domain x %D is smaller than stencil width s %D",x,s); 30947c6ae99SBarry Smith 31047c6ae99SBarry Smith if (!ly) { 31147c6ae99SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr); 31247c6ae99SBarry Smith ly = dd->ly; 31347c6ae99SBarry Smith for (i=0; i<n; i++) { 31447c6ae99SBarry Smith ly[i] = N/n + ((N % n) > (i % n)); 31547c6ae99SBarry Smith } 31647c6ae99SBarry Smith } 31747c6ae99SBarry Smith y = ly[(rank % (m*n))/m]; 318*d9c9ebe5SPeter Brune if ((y < s) && ((n > 1) || (by == DMDA_BOUNDARY_PERIODIC))) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local y-width of domain y %D is smaller than stencil width s %D",y,s); 31930729d88SBarry Smith 32047c6ae99SBarry Smith ys = 0; 32147c6ae99SBarry Smith for (i=0; i<(rank % (m*n))/m; i++) { ys += ly[i];} 32247c6ae99SBarry Smith 32347c6ae99SBarry Smith if (!lz) { 32447c6ae99SBarry Smith ierr = PetscMalloc(p*sizeof(PetscInt), &dd->lz);CHKERRQ(ierr); 32547c6ae99SBarry Smith lz = dd->lz; 32647c6ae99SBarry Smith for (i=0; i<p; i++) { 32747c6ae99SBarry Smith lz[i] = P/p + ((P % p) > (i % p)); 32847c6ae99SBarry Smith } 32947c6ae99SBarry Smith } 33047c6ae99SBarry Smith z = lz[rank/(m*n)]; 331bcea557cSEthan Coon 332fdc81ce8SEthan Coon /* note this is different than x- and y-, as we will handle as an important special 333fdc81ce8SEthan Coon case when p=P=1 and DMDA_BOUNDARY_PERIODIC and s > z. This is to deal with 2D problems 334fdc81ce8SEthan Coon in a 3D code. Additional code for this case is noted with "2d case" comments */ 3356f951b95Secoon twod = PETSC_FALSE; 3366f951b95Secoon if (P == 1) { 3376f951b95Secoon twod = PETSC_TRUE; 338*d9c9ebe5SPeter Brune } else if ((z < s) && ((p > 1) || (bz == DMDA_BOUNDARY_PERIODIC))) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local z-width of domain z %D is smaller than stencil width s %D",z,s); 33947c6ae99SBarry Smith zs = 0; 34047c6ae99SBarry Smith for (i=0; i<(rank/(m*n)); i++) { zs += lz[i];} 34147c6ae99SBarry Smith ye = ys + y; 34247c6ae99SBarry Smith xe = xs + x; 34347c6ae99SBarry Smith ze = zs + z; 34447c6ae99SBarry Smith 34588661749SPeter Brune /* determine ghost region (Xs) and region scattered into (IXs) */ 346*d9c9ebe5SPeter Brune if (xs-s > 0) { 347*d9c9ebe5SPeter Brune Xs = xs - s; IXs = xs - s; 34888661749SPeter Brune } else { 34988661749SPeter Brune if (bx) { 35088661749SPeter Brune Xs = xs - s; 35188661749SPeter Brune } else { 35288661749SPeter Brune Xs = 0; 35388661749SPeter Brune } 35488661749SPeter Brune IXs = 0; 35588661749SPeter Brune } 356*d9c9ebe5SPeter Brune if (xe+s <= M) { 357*d9c9ebe5SPeter Brune Xe = xe + s; IXe = xe + s; 35888661749SPeter Brune } else { 35988661749SPeter Brune if (bx) { 360*d9c9ebe5SPeter Brune Xs = xs - s; Xe = xe + s; 36188661749SPeter Brune } else { 36288661749SPeter Brune Xe = M; 36388661749SPeter Brune } 36488661749SPeter Brune IXe = M; 36588661749SPeter Brune } 36647c6ae99SBarry Smith 367c2859e5eSBarry Smith if (bx == DMDA_BOUNDARY_PERIODIC || bx == DMDA_BOUNDARY_MIRROR) { 368*d9c9ebe5SPeter Brune IXs = xs - s; 369*d9c9ebe5SPeter Brune IXe = xe + s; 370*d9c9ebe5SPeter Brune Xs = xs - s; 371*d9c9ebe5SPeter Brune Xe = xe + s; 37288661749SPeter Brune } 37388661749SPeter Brune 374*d9c9ebe5SPeter Brune if (ys-s > 0) { 375*d9c9ebe5SPeter Brune Ys = ys - s; IYs = ys - s; 37688661749SPeter Brune } else { 37788661749SPeter Brune if (by) { 37888661749SPeter Brune Ys = ys - s; 37988661749SPeter Brune } else { 38088661749SPeter Brune Ys = 0; 38188661749SPeter Brune } 38288661749SPeter Brune IYs = 0; 38388661749SPeter Brune } 384*d9c9ebe5SPeter Brune if (ye+s <= N) { 385*d9c9ebe5SPeter Brune Ye = ye + s; IYe = ye + s; 38688661749SPeter Brune } else { 38788661749SPeter Brune if (by) { 38888661749SPeter Brune Ye = ye + s; 38988661749SPeter Brune } else { 39088661749SPeter Brune Ye = N; 39188661749SPeter Brune } 39288661749SPeter Brune IYe = N; 39388661749SPeter Brune } 39488661749SPeter Brune 395c2859e5eSBarry Smith if (by == DMDA_BOUNDARY_PERIODIC || by == DMDA_BOUNDARY_MIRROR) { 396*d9c9ebe5SPeter Brune IYs = ys - s; 397*d9c9ebe5SPeter Brune IYe = ye + s; 398*d9c9ebe5SPeter Brune Ys = ys - s; 399*d9c9ebe5SPeter Brune Ye = ye + s; 40088661749SPeter Brune } 40188661749SPeter Brune 402*d9c9ebe5SPeter Brune if (zs-s > 0) { 403*d9c9ebe5SPeter Brune Zs = zs - s; IZs = zs - s; 40488661749SPeter Brune } else { 40588661749SPeter Brune if (bz) { 40688661749SPeter Brune Zs = zs - s; 40788661749SPeter Brune } else { 40888661749SPeter Brune Zs = 0; 40988661749SPeter Brune } 41088661749SPeter Brune IZs = 0; 41188661749SPeter Brune } 412*d9c9ebe5SPeter Brune if (ze+s <= P) { 413*d9c9ebe5SPeter Brune Ze = ze + s; IZe = ze + s; 41488661749SPeter Brune } else { 41588661749SPeter Brune if (bz) { 41688661749SPeter Brune Ze = ze + s; 41788661749SPeter Brune } else { 41888661749SPeter Brune Ze = P; 41988661749SPeter Brune } 42088661749SPeter Brune IZe = P; 42188661749SPeter Brune } 42288661749SPeter Brune 423c2859e5eSBarry Smith if (bz == DMDA_BOUNDARY_PERIODIC || bz == DMDA_BOUNDARY_MIRROR) { 424*d9c9ebe5SPeter Brune IZs = zs - s; 425*d9c9ebe5SPeter Brune IZe = ze + s; 426*d9c9ebe5SPeter Brune Zs = zs - s; 427*d9c9ebe5SPeter Brune Ze = ze + s; 42888661749SPeter Brune } 42947c6ae99SBarry Smith 43047c6ae99SBarry Smith /* Resize all X parameters to reflect w */ 431*d9c9ebe5SPeter Brune s_x = s; 432*d9c9ebe5SPeter Brune s_y = s; 433*d9c9ebe5SPeter Brune s_z = s; 43447c6ae99SBarry Smith 43547c6ae99SBarry Smith /* determine starting point of each processor */ 43647c6ae99SBarry Smith nn = x*y*z; 43747c6ae99SBarry Smith ierr = PetscMalloc2(size+1,PetscInt,&bases,size,PetscInt,&ldims);CHKERRQ(ierr); 43847c6ae99SBarry Smith ierr = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr); 43947c6ae99SBarry Smith bases[0] = 0; 44047c6ae99SBarry Smith for (i=1; i<=size; i++) { 44147c6ae99SBarry Smith bases[i] = ldims[i-1]; 44247c6ae99SBarry Smith } 44347c6ae99SBarry Smith for (i=1; i<=size; i++) { 44447c6ae99SBarry Smith bases[i] += bases[i-1]; 44547c6ae99SBarry Smith } 446ce00eea3SSatish Balay base = bases[rank]*dof; 44747c6ae99SBarry Smith 44847c6ae99SBarry Smith /* allocate the base parallel and sequential vectors */ 449ce00eea3SSatish Balay dd->Nlocal = x*y*z*dof; 450778a2246SBarry Smith ierr = VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,0,&global);CHKERRQ(ierr); 451ce00eea3SSatish Balay dd->nlocal = (Xe-Xs)*(Ye-Ys)*(Ze-Zs)*dof; 452778a2246SBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,0,&local);CHKERRQ(ierr); 45347c6ae99SBarry Smith 45447c6ae99SBarry Smith /* generate appropriate vector scatters */ 45547c6ae99SBarry Smith /* local to global inserts non-ghost point region into global */ 45647c6ae99SBarry Smith ierr = VecGetOwnershipRange(global,&start,&end);CHKERRQ(ierr); 457ce00eea3SSatish Balay ierr = ISCreateStride(comm,x*y*z*dof,start,1,&to);CHKERRQ(ierr); 45847c6ae99SBarry Smith 459ce00eea3SSatish Balay ierr = PetscMalloc(x*y*z*sizeof(PetscInt),&idx);CHKERRQ(ierr); 460ce00eea3SSatish Balay left = xs - Xs; right = left + x; 46147c6ae99SBarry Smith bottom = ys - Ys; top = bottom + y; 46247c6ae99SBarry Smith down = zs - Zs; up = down + z; 46347c6ae99SBarry Smith count = 0; 46447c6ae99SBarry Smith for (i=down; i<up; i++) { 46547c6ae99SBarry Smith for (j=bottom; j<top; j++) { 466ce00eea3SSatish Balay for (k=left; k<right; k++) { 467ce00eea3SSatish Balay idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 46847c6ae99SBarry Smith } 46947c6ae99SBarry Smith } 47047c6ae99SBarry Smith } 47147c6ae99SBarry Smith 47247c6ae99SBarry Smith ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&from);CHKERRQ(ierr); 47347c6ae99SBarry Smith ierr = VecScatterCreate(local,from,global,to,<og);CHKERRQ(ierr); 47447c6ae99SBarry Smith ierr = PetscLogObjectParent(da,ltog);CHKERRQ(ierr); 475fcfd50ebSBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 476fcfd50ebSBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 47747c6ae99SBarry Smith 478ce00eea3SSatish Balay /* global to local must include ghost points within the domain, 479ce00eea3SSatish Balay but not ghost points outside the domain that aren't periodic */ 480*d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_BOX) { 481db87c5ecSEthan Coon count = (IXe-IXs)*(IYe-IYs)*(IZe-IZs); 482db87c5ecSEthan Coon ierr = PetscMalloc(count*sizeof(PetscInt),&idx);CHKERRQ(ierr); 483ce00eea3SSatish Balay 484ce00eea3SSatish Balay left = IXs - Xs; right = left + (IXe-IXs); 485ce00eea3SSatish Balay bottom = IYs - Ys; top = bottom + (IYe-IYs); 486ce00eea3SSatish Balay down = IZs - Zs; up = down + (IZe-IZs); 487ce00eea3SSatish Balay count = 0; 488ce00eea3SSatish Balay for (i=down; i<up; i++) { 489ce00eea3SSatish Balay for (j=bottom; j<top; j++) { 490ce00eea3SSatish Balay for (k=left; k<right; k++) { 491ce00eea3SSatish Balay idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 492ce00eea3SSatish Balay } 493ce00eea3SSatish Balay } 494ce00eea3SSatish Balay } 495ce00eea3SSatish Balay ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr); 496ce00eea3SSatish Balay 49747c6ae99SBarry Smith } else { 49847c6ae99SBarry Smith /* This is way ugly! We need to list the funny cross type region */ 499db87c5ecSEthan Coon count = ((ys-IYs) + (IYe-ye))*x*z + ((xs-IXs) + (IXe-xe))*y*z + ((zs-IZs) + (IZe-ze))*x*y + x*y*z; 500db87c5ecSEthan Coon ierr = PetscMalloc(count*sizeof(PetscInt),&idx);CHKERRQ(ierr); 501ce00eea3SSatish Balay 502ce00eea3SSatish Balay left = xs - Xs; right = left + x; 50347c6ae99SBarry Smith bottom = ys - Ys; top = bottom + y; 50447c6ae99SBarry Smith down = zs - Zs; up = down + z; 50547c6ae99SBarry Smith count = 0; 506ce00eea3SSatish Balay /* the bottom chunck */ 507ce00eea3SSatish Balay for (i=(IZs-Zs); i<down; i++) { 50847c6ae99SBarry Smith for (j=bottom; j<top; j++) { 509ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 51047c6ae99SBarry Smith } 51147c6ae99SBarry Smith } 51247c6ae99SBarry Smith /* the middle piece */ 51347c6ae99SBarry Smith for (i=down; i<up; i++) { 51447c6ae99SBarry Smith /* front */ 515ce00eea3SSatish Balay for (j=(IYs-Ys); j<bottom; j++) { 516ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 51747c6ae99SBarry Smith } 51847c6ae99SBarry Smith /* middle */ 51947c6ae99SBarry Smith for (j=bottom; j<top; j++) { 520ce00eea3SSatish Balay for (k=IXs-Xs; k<IXe-Xs; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 52147c6ae99SBarry Smith } 52247c6ae99SBarry Smith /* back */ 523ce00eea3SSatish Balay for (j=top; j<top+IYe-ye; j++) { 524ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 52547c6ae99SBarry Smith } 52647c6ae99SBarry Smith } 52747c6ae99SBarry Smith /* the top piece */ 528ce00eea3SSatish Balay for (i=up; i<up+IZe-ze; i++) { 52947c6ae99SBarry Smith for (j=bottom; j<top; j++) { 530ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 53147c6ae99SBarry Smith } 53247c6ae99SBarry Smith } 53347c6ae99SBarry Smith ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr); 53447c6ae99SBarry Smith } 53547c6ae99SBarry Smith 53647c6ae99SBarry Smith /* determine who lies on each side of use stored in n24 n25 n26 53747c6ae99SBarry Smith n21 n22 n23 53847c6ae99SBarry Smith n18 n19 n20 53947c6ae99SBarry Smith 54047c6ae99SBarry Smith n15 n16 n17 54147c6ae99SBarry Smith n12 n14 54247c6ae99SBarry Smith n9 n10 n11 54347c6ae99SBarry Smith 54447c6ae99SBarry Smith n6 n7 n8 54547c6ae99SBarry Smith n3 n4 n5 54647c6ae99SBarry Smith n0 n1 n2 54747c6ae99SBarry Smith */ 54847c6ae99SBarry Smith 54947c6ae99SBarry Smith /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */ 55047c6ae99SBarry Smith /* Assume Nodes are Internal to the Cube */ 55147c6ae99SBarry Smith n0 = rank - m*n - m - 1; 55247c6ae99SBarry Smith n1 = rank - m*n - m; 55347c6ae99SBarry Smith n2 = rank - m*n - m + 1; 55447c6ae99SBarry Smith n3 = rank - m*n -1; 55547c6ae99SBarry Smith n4 = rank - m*n; 55647c6ae99SBarry Smith n5 = rank - m*n + 1; 55747c6ae99SBarry Smith n6 = rank - m*n + m - 1; 55847c6ae99SBarry Smith n7 = rank - m*n + m; 55947c6ae99SBarry Smith n8 = rank - m*n + m + 1; 56047c6ae99SBarry Smith 56147c6ae99SBarry Smith n9 = rank - m - 1; 56247c6ae99SBarry Smith n10 = rank - m; 56347c6ae99SBarry Smith n11 = rank - m + 1; 56447c6ae99SBarry Smith n12 = rank - 1; 56547c6ae99SBarry Smith n14 = rank + 1; 56647c6ae99SBarry Smith n15 = rank + m - 1; 56747c6ae99SBarry Smith n16 = rank + m; 56847c6ae99SBarry Smith n17 = rank + m + 1; 56947c6ae99SBarry Smith 57047c6ae99SBarry Smith n18 = rank + m*n - m - 1; 57147c6ae99SBarry Smith n19 = rank + m*n - m; 57247c6ae99SBarry Smith n20 = rank + m*n - m + 1; 57347c6ae99SBarry Smith n21 = rank + m*n - 1; 57447c6ae99SBarry Smith n22 = rank + m*n; 57547c6ae99SBarry Smith n23 = rank + m*n + 1; 57647c6ae99SBarry Smith n24 = rank + m*n + m - 1; 57747c6ae99SBarry Smith n25 = rank + m*n + m; 57847c6ae99SBarry Smith n26 = rank + m*n + m + 1; 57947c6ae99SBarry Smith 58047c6ae99SBarry Smith /* Assume Pieces are on Faces of Cube */ 58147c6ae99SBarry Smith 58247c6ae99SBarry Smith if (xs == 0) { /* First assume not corner or edge */ 58347c6ae99SBarry Smith n0 = rank -1 - (m*n); 58447c6ae99SBarry Smith n3 = rank + m -1 - (m*n); 58547c6ae99SBarry Smith n6 = rank + 2*m -1 - (m*n); 58647c6ae99SBarry Smith n9 = rank -1; 58747c6ae99SBarry Smith n12 = rank + m -1; 58847c6ae99SBarry Smith n15 = rank + 2*m -1; 58947c6ae99SBarry Smith n18 = rank -1 + (m*n); 59047c6ae99SBarry Smith n21 = rank + m -1 + (m*n); 59147c6ae99SBarry Smith n24 = rank + 2*m -1 + (m*n); 59247c6ae99SBarry Smith } 59347c6ae99SBarry Smith 594ce00eea3SSatish Balay if (xe == M) { /* First assume not corner or edge */ 59547c6ae99SBarry Smith n2 = rank -2*m +1 - (m*n); 59647c6ae99SBarry Smith n5 = rank - m +1 - (m*n); 59747c6ae99SBarry Smith n8 = rank +1 - (m*n); 59847c6ae99SBarry Smith n11 = rank -2*m +1; 59947c6ae99SBarry Smith n14 = rank - m +1; 60047c6ae99SBarry Smith n17 = rank +1; 60147c6ae99SBarry Smith n20 = rank -2*m +1 + (m*n); 60247c6ae99SBarry Smith n23 = rank - m +1 + (m*n); 60347c6ae99SBarry Smith n26 = rank +1 + (m*n); 60447c6ae99SBarry Smith } 60547c6ae99SBarry Smith 60647c6ae99SBarry Smith if (ys==0) { /* First assume not corner or edge */ 60747c6ae99SBarry Smith n0 = rank + m * (n-1) -1 - (m*n); 60847c6ae99SBarry Smith n1 = rank + m * (n-1) - (m*n); 60947c6ae99SBarry Smith n2 = rank + m * (n-1) +1 - (m*n); 61047c6ae99SBarry Smith n9 = rank + m * (n-1) -1; 61147c6ae99SBarry Smith n10 = rank + m * (n-1); 61247c6ae99SBarry Smith n11 = rank + m * (n-1) +1; 61347c6ae99SBarry Smith n18 = rank + m * (n-1) -1 + (m*n); 61447c6ae99SBarry Smith n19 = rank + m * (n-1) + (m*n); 61547c6ae99SBarry Smith n20 = rank + m * (n-1) +1 + (m*n); 61647c6ae99SBarry Smith } 61747c6ae99SBarry Smith 61847c6ae99SBarry Smith if (ye == N) { /* First assume not corner or edge */ 61947c6ae99SBarry Smith n6 = rank - m * (n-1) -1 - (m*n); 62047c6ae99SBarry Smith n7 = rank - m * (n-1) - (m*n); 62147c6ae99SBarry Smith n8 = rank - m * (n-1) +1 - (m*n); 62247c6ae99SBarry Smith n15 = rank - m * (n-1) -1; 62347c6ae99SBarry Smith n16 = rank - m * (n-1); 62447c6ae99SBarry Smith n17 = rank - m * (n-1) +1; 62547c6ae99SBarry Smith n24 = rank - m * (n-1) -1 + (m*n); 62647c6ae99SBarry Smith n25 = rank - m * (n-1) + (m*n); 62747c6ae99SBarry Smith n26 = rank - m * (n-1) +1 + (m*n); 62847c6ae99SBarry Smith } 62947c6ae99SBarry Smith 63047c6ae99SBarry Smith if (zs == 0) { /* First assume not corner or edge */ 63147c6ae99SBarry Smith n0 = size - (m*n) + rank - m - 1; 63247c6ae99SBarry Smith n1 = size - (m*n) + rank - m; 63347c6ae99SBarry Smith n2 = size - (m*n) + rank - m + 1; 63447c6ae99SBarry Smith n3 = size - (m*n) + rank - 1; 63547c6ae99SBarry Smith n4 = size - (m*n) + rank; 63647c6ae99SBarry Smith n5 = size - (m*n) + rank + 1; 63747c6ae99SBarry Smith n6 = size - (m*n) + rank + m - 1; 63847c6ae99SBarry Smith n7 = size - (m*n) + rank + m ; 63947c6ae99SBarry Smith n8 = size - (m*n) + rank + m + 1; 64047c6ae99SBarry Smith } 64147c6ae99SBarry Smith 64247c6ae99SBarry Smith if (ze == P) { /* First assume not corner or edge */ 64347c6ae99SBarry Smith n18 = (m*n) - (size-rank) - m - 1; 64447c6ae99SBarry Smith n19 = (m*n) - (size-rank) - m; 64547c6ae99SBarry Smith n20 = (m*n) - (size-rank) - m + 1; 64647c6ae99SBarry Smith n21 = (m*n) - (size-rank) - 1; 64747c6ae99SBarry Smith n22 = (m*n) - (size-rank); 64847c6ae99SBarry Smith n23 = (m*n) - (size-rank) + 1; 64947c6ae99SBarry Smith n24 = (m*n) - (size-rank) + m - 1; 65047c6ae99SBarry Smith n25 = (m*n) - (size-rank) + m; 65147c6ae99SBarry Smith n26 = (m*n) - (size-rank) + m + 1; 65247c6ae99SBarry Smith } 65347c6ae99SBarry Smith 65447c6ae99SBarry Smith if ((xs==0) && (zs==0)) { /* Assume an edge, not corner */ 65547c6ae99SBarry Smith n0 = size - m*n + rank + m-1 - m; 65647c6ae99SBarry Smith n3 = size - m*n + rank + m-1; 65747c6ae99SBarry Smith n6 = size - m*n + rank + m-1 + m; 65847c6ae99SBarry Smith } 65947c6ae99SBarry Smith 66047c6ae99SBarry Smith if ((xs==0) && (ze==P)) { /* Assume an edge, not corner */ 66147c6ae99SBarry Smith n18 = m*n - (size - rank) + m-1 - m; 66247c6ae99SBarry Smith n21 = m*n - (size - rank) + m-1; 66347c6ae99SBarry Smith n24 = m*n - (size - rank) + m-1 + m; 66447c6ae99SBarry Smith } 66547c6ae99SBarry Smith 66647c6ae99SBarry Smith if ((xs==0) && (ys==0)) { /* Assume an edge, not corner */ 66747c6ae99SBarry Smith n0 = rank + m*n -1 - m*n; 66847c6ae99SBarry Smith n9 = rank + m*n -1; 66947c6ae99SBarry Smith n18 = rank + m*n -1 + m*n; 67047c6ae99SBarry Smith } 67147c6ae99SBarry Smith 67247c6ae99SBarry Smith if ((xs==0) && (ye==N)) { /* Assume an edge, not corner */ 67347c6ae99SBarry Smith n6 = rank - m*(n-1) + m-1 - m*n; 67447c6ae99SBarry Smith n15 = rank - m*(n-1) + m-1; 67547c6ae99SBarry Smith n24 = rank - m*(n-1) + m-1 + m*n; 67647c6ae99SBarry Smith } 67747c6ae99SBarry Smith 678ce00eea3SSatish Balay if ((xe==M) && (zs==0)) { /* Assume an edge, not corner */ 67947c6ae99SBarry Smith n2 = size - (m*n-rank) - (m-1) - m; 68047c6ae99SBarry Smith n5 = size - (m*n-rank) - (m-1); 68147c6ae99SBarry Smith n8 = size - (m*n-rank) - (m-1) + m; 68247c6ae99SBarry Smith } 68347c6ae99SBarry Smith 684ce00eea3SSatish Balay if ((xe==M) && (ze==P)) { /* Assume an edge, not corner */ 68547c6ae99SBarry Smith n20 = m*n - (size - rank) - (m-1) - m; 68647c6ae99SBarry Smith n23 = m*n - (size - rank) - (m-1); 68747c6ae99SBarry Smith n26 = m*n - (size - rank) - (m-1) + m; 68847c6ae99SBarry Smith } 68947c6ae99SBarry Smith 690ce00eea3SSatish Balay if ((xe==M) && (ys==0)) { /* Assume an edge, not corner */ 69147c6ae99SBarry Smith n2 = rank + m*(n-1) - (m-1) - m*n; 69247c6ae99SBarry Smith n11 = rank + m*(n-1) - (m-1); 69347c6ae99SBarry Smith n20 = rank + m*(n-1) - (m-1) + m*n; 69447c6ae99SBarry Smith } 69547c6ae99SBarry Smith 696ce00eea3SSatish Balay if ((xe==M) && (ye==N)) { /* Assume an edge, not corner */ 69747c6ae99SBarry Smith n8 = rank - m*n +1 - m*n; 69847c6ae99SBarry Smith n17 = rank - m*n +1; 69947c6ae99SBarry Smith n26 = rank - m*n +1 + m*n; 70047c6ae99SBarry Smith } 70147c6ae99SBarry Smith 70247c6ae99SBarry Smith if ((ys==0) && (zs==0)) { /* Assume an edge, not corner */ 70347c6ae99SBarry Smith n0 = size - m + rank -1; 70447c6ae99SBarry Smith n1 = size - m + rank; 70547c6ae99SBarry Smith n2 = size - m + rank +1; 70647c6ae99SBarry Smith } 70747c6ae99SBarry Smith 70847c6ae99SBarry Smith if ((ys==0) && (ze==P)) { /* Assume an edge, not corner */ 70947c6ae99SBarry Smith n18 = m*n - (size - rank) + m*(n-1) -1; 71047c6ae99SBarry Smith n19 = m*n - (size - rank) + m*(n-1); 71147c6ae99SBarry Smith n20 = m*n - (size - rank) + m*(n-1) +1; 71247c6ae99SBarry Smith } 71347c6ae99SBarry Smith 71447c6ae99SBarry Smith if ((ye==N) && (zs==0)) { /* Assume an edge, not corner */ 71547c6ae99SBarry Smith n6 = size - (m*n-rank) - m * (n-1) -1; 71647c6ae99SBarry Smith n7 = size - (m*n-rank) - m * (n-1); 71747c6ae99SBarry Smith n8 = size - (m*n-rank) - m * (n-1) +1; 71847c6ae99SBarry Smith } 71947c6ae99SBarry Smith 72047c6ae99SBarry Smith if ((ye==N) && (ze==P)) { /* Assume an edge, not corner */ 72147c6ae99SBarry Smith n24 = rank - (size-m) -1; 72247c6ae99SBarry Smith n25 = rank - (size-m); 72347c6ae99SBarry Smith n26 = rank - (size-m) +1; 72447c6ae99SBarry Smith } 72547c6ae99SBarry Smith 72647c6ae99SBarry Smith /* Check for Corners */ 72747c6ae99SBarry Smith if ((xs==0) && (ys==0) && (zs==0)) { n0 = size -1;} 72847c6ae99SBarry Smith if ((xs==0) && (ys==0) && (ze==P)) { n18 = m*n-1;} 72947c6ae99SBarry Smith if ((xs==0) && (ye==N) && (zs==0)) { n6 = (size-1)-m*(n-1);} 73047c6ae99SBarry Smith if ((xs==0) && (ye==N) && (ze==P)) { n24 = m-1;} 731ce00eea3SSatish Balay if ((xe==M) && (ys==0) && (zs==0)) { n2 = size-m;} 732ce00eea3SSatish Balay if ((xe==M) && (ys==0) && (ze==P)) { n20 = m*n-m;} 733ce00eea3SSatish Balay if ((xe==M) && (ye==N) && (zs==0)) { n8 = size-m*n;} 734ce00eea3SSatish Balay if ((xe==M) && (ye==N) && (ze==P)) { n26 = 0;} 73547c6ae99SBarry Smith 73647c6ae99SBarry Smith /* Check for when not X,Y, and Z Periodic */ 73747c6ae99SBarry Smith 73847c6ae99SBarry Smith /* If not X periodic */ 7391321219cSEthan Coon if (bx != DMDA_BOUNDARY_PERIODIC) { 74047c6ae99SBarry Smith if (xs==0) {n0 = n3 = n6 = n9 = n12 = n15 = n18 = n21 = n24 = -2;} 741ce00eea3SSatish Balay if (xe==M) {n2 = n5 = n8 = n11 = n14 = n17 = n20 = n23 = n26 = -2;} 74247c6ae99SBarry Smith } 74347c6ae99SBarry Smith 74447c6ae99SBarry Smith /* If not Y periodic */ 7451321219cSEthan Coon if (by != DMDA_BOUNDARY_PERIODIC) { 74647c6ae99SBarry Smith if (ys==0) {n0 = n1 = n2 = n9 = n10 = n11 = n18 = n19 = n20 = -2;} 74747c6ae99SBarry Smith if (ye==N) {n6 = n7 = n8 = n15 = n16 = n17 = n24 = n25 = n26 = -2;} 74847c6ae99SBarry Smith } 74947c6ae99SBarry Smith 75047c6ae99SBarry Smith /* If not Z periodic */ 7511321219cSEthan Coon if (bz != DMDA_BOUNDARY_PERIODIC) { 75247c6ae99SBarry Smith if (zs==0) {n0 = n1 = n2 = n3 = n4 = n5 = n6 = n7 = n8 = -2;} 75347c6ae99SBarry Smith if (ze==P) {n18 = n19 = n20 = n21 = n22 = n23 = n24 = n25 = n26 = -2;} 75447c6ae99SBarry Smith } 75547c6ae99SBarry Smith 75647c6ae99SBarry Smith ierr = PetscMalloc(27*sizeof(PetscInt),&dd->neighbors);CHKERRQ(ierr); 75747c6ae99SBarry Smith dd->neighbors[0] = n0; 75847c6ae99SBarry Smith dd->neighbors[1] = n1; 75947c6ae99SBarry Smith dd->neighbors[2] = n2; 76047c6ae99SBarry Smith dd->neighbors[3] = n3; 76147c6ae99SBarry Smith dd->neighbors[4] = n4; 76247c6ae99SBarry Smith dd->neighbors[5] = n5; 76347c6ae99SBarry Smith dd->neighbors[6] = n6; 76447c6ae99SBarry Smith dd->neighbors[7] = n7; 76547c6ae99SBarry Smith dd->neighbors[8] = n8; 76647c6ae99SBarry Smith dd->neighbors[9] = n9; 76747c6ae99SBarry Smith dd->neighbors[10] = n10; 76847c6ae99SBarry Smith dd->neighbors[11] = n11; 76947c6ae99SBarry Smith dd->neighbors[12] = n12; 77047c6ae99SBarry Smith dd->neighbors[13] = rank; 77147c6ae99SBarry Smith dd->neighbors[14] = n14; 77247c6ae99SBarry Smith dd->neighbors[15] = n15; 77347c6ae99SBarry Smith dd->neighbors[16] = n16; 77447c6ae99SBarry Smith dd->neighbors[17] = n17; 77547c6ae99SBarry Smith dd->neighbors[18] = n18; 77647c6ae99SBarry Smith dd->neighbors[19] = n19; 77747c6ae99SBarry Smith dd->neighbors[20] = n20; 77847c6ae99SBarry Smith dd->neighbors[21] = n21; 77947c6ae99SBarry Smith dd->neighbors[22] = n22; 78047c6ae99SBarry Smith dd->neighbors[23] = n23; 78147c6ae99SBarry Smith dd->neighbors[24] = n24; 78247c6ae99SBarry Smith dd->neighbors[25] = n25; 78347c6ae99SBarry Smith dd->neighbors[26] = n26; 78447c6ae99SBarry Smith 78547c6ae99SBarry Smith /* If star stencil then delete the corner neighbors */ 786*d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 78747c6ae99SBarry Smith /* save information about corner neighbors */ 78847c6ae99SBarry Smith sn0 = n0; sn1 = n1; sn2 = n2; sn3 = n3; sn5 = n5; sn6 = n6; sn7 = n7; 78947c6ae99SBarry Smith sn8 = n8; sn9 = n9; sn11 = n11; sn15 = n15; sn17 = n17; sn18 = n18; 79047c6ae99SBarry Smith sn19 = n19; sn20 = n20; sn21 = n21; sn23 = n23; sn24 = n24; sn25 = n25; 79147c6ae99SBarry Smith sn26 = n26; 79247c6ae99SBarry Smith n0 = n1 = n2 = n3 = n5 = n6 = n7 = n8 = n9 = n11 = 79347c6ae99SBarry Smith n15 = n17 = n18 = n19 = n20 = n21 = n23 = n24 = n25 = n26 = -1; 79447c6ae99SBarry Smith } 79547c6ae99SBarry Smith 79647c6ae99SBarry Smith 79747c6ae99SBarry Smith ierr = PetscMalloc((Xe-Xs)*(Ye-Ys)*(Ze-Zs)*sizeof(PetscInt),&idx);CHKERRQ(ierr); 79847c6ae99SBarry Smith ierr = PetscLogObjectMemory(da,(Xe-Xs)*(Ye-Ys)*(Ze-Zs)*sizeof(PetscInt));CHKERRQ(ierr); 79947c6ae99SBarry Smith 80047c6ae99SBarry Smith nn = 0; 80147c6ae99SBarry Smith /* Bottom Level */ 80247c6ae99SBarry Smith for (k=0; k<s_z; k++) { 80347c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 80447c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 805ce00eea3SSatish Balay x_t = lx[n0 % m]; 80647c6ae99SBarry Smith y_t = ly[(n0 % (m*n))/m]; 80747c6ae99SBarry Smith z_t = lz[n0 / (m*n)]; 80847c6ae99SBarry 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; 8096f951b95Secoon if (twod && (s_t < 0)) {s_t = bases[n0] + x_t*y_t*z_t - (s_y-i)*x_t - s_x;} /* 2D case */ 81047c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 81147c6ae99SBarry Smith } 81247c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 81347c6ae99SBarry Smith x_t = x; 81447c6ae99SBarry Smith y_t = ly[(n1 % (m*n))/m]; 81547c6ae99SBarry Smith z_t = lz[n1 / (m*n)]; 81647c6ae99SBarry 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; 8176f951b95Secoon if (twod && (s_t < 0)) {s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t;} /* 2D case */ 81847c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 81947c6ae99SBarry Smith } 82047c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 821ce00eea3SSatish Balay x_t = lx[n2 % m]; 82247c6ae99SBarry Smith y_t = ly[(n2 % (m*n))/m]; 82347c6ae99SBarry Smith z_t = lz[n2 / (m*n)]; 82447c6ae99SBarry 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; 8256f951b95Secoon if (twod && (s_t < 0)) {s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t;} /* 2D case */ 82647c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 82747c6ae99SBarry Smith } 82847c6ae99SBarry Smith } 82947c6ae99SBarry Smith 83047c6ae99SBarry Smith for (i=0; i<y; i++) { 83147c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 832ce00eea3SSatish Balay x_t = lx[n3 % m]; 83347c6ae99SBarry Smith y_t = y; 83447c6ae99SBarry Smith z_t = lz[n3 / (m*n)]; 83547c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8366f951b95Secoon 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 */ 83747c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 83847c6ae99SBarry Smith } 83947c6ae99SBarry Smith 84047c6ae99SBarry Smith if (n4 >= 0) { /* middle */ 84147c6ae99SBarry Smith x_t = x; 84247c6ae99SBarry Smith y_t = y; 84347c6ae99SBarry Smith z_t = lz[n4 / (m*n)]; 84447c6ae99SBarry Smith s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8456f951b95Secoon if (twod && (s_t < 0)) {s_t = bases[n4] + i*x_t + x_t*y_t*z_t - x_t*y_t;} /* 2D case */ 84647c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 847c2859e5eSBarry Smith } else if (bz == DMDA_BOUNDARY_MIRROR) { 848c2859e5eSBarry Smith for (j=0; j<x; j++) { idx[nn++] = 0;} 84947c6ae99SBarry Smith } 85047c6ae99SBarry Smith 85147c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 852ce00eea3SSatish Balay x_t = lx[n5 % m]; 85347c6ae99SBarry Smith y_t = y; 85447c6ae99SBarry Smith z_t = lz[n5 / (m*n)]; 85547c6ae99SBarry Smith s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8566f951b95Secoon if (twod && (s_t < 0)) {s_t = bases[n5] + i*x_t + x_t*y_t*z_t - x_t*y_t;} /* 2D case */ 85747c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 85847c6ae99SBarry Smith } 85947c6ae99SBarry Smith } 86047c6ae99SBarry Smith 86147c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 86247c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 863ce00eea3SSatish Balay x_t = lx[n6 % m]; 86447c6ae99SBarry Smith y_t = ly[(n6 % (m*n))/m]; 86547c6ae99SBarry Smith z_t = lz[n6 / (m*n)]; 86647c6ae99SBarry Smith s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8676f951b95Secoon 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 */ 86847c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 86947c6ae99SBarry Smith } 87047c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 87147c6ae99SBarry Smith x_t = x; 87247c6ae99SBarry Smith y_t = ly[(n7 % (m*n))/m]; 87347c6ae99SBarry Smith z_t = lz[n7 / (m*n)]; 87447c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8756f951b95Secoon 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 */ 87647c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 87747c6ae99SBarry Smith } 87847c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 879ce00eea3SSatish Balay x_t = lx[n8 % m]; 88047c6ae99SBarry Smith y_t = ly[(n8 % (m*n))/m]; 88147c6ae99SBarry Smith z_t = lz[n8 / (m*n)]; 88247c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8836f951b95Secoon 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 */ 88447c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 88547c6ae99SBarry Smith } 88647c6ae99SBarry Smith } 88747c6ae99SBarry Smith } 88847c6ae99SBarry Smith 88947c6ae99SBarry Smith /* Middle Level */ 89047c6ae99SBarry Smith for (k=0; k<z; k++) { 89147c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 89247c6ae99SBarry Smith if (n9 >= 0) { /* left below */ 893ce00eea3SSatish Balay x_t = lx[n9 % m]; 89447c6ae99SBarry Smith y_t = ly[(n9 % (m*n))/m]; 89547c6ae99SBarry Smith /* z_t = z; */ 89647c6ae99SBarry Smith s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 89747c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 89847c6ae99SBarry Smith } 89947c6ae99SBarry Smith if (n10 >= 0) { /* directly below */ 90047c6ae99SBarry Smith x_t = x; 90147c6ae99SBarry Smith y_t = ly[(n10 % (m*n))/m]; 90247c6ae99SBarry Smith /* z_t = z; */ 90347c6ae99SBarry Smith s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 90447c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 905c2859e5eSBarry Smith } else if (by == DMDA_BOUNDARY_MIRROR) { 906c2859e5eSBarry Smith for (j=0; j<x; j++) { idx[nn++] = 0;} 90747c6ae99SBarry Smith } 90847c6ae99SBarry Smith if (n11 >= 0) { /* right below */ 909ce00eea3SSatish Balay x_t = lx[n11 % m]; 91047c6ae99SBarry Smith y_t = ly[(n11 % (m*n))/m]; 91147c6ae99SBarry Smith /* z_t = z; */ 91247c6ae99SBarry Smith s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 91347c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 91447c6ae99SBarry Smith } 91547c6ae99SBarry Smith } 91647c6ae99SBarry Smith 91747c6ae99SBarry Smith for (i=0; i<y; i++) { 91847c6ae99SBarry Smith if (n12 >= 0) { /* directly left */ 919ce00eea3SSatish Balay x_t = lx[n12 % m]; 92047c6ae99SBarry Smith y_t = y; 92147c6ae99SBarry Smith /* z_t = z; */ 92247c6ae99SBarry Smith s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t; 92347c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 924c2859e5eSBarry Smith } else if (bx == DMDA_BOUNDARY_MIRROR) { 925c2859e5eSBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = 0;} 92647c6ae99SBarry Smith } 92747c6ae99SBarry Smith 92847c6ae99SBarry Smith /* Interior */ 92947c6ae99SBarry Smith s_t = bases[rank] + i*x + k*x*y; 93047c6ae99SBarry Smith for (j=0; j<x; j++) { idx[nn++] = s_t++;} 93147c6ae99SBarry Smith 93247c6ae99SBarry Smith if (n14 >= 0) { /* directly right */ 933ce00eea3SSatish Balay x_t = lx[n14 % m]; 93447c6ae99SBarry Smith y_t = y; 93547c6ae99SBarry Smith /* z_t = z; */ 93647c6ae99SBarry Smith s_t = bases[n14] + i*x_t + k*x_t*y_t; 93747c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 938c2859e5eSBarry Smith } else if (bx == DMDA_BOUNDARY_MIRROR) { 939c2859e5eSBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = 0;} 94047c6ae99SBarry Smith } 94147c6ae99SBarry Smith } 94247c6ae99SBarry Smith 94347c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 94447c6ae99SBarry Smith if (n15 >= 0) { /* left above */ 945ce00eea3SSatish Balay x_t = lx[n15 % m]; 94647c6ae99SBarry Smith y_t = ly[(n15 % (m*n))/m]; 94747c6ae99SBarry Smith /* z_t = z; */ 94847c6ae99SBarry Smith s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t; 94947c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 95047c6ae99SBarry Smith } 95147c6ae99SBarry Smith if (n16 >= 0) { /* directly above */ 95247c6ae99SBarry Smith x_t = x; 95347c6ae99SBarry Smith y_t = ly[(n16 % (m*n))/m]; 95447c6ae99SBarry Smith /* z_t = z; */ 95547c6ae99SBarry Smith s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t; 95647c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 957c2859e5eSBarry Smith } else if (by == DMDA_BOUNDARY_MIRROR) { 958c2859e5eSBarry Smith for (j=0; j<x; j++) { idx[nn++] = 0;} 95947c6ae99SBarry Smith } 96047c6ae99SBarry Smith if (n17 >= 0) { /* right above */ 961ce00eea3SSatish Balay x_t = lx[n17 % m]; 96247c6ae99SBarry Smith y_t = ly[(n17 % (m*n))/m]; 96347c6ae99SBarry Smith /* z_t = z; */ 96447c6ae99SBarry Smith s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t; 96547c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 96647c6ae99SBarry Smith } 96747c6ae99SBarry Smith } 96847c6ae99SBarry Smith } 96947c6ae99SBarry Smith 97047c6ae99SBarry Smith /* Upper Level */ 97147c6ae99SBarry Smith for (k=0; k<s_z; k++) { 97247c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 97347c6ae99SBarry Smith if (n18 >= 0) { /* left below */ 974ce00eea3SSatish Balay x_t = lx[n18 % m]; 97547c6ae99SBarry Smith y_t = ly[(n18 % (m*n))/m]; 97647c6ae99SBarry Smith /* z_t = lz[n18 / (m*n)]; */ 97747c6ae99SBarry Smith s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 9786f951b95Secoon if (twod && (s_t >= M*N*P)) {s_t = bases[n18] - (s_y-i)*x_t -s_x + x_t*y_t;} /* 2d case */ 97947c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 98047c6ae99SBarry Smith } 98147c6ae99SBarry Smith if (n19 >= 0) { /* directly below */ 98247c6ae99SBarry Smith x_t = x; 98347c6ae99SBarry Smith y_t = ly[(n19 % (m*n))/m]; 98447c6ae99SBarry Smith /* z_t = lz[n19 / (m*n)]; */ 98547c6ae99SBarry Smith s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 9866f951b95Secoon if (twod && (s_t >= M*N*P)) {s_t = bases[n19] - (s_y+1-i)*x_t + x_t*y_t;} /* 2d case */ 98747c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 98847c6ae99SBarry Smith } 98947c6ae99SBarry Smith if (n20 >= 0) { /* right below */ 990ce00eea3SSatish Balay x_t = lx[n20 % m]; 99147c6ae99SBarry Smith y_t = ly[(n20 % (m*n))/m]; 99247c6ae99SBarry Smith /* z_t = lz[n20 / (m*n)]; */ 99347c6ae99SBarry Smith s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 9946f951b95Secoon if (twod && (s_t >= M*N*P)) {s_t = bases[n20] - (s_y+1-i)*x_t + x_t*y_t;} /* 2d case */ 99547c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 99647c6ae99SBarry Smith } 99747c6ae99SBarry Smith } 99847c6ae99SBarry Smith 99947c6ae99SBarry Smith for (i=0; i<y; i++) { 100047c6ae99SBarry Smith if (n21 >= 0) { /* directly left */ 1001ce00eea3SSatish Balay x_t = lx[n21 % m]; 100247c6ae99SBarry Smith y_t = y; 100347c6ae99SBarry Smith /* z_t = lz[n21 / (m*n)]; */ 100447c6ae99SBarry Smith s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t; 10056f951b95Secoon if (twod && (s_t >= M*N*P)) {s_t = bases[n21] + (i+1)*x_t - s_x;} /* 2d case */ 100647c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 100747c6ae99SBarry Smith } 100847c6ae99SBarry Smith 100947c6ae99SBarry Smith if (n22 >= 0) { /* middle */ 101047c6ae99SBarry Smith x_t = x; 101147c6ae99SBarry Smith y_t = y; 101247c6ae99SBarry Smith /* z_t = lz[n22 / (m*n)]; */ 101347c6ae99SBarry Smith s_t = bases[n22] + i*x_t + k*x_t*y_t; 10146f951b95Secoon if (twod && (s_t >= M*N*P)) {s_t = bases[n22] + i*x_t;} /* 2d case */ 101547c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 1016c2859e5eSBarry Smith } else if (bz == DMDA_BOUNDARY_MIRROR) { 1017c2859e5eSBarry Smith for (j=0; j<x; j++) { idx[nn++] = 0;} 101847c6ae99SBarry Smith } 101947c6ae99SBarry Smith 102047c6ae99SBarry Smith if (n23 >= 0) { /* directly right */ 1021ce00eea3SSatish Balay x_t = lx[n23 % m]; 102247c6ae99SBarry Smith y_t = y; 102347c6ae99SBarry Smith /* z_t = lz[n23 / (m*n)]; */ 102447c6ae99SBarry Smith s_t = bases[n23] + i*x_t + k*x_t*y_t; 10256f951b95Secoon if (twod && (s_t >= M*N*P)) {s_t = bases[n23] + i*x_t;} /* 2d case */ 102647c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 102747c6ae99SBarry Smith } 102847c6ae99SBarry Smith } 102947c6ae99SBarry Smith 103047c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 103147c6ae99SBarry Smith if (n24 >= 0) { /* left above */ 1032ce00eea3SSatish Balay x_t = lx[n24 % m]; 103347c6ae99SBarry Smith y_t = ly[(n24 % (m*n))/m]; 103447c6ae99SBarry Smith /* z_t = lz[n24 / (m*n)]; */ 103547c6ae99SBarry Smith s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t; 10366f951b95Secoon if (twod && (s_t >= M*N*P)) {s_t = bases[n24] + i*x_t - s_x;} /* 2d case */ 103747c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 103847c6ae99SBarry Smith } 103947c6ae99SBarry Smith if (n25 >= 0) { /* directly above */ 104047c6ae99SBarry Smith x_t = x; 104147c6ae99SBarry Smith y_t = ly[(n25 % (m*n))/m]; 104247c6ae99SBarry Smith /* z_t = lz[n25 / (m*n)]; */ 104347c6ae99SBarry Smith s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t; 10446f951b95Secoon if (twod && (s_t >= M*N*P)) {s_t = bases[n25] + (i-1)*x_t;} /* 2d case */ 104547c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 104647c6ae99SBarry Smith } 104747c6ae99SBarry Smith if (n26 >= 0) { /* right above */ 1048ce00eea3SSatish Balay x_t = lx[n26 % m]; 104947c6ae99SBarry Smith y_t = ly[(n26 % (m*n))/m]; 105047c6ae99SBarry Smith /* z_t = lz[n26 / (m*n)]; */ 105147c6ae99SBarry Smith s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t; 10526f951b95Secoon if (twod && (s_t >= M*N*P)) {s_t = bases[n26] + (i-1)*x_t;} /* 2d case */ 105347c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 105447c6ae99SBarry Smith } 105547c6ae99SBarry Smith } 105647c6ae99SBarry Smith } 1057ce00eea3SSatish Balay 1058ce00eea3SSatish Balay ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_COPY_VALUES,&from);CHKERRQ(ierr); 105947c6ae99SBarry Smith ierr = VecScatterCreate(global,from,local,to,>ol);CHKERRQ(ierr); 106047c6ae99SBarry Smith ierr = PetscLogObjectParent(da,gtol);CHKERRQ(ierr); 1061fcfd50ebSBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 1062fcfd50ebSBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 106347c6ae99SBarry Smith 1064*d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 106547c6ae99SBarry Smith n0 = sn0; n1 = sn1; n2 = sn2; n3 = sn3; n5 = sn5; n6 = sn6; n7 = sn7; 106647c6ae99SBarry Smith n8 = sn8; n9 = sn9; n11 = sn11; n15 = sn15; n17 = sn17; n18 = sn18; 106747c6ae99SBarry Smith n19 = sn19; n20 = sn20; n21 = sn21; n23 = sn23; n24 = sn24; n25 = sn25; 106847c6ae99SBarry Smith n26 = sn26; 1069ce00eea3SSatish Balay } 107047c6ae99SBarry Smith 107188661749SPeter Brune if (((stencil_type == DMDA_STENCIL_STAR) || 10721321219cSEthan Coon (bx != DMDA_BOUNDARY_PERIODIC && bx) || 10731321219cSEthan Coon (by != DMDA_BOUNDARY_PERIODIC && by) || 1074*d9c9ebe5SPeter Brune (bz != DMDA_BOUNDARY_PERIODIC && bz))) { 1075ce00eea3SSatish Balay /* 1076ce00eea3SSatish Balay Recompute the local to global mappings, this time keeping the 1077ce00eea3SSatish Balay information about the cross corner processor numbers. 1078ce00eea3SSatish Balay */ 107947c6ae99SBarry Smith nn = 0; 108047c6ae99SBarry Smith /* Bottom Level */ 108147c6ae99SBarry Smith for (k=0; k<s_z; k++) { 108247c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 108347c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 1084ce00eea3SSatish Balay x_t = lx[n0 % m]; 108547c6ae99SBarry Smith y_t = ly[(n0 % (m*n))/m]; 108647c6ae99SBarry Smith z_t = lz[n0 / (m*n)]; 108747c6ae99SBarry 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; 108847c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1089ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0 && Zs-zs < 0) { 1090ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 109147c6ae99SBarry Smith } 109247c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 109347c6ae99SBarry Smith x_t = x; 109447c6ae99SBarry Smith y_t = ly[(n1 % (m*n))/m]; 109547c6ae99SBarry Smith z_t = lz[n1 / (m*n)]; 109647c6ae99SBarry 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; 109747c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 1098ce00eea3SSatish Balay } else if (Ys-ys < 0 && Zs-zs < 0) { 1099ce00eea3SSatish Balay for (j=0; j<x; j++) { idx[nn++] = -1;} 110047c6ae99SBarry Smith } 110147c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 1102ce00eea3SSatish Balay x_t = lx[n2 % m]; 110347c6ae99SBarry Smith y_t = ly[(n2 % (m*n))/m]; 110447c6ae99SBarry Smith z_t = lz[n2 / (m*n)]; 110547c6ae99SBarry 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; 110647c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1107ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0 && Zs-zs < 0) { 1108ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 110947c6ae99SBarry Smith } 111047c6ae99SBarry Smith } 111147c6ae99SBarry Smith 111247c6ae99SBarry Smith for (i=0; i<y; i++) { 111347c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 1114ce00eea3SSatish Balay x_t = lx[n3 % m]; 111547c6ae99SBarry Smith y_t = y; 111647c6ae99SBarry Smith z_t = lz[n3 / (m*n)]; 111747c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 111847c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1119ce00eea3SSatish Balay } else if (Xs-xs < 0 && Zs-zs < 0) { 1120ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 112147c6ae99SBarry Smith } 112247c6ae99SBarry Smith 112347c6ae99SBarry Smith if (n4 >= 0) { /* middle */ 112447c6ae99SBarry Smith x_t = x; 112547c6ae99SBarry Smith y_t = y; 112647c6ae99SBarry Smith z_t = lz[n4 / (m*n)]; 112747c6ae99SBarry Smith s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 112847c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 1129ce00eea3SSatish Balay } else if (Zs-zs < 0) { 1130c2859e5eSBarry Smith if (bz == DMDA_BOUNDARY_MIRROR) { 1131c2859e5eSBarry Smith for (j=0; j<x; j++) { idx[nn++] = 0;} 1132c2859e5eSBarry Smith } else { 1133ce00eea3SSatish Balay for (j=0; j<x; j++) { idx[nn++] = -1;} 113447c6ae99SBarry Smith } 1135c2859e5eSBarry Smith } 113647c6ae99SBarry Smith 113747c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 1138ce00eea3SSatish Balay x_t = lx[n5 % m]; 113947c6ae99SBarry Smith y_t = y; 114047c6ae99SBarry Smith z_t = lz[n5 / (m*n)]; 114147c6ae99SBarry Smith s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 114247c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1143ce00eea3SSatish Balay } else if (xe-Xe < 0 && Zs-zs < 0) { 1144ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 114547c6ae99SBarry Smith } 114647c6ae99SBarry Smith } 114747c6ae99SBarry Smith 114847c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 114947c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 1150ce00eea3SSatish Balay x_t = lx[n6 % m]; 115147c6ae99SBarry Smith y_t = ly[(n6 % (m*n))/m]; 115247c6ae99SBarry Smith z_t = lz[n6 / (m*n)]; 115347c6ae99SBarry Smith s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 115447c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1155ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0 && Zs-zs < 0) { 1156ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 115747c6ae99SBarry Smith } 115847c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 115947c6ae99SBarry Smith x_t = x; 116047c6ae99SBarry Smith y_t = ly[(n7 % (m*n))/m]; 116147c6ae99SBarry Smith z_t = lz[n7 / (m*n)]; 116247c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 116347c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 1164ce00eea3SSatish Balay } else if (ye-Ye < 0 && Zs-zs < 0) { 1165ce00eea3SSatish Balay for (j=0; j<x; j++) { idx[nn++] = -1;} 116647c6ae99SBarry Smith } 116747c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 1168ce00eea3SSatish Balay x_t = lx[n8 % m]; 116947c6ae99SBarry Smith y_t = ly[(n8 % (m*n))/m]; 117047c6ae99SBarry Smith z_t = lz[n8 / (m*n)]; 117147c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 117247c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1173ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0 && Zs-zs < 0) { 1174ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 117547c6ae99SBarry Smith } 117647c6ae99SBarry Smith } 117747c6ae99SBarry Smith } 117847c6ae99SBarry Smith 117947c6ae99SBarry Smith /* Middle Level */ 118047c6ae99SBarry Smith for (k=0; k<z; k++) { 118147c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 118247c6ae99SBarry Smith if (n9 >= 0) { /* left below */ 1183ce00eea3SSatish Balay x_t = lx[n9 % m]; 118447c6ae99SBarry Smith y_t = ly[(n9 % (m*n))/m]; 118547c6ae99SBarry Smith /* z_t = z; */ 118647c6ae99SBarry Smith s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 118747c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1188ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0) { 1189ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 119047c6ae99SBarry Smith } 119147c6ae99SBarry Smith if (n10 >= 0) { /* directly below */ 119247c6ae99SBarry Smith x_t = x; 119347c6ae99SBarry Smith y_t = ly[(n10 % (m*n))/m]; 119447c6ae99SBarry Smith /* z_t = z; */ 119547c6ae99SBarry Smith s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 119647c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 1197ce00eea3SSatish Balay } else if (Ys-ys < 0) { 1198c2859e5eSBarry Smith if (by == DMDA_BOUNDARY_MIRROR) { 1199ce00eea3SSatish Balay for (j=0; j<x; j++) { idx[nn++] = -1;} 1200c2859e5eSBarry Smith } else { 1201c2859e5eSBarry Smith for (j=0; j<x; j++) { idx[nn++] = -1;} 1202c2859e5eSBarry Smith } 120347c6ae99SBarry Smith } 120447c6ae99SBarry Smith if (n11 >= 0) { /* right below */ 1205ce00eea3SSatish Balay x_t = lx[n11 % m]; 120647c6ae99SBarry Smith y_t = ly[(n11 % (m*n))/m]; 120747c6ae99SBarry Smith /* z_t = z; */ 120847c6ae99SBarry Smith s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 120947c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1210ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0) { 1211ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 121247c6ae99SBarry Smith } 121347c6ae99SBarry Smith } 121447c6ae99SBarry Smith 121547c6ae99SBarry Smith for (i=0; i<y; i++) { 121647c6ae99SBarry Smith if (n12 >= 0) { /* directly left */ 1217ce00eea3SSatish Balay x_t = lx[n12 % m]; 121847c6ae99SBarry Smith y_t = y; 121947c6ae99SBarry Smith /* z_t = z; */ 122047c6ae99SBarry Smith s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t; 122147c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1222ce00eea3SSatish Balay } else if (Xs-xs < 0) { 1223c2859e5eSBarry Smith if (bx == DMDA_BOUNDARY_MIRROR) { 1224c2859e5eSBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = 0;} 1225c2859e5eSBarry Smith } else { 1226ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 122747c6ae99SBarry Smith } 1228c2859e5eSBarry Smith } 122947c6ae99SBarry Smith 123047c6ae99SBarry Smith /* Interior */ 123147c6ae99SBarry Smith s_t = bases[rank] + i*x + k*x*y; 123247c6ae99SBarry Smith for (j=0; j<x; j++) { idx[nn++] = s_t++;} 123347c6ae99SBarry Smith 123447c6ae99SBarry Smith if (n14 >= 0) { /* directly right */ 1235ce00eea3SSatish Balay x_t = lx[n14 % m]; 123647c6ae99SBarry Smith y_t = y; 123747c6ae99SBarry Smith /* z_t = z; */ 123847c6ae99SBarry Smith s_t = bases[n14] + i*x_t + k*x_t*y_t; 123947c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1240ce00eea3SSatish Balay } else if (xe-Xe < 0) { 1241c2859e5eSBarry Smith if (bx == DMDA_BOUNDARY_MIRROR) { 1242c2859e5eSBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = 0;} 1243c2859e5eSBarry Smith } else { 1244ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 124547c6ae99SBarry Smith } 124647c6ae99SBarry Smith } 1247c2859e5eSBarry Smith } 124847c6ae99SBarry Smith 124947c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 125047c6ae99SBarry Smith if (n15 >= 0) { /* left above */ 1251ce00eea3SSatish Balay x_t = lx[n15 % m]; 125247c6ae99SBarry Smith y_t = ly[(n15 % (m*n))/m]; 125347c6ae99SBarry Smith /* z_t = z; */ 125447c6ae99SBarry Smith s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t; 125547c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1256ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0) { 1257ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 125847c6ae99SBarry Smith } 125947c6ae99SBarry Smith if (n16 >= 0) { /* directly above */ 126047c6ae99SBarry Smith x_t = x; 126147c6ae99SBarry Smith y_t = ly[(n16 % (m*n))/m]; 126247c6ae99SBarry Smith /* z_t = z; */ 126347c6ae99SBarry Smith s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t; 126447c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 1265ce00eea3SSatish Balay } else if (ye-Ye < 0) { 1266c2859e5eSBarry Smith if (by == DMDA_BOUNDARY_MIRROR) { 1267c2859e5eSBarry Smith for (j=0; j<x; j++) { idx[nn++] = 0;} 1268c2859e5eSBarry Smith } else { 1269ce00eea3SSatish Balay for (j=0; j<x; j++) { idx[nn++] = -1;} 127047c6ae99SBarry Smith } 1271c2859e5eSBarry Smith } 127247c6ae99SBarry Smith if (n17 >= 0) { /* right above */ 1273ce00eea3SSatish Balay x_t = lx[n17 % m]; 127447c6ae99SBarry Smith y_t = ly[(n17 % (m*n))/m]; 127547c6ae99SBarry Smith /* z_t = z; */ 127647c6ae99SBarry Smith s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t; 127747c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1278ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0) { 1279ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 128047c6ae99SBarry Smith } 128147c6ae99SBarry Smith } 128247c6ae99SBarry Smith } 128347c6ae99SBarry Smith 128447c6ae99SBarry Smith /* Upper Level */ 128547c6ae99SBarry Smith for (k=0; k<s_z; k++) { 128647c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 128747c6ae99SBarry Smith if (n18 >= 0) { /* left below */ 1288ce00eea3SSatish Balay x_t = lx[n18 % m]; 128947c6ae99SBarry Smith y_t = ly[(n18 % (m*n))/m]; 129047c6ae99SBarry Smith /* z_t = lz[n18 / (m*n)]; */ 129147c6ae99SBarry Smith s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 129247c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1293ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0 && ze-Ze < 0) { 1294ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 129547c6ae99SBarry Smith } 129647c6ae99SBarry Smith if (n19 >= 0) { /* directly below */ 129747c6ae99SBarry Smith x_t = x; 129847c6ae99SBarry Smith y_t = ly[(n19 % (m*n))/m]; 129947c6ae99SBarry Smith /* z_t = lz[n19 / (m*n)]; */ 130047c6ae99SBarry Smith s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 130147c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 1302ce00eea3SSatish Balay } else if (Ys-ys < 0 && ze-Ze < 0) { 1303ce00eea3SSatish Balay for (j=0; j<x; j++) { idx[nn++] = -1;} 130447c6ae99SBarry Smith } 130547c6ae99SBarry Smith if (n20 >= 0) { /* right below */ 1306ce00eea3SSatish Balay x_t = lx[n20 % m]; 130747c6ae99SBarry Smith y_t = ly[(n20 % (m*n))/m]; 130847c6ae99SBarry Smith /* z_t = lz[n20 / (m*n)]; */ 130947c6ae99SBarry Smith s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 131047c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1311ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0 && ze-Ze < 0) { 1312ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 131347c6ae99SBarry Smith } 131447c6ae99SBarry Smith } 131547c6ae99SBarry Smith 131647c6ae99SBarry Smith for (i=0; i<y; i++) { 131747c6ae99SBarry Smith if (n21 >= 0) { /* directly left */ 1318ce00eea3SSatish Balay x_t = lx[n21 % m]; 131947c6ae99SBarry Smith y_t = y; 132047c6ae99SBarry Smith /* z_t = lz[n21 / (m*n)]; */ 132147c6ae99SBarry Smith s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t; 132247c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1323ce00eea3SSatish Balay } else if (Xs-xs < 0 && ze-Ze < 0) { 1324ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 132547c6ae99SBarry Smith } 132647c6ae99SBarry Smith 132747c6ae99SBarry Smith if (n22 >= 0) { /* middle */ 132847c6ae99SBarry Smith x_t = x; 132947c6ae99SBarry Smith y_t = y; 133047c6ae99SBarry Smith /* z_t = lz[n22 / (m*n)]; */ 133147c6ae99SBarry Smith s_t = bases[n22] + i*x_t + k*x_t*y_t; 133247c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 1333ce00eea3SSatish Balay } else if (ze-Ze < 0) { 1334c2859e5eSBarry Smith if (bz == DMDA_BOUNDARY_MIRROR) { 1335c2859e5eSBarry Smith for (j=0; j<x; j++) { idx[nn++] = 0;} 1336c2859e5eSBarry Smith } else { 1337ce00eea3SSatish Balay for (j=0; j<x; j++) { idx[nn++] = -1;} 133847c6ae99SBarry Smith } 1339c2859e5eSBarry Smith } 134047c6ae99SBarry Smith 134147c6ae99SBarry Smith if (n23 >= 0) { /* directly right */ 1342ce00eea3SSatish Balay x_t = lx[n23 % m]; 134347c6ae99SBarry Smith y_t = y; 134447c6ae99SBarry Smith /* z_t = lz[n23 / (m*n)]; */ 134547c6ae99SBarry Smith s_t = bases[n23] + i*x_t + k*x_t*y_t; 134647c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1347ce00eea3SSatish Balay } else if (xe-Xe < 0 && ze-Ze < 0) { 1348ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 134947c6ae99SBarry Smith } 135047c6ae99SBarry Smith } 135147c6ae99SBarry Smith 135247c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 135347c6ae99SBarry Smith if (n24 >= 0) { /* left above */ 1354ce00eea3SSatish Balay x_t = lx[n24 % m]; 135547c6ae99SBarry Smith y_t = ly[(n24 % (m*n))/m]; 135647c6ae99SBarry Smith /* z_t = lz[n24 / (m*n)]; */ 135747c6ae99SBarry Smith s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t; 135847c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1359ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0 && ze-Ze < 0) { 1360ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 136147c6ae99SBarry Smith } 136247c6ae99SBarry Smith if (n25 >= 0) { /* directly above */ 136347c6ae99SBarry Smith x_t = x; 136447c6ae99SBarry Smith y_t = ly[(n25 % (m*n))/m]; 136547c6ae99SBarry Smith /* z_t = lz[n25 / (m*n)]; */ 136647c6ae99SBarry Smith s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t; 136747c6ae99SBarry Smith for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 1368ce00eea3SSatish Balay } else if (ye-Ye < 0 && ze-Ze < 0) { 1369ce00eea3SSatish Balay for (j=0; j<x; j++) { idx[nn++] = -1;} 137047c6ae99SBarry Smith } 137147c6ae99SBarry Smith if (n26 >= 0) { /* right above */ 1372ce00eea3SSatish Balay x_t = lx[n26 % m]; 137347c6ae99SBarry Smith y_t = ly[(n26 % (m*n))/m]; 137447c6ae99SBarry Smith /* z_t = lz[n26 / (m*n)]; */ 137547c6ae99SBarry Smith s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t; 137647c6ae99SBarry Smith for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1377ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0 && ze-Ze < 0) { 1378ce00eea3SSatish Balay for (j=0; j<s_x; j++) { idx[nn++] = -1;} 137947c6ae99SBarry Smith } 138047c6ae99SBarry Smith } 138147c6ae99SBarry Smith } 138247c6ae99SBarry Smith } 138347c6ae99SBarry Smith /* 138447c6ae99SBarry Smith Set the local to global ordering in the global vector, this allows use 138547c6ae99SBarry Smith of VecSetValuesLocal(). 138647c6ae99SBarry Smith */ 1387ce00eea3SSatish Balay ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_OWN_POINTER,<ogis);CHKERRQ(ierr); 1388ce00eea3SSatish Balay ierr = PetscMalloc(nn*dof*sizeof(PetscInt),&idx_cpy);CHKERRQ(ierr); 1389db87c5ecSEthan Coon ierr = PetscLogObjectMemory(da,nn*dof*sizeof(PetscInt));CHKERRQ(ierr); 13903bf036e2SBarry Smith ierr = ISGetIndices(ltogis, &idx_full);CHKERRQ(ierr); 1391ce00eea3SSatish Balay ierr = PetscMemcpy(idx_cpy,idx_full,nn*dof*sizeof(PetscInt));CHKERRQ(ierr); 13923bf036e2SBarry Smith ierr = ISRestoreIndices(ltogis, &idx_full);CHKERRQ(ierr); 1393ce00eea3SSatish Balay ierr = ISLocalToGlobalMappingCreateIS(ltogis,&da->ltogmap);CHKERRQ(ierr); 1394ce00eea3SSatish Balay ierr = PetscLogObjectParent(da,da->ltogmap);CHKERRQ(ierr); 1395fcfd50ebSBarry Smith ierr = ISDestroy(<ogis);CHKERRQ(ierr); 13961411c6eeSJed Brown ierr = ISLocalToGlobalMappingBlock(da->ltogmap,dd->w,&da->ltogmapb);CHKERRQ(ierr); 13971411c6eeSJed Brown ierr = PetscLogObjectParent(da,da->ltogmap);CHKERRQ(ierr); 139847c6ae99SBarry Smith 1399ce00eea3SSatish Balay ierr = PetscFree2(bases,ldims);CHKERRQ(ierr); 1400ce00eea3SSatish Balay dd->m = m; dd->n = n; dd->p = p; 1401ce00eea3SSatish Balay /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */ 1402ce00eea3SSatish Balay dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = zs; dd->ze = ze; 1403ce00eea3SSatish Balay dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = Zs; dd->Ze = Ze; 1404ce00eea3SSatish Balay 1405fcfd50ebSBarry Smith ierr = VecDestroy(&local);CHKERRQ(ierr); 1406fcfd50ebSBarry Smith ierr = VecDestroy(&global);CHKERRQ(ierr); 1407ce00eea3SSatish Balay 1408ce00eea3SSatish Balay dd->gtol = gtol; 1409ce00eea3SSatish Balay dd->ltog = ltog; 1410ce00eea3SSatish Balay dd->idx = idx_cpy; 1411ce00eea3SSatish Balay dd->Nl = nn*dof; 1412ce00eea3SSatish Balay dd->base = base; 1413ce00eea3SSatish Balay da->ops->view = DMView_DA_3d; 141447c6ae99SBarry Smith dd->ltol = PETSC_NULL; 141547c6ae99SBarry Smith dd->ao = PETSC_NULL; 1416ce00eea3SSatish Balay 141747c6ae99SBarry Smith PetscFunctionReturn(0); 141847c6ae99SBarry Smith } 1419564755cdSBarry Smith 142047c6ae99SBarry Smith 142147c6ae99SBarry Smith #undef __FUNCT__ 1422aa219208SBarry Smith #define __FUNCT__ "DMDACreate3d" 142347c6ae99SBarry Smith /*@C 1424aa219208SBarry Smith DMDACreate3d - Creates an object that will manage the communication of three-dimensional 142547c6ae99SBarry Smith regular array data that is distributed across some processors. 142647c6ae99SBarry Smith 142747c6ae99SBarry Smith Collective on MPI_Comm 142847c6ae99SBarry Smith 142947c6ae99SBarry Smith Input Parameters: 143047c6ae99SBarry Smith + comm - MPI communicator 14311321219cSEthan Coon . bx,by,bz - type of ghost nodes the array have. 14321321219cSEthan Coon Use one of DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_PERIODIC. 1433aa219208SBarry Smith . stencil_type - Type of stencil (DMDA_STENCIL_STAR or DMDA_STENCIL_BOX) 143447c6ae99SBarry 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 143547c6ae99SBarry Smith from the command line with -da_grid_x <M> -da_grid_y <N> -da_grid_z <P>) 143647c6ae99SBarry Smith . m,n,p - corresponding number of processors in each dimension 143747c6ae99SBarry Smith (or PETSC_DECIDE to have calculated) 143847c6ae99SBarry Smith . dof - number of degrees of freedom per node 143910d7c030SMatthew G Knepley . s - stencil width 144010d7c030SMatthew G Knepley - lx, ly, lz - arrays containing the number of nodes in each cell along 144147c6ae99SBarry Smith the x, y, and z coordinates, or PETSC_NULL. If non-null, these 144247c6ae99SBarry Smith must be of length as m,n,p and the corresponding 144347c6ae99SBarry Smith m,n, or p cannot be PETSC_DECIDE. Sum of the lx[] entries must be M, sum of 144447c6ae99SBarry Smith the ly[] must N, sum of the lz[] must be P 144547c6ae99SBarry Smith 144647c6ae99SBarry Smith Output Parameter: 144747c6ae99SBarry Smith . da - the resulting distributed array object 144847c6ae99SBarry Smith 144947c6ae99SBarry Smith Options Database Key: 1450706a11cbSBarry Smith + -dm_view - Calls DMView() at the conclusion of DMDACreate3d() 145147c6ae99SBarry Smith . -da_grid_x <nx> - number of grid points in x direction, if M < 0 145247c6ae99SBarry Smith . -da_grid_y <ny> - number of grid points in y direction, if N < 0 145347c6ae99SBarry Smith . -da_grid_z <nz> - number of grid points in z direction, if P < 0 145447c6ae99SBarry Smith . -da_processors_x <MX> - number of processors in x direction 145547c6ae99SBarry Smith . -da_processors_y <MY> - number of processors in y direction 145647c6ae99SBarry Smith . -da_processors_z <MZ> - number of processors in z direction 1457e0f5d30fSBarry Smith . -da_refine_x <rx> - refinement ratio in x direction 1458e0f5d30fSBarry Smith . -da_refine_y <ry> - refinement ratio in y direction 1459e0f5d30fSBarry Smith . -da_refine_z <rz>- refinement ratio in z directio 1460e0f5d30fSBarry Smith - -da_refine <n> - refine the DMDA n times before creating it, , if M, N, or P < 0 146147c6ae99SBarry Smith 146247c6ae99SBarry Smith Level: beginner 146347c6ae99SBarry Smith 146447c6ae99SBarry Smith Notes: 1465aa219208SBarry Smith The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the 1466aa219208SBarry Smith standard 7-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes 146747c6ae99SBarry Smith the standard 27-pt stencil. 146847c6ae99SBarry Smith 1469aa219208SBarry Smith The array data itself is NOT stored in the DMDA, it is stored in Vec objects; 1470564755cdSBarry Smith The appropriate vector objects can be obtained with calls to DMCreateGlobalVector() 1471564755cdSBarry Smith and DMCreateLocalVector() and calls to VecDuplicate() if more are needed. 147247c6ae99SBarry Smith 147347c6ae99SBarry Smith .keywords: distributed array, create, three-dimensional 147447c6ae99SBarry Smith 1475aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate2d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(), 1476aa219208SBarry Smith DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMDALocalToLocalBegin(), DMDALocalToLocalEnd(), DMDASetRefinementFactor(), 1477d461ba97SBarry Smith DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges() 147847c6ae99SBarry Smith 147947c6ae99SBarry Smith @*/ 14801321219cSEthan Coon PetscErrorCode DMDACreate3d(MPI_Comm comm,DMDABoundaryType bx,DMDABoundaryType by,DMDABoundaryType bz,DMDAStencilType stencil_type,PetscInt M, 14819a42bb27SBarry 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) 148247c6ae99SBarry Smith { 148347c6ae99SBarry Smith PetscErrorCode ierr; 148447c6ae99SBarry Smith 148547c6ae99SBarry Smith PetscFunctionBegin; 1486aa219208SBarry Smith ierr = DMDACreate(comm, da);CHKERRQ(ierr); 1487aa219208SBarry Smith ierr = DMDASetDim(*da, 3);CHKERRQ(ierr); 1488aa219208SBarry Smith ierr = DMDASetSizes(*da, M, N, P);CHKERRQ(ierr); 1489aa219208SBarry Smith ierr = DMDASetNumProcs(*da, m, n, p);CHKERRQ(ierr); 14901321219cSEthan Coon ierr = DMDASetBoundaryType(*da, bx, by, bz);CHKERRQ(ierr); 1491aa219208SBarry Smith ierr = DMDASetDof(*da, dof);CHKERRQ(ierr); 1492aa219208SBarry Smith ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr); 1493aa219208SBarry Smith ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr); 1494aa219208SBarry Smith ierr = DMDASetOwnershipRanges(*da, lx, ly, lz);CHKERRQ(ierr); 149547c6ae99SBarry Smith /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */ 14969a42bb27SBarry Smith ierr = DMSetFromOptions(*da);CHKERRQ(ierr); 14979a42bb27SBarry Smith ierr = DMSetUp(*da);CHKERRQ(ierr); 1498ca266f36SBarry Smith ierr = DMViewFromOptions(*da,"-dm_view");CHKERRQ(ierr); 149947c6ae99SBarry Smith PetscFunctionReturn(0); 150047c6ae99SBarry Smith } 1501