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 7af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 847c6ae99SBarry Smith 99804daf3SBarry Smith #include <petscdraw.h> 10e0877f53SBarry Smith static PetscErrorCode DMView_DA_3d(DM da,PetscViewer viewer) 1147c6ae99SBarry Smith { 1247c6ae99SBarry Smith PetscErrorCode ierr; 1347c6ae99SBarry Smith PetscMPIInt rank; 14c9493c35SStefano Zampini PetscBool iascii,isdraw,isglvis,isbinary; 1547c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 169a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 179a42bb27SBarry Smith PetscBool ismatlab; 189a42bb27SBarry Smith #endif 1947c6ae99SBarry Smith 2047c6ae99SBarry Smith PetscFunctionBegin; 21ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank);CHKERRQ(ierr); 2247c6ae99SBarry Smith 23251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 24251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 25c9493c35SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis);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 331575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 3447c6ae99SBarry Smith ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 35*76a8abe0SBarry Smith if (format == PETSC_VIEWER_LOAD_BALANCE) { 36*76a8abe0SBarry Smith PetscInt i,nmax = 0,nmin = PETSC_MAX_INT,navg = 0,*nz,nzlocal; 37*76a8abe0SBarry Smith DMDALocalInfo info; 38*76a8abe0SBarry Smith PetscMPIInt size; 39*76a8abe0SBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)da),&size);CHKERRQ(ierr); 40*76a8abe0SBarry Smith ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr); 41*76a8abe0SBarry Smith nzlocal = info.xm*info.ym*info.zm; 42*76a8abe0SBarry Smith ierr = PetscMalloc1(size,&nz);CHKERRQ(ierr); 43*76a8abe0SBarry Smith ierr = MPI_Allgather(&nzlocal,1,MPIU_INT,nz,1,MPIU_INT,PetscObjectComm((PetscObject)da));CHKERRQ(ierr); 44*76a8abe0SBarry Smith for (i=0; i<(PetscInt)size; i++) { 45*76a8abe0SBarry Smith nmax = PetscMax(nmax,nz[i]); 46*76a8abe0SBarry Smith nmin = PetscMin(nmin,nz[i]); 47*76a8abe0SBarry Smith navg += nz[i]; 48*76a8abe0SBarry Smith } 49*76a8abe0SBarry Smith ierr = PetscFree(nz);CHKERRQ(ierr); 50*76a8abe0SBarry Smith navg = navg/size; 51*76a8abe0SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Load Balance - Grid Points: Min %D avg %D max %D\n",nmin,navg,nmax);CHKERRQ(ierr); 52*76a8abe0SBarry Smith PetscFunctionReturn(0); 53*76a8abe0SBarry Smith } 548135c375SStefano Zampini if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL && format != PETSC_VIEWER_ASCII_GLVIS) { 55aa219208SBarry Smith DMDALocalInfo info; 56aa219208SBarry Smith ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr); 5747c6ae99SBarry 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); 5847c6ae99SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"X range of indices: %D %D, Y range of indices: %D %D, Z range of indices: %D %D\n", 5947c6ae99SBarry Smith info.xs,info.xs+info.xm,info.ys,info.ys+info.ym,info.zs,info.zs+info.zm);CHKERRQ(ierr); 6047c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX) 616636e97aSMatthew G Knepley if (da->coordinates) { 6247c6ae99SBarry Smith PetscInt last; 6347c6ae99SBarry Smith const PetscReal *coors; 646636e97aSMatthew G Knepley ierr = VecGetArrayRead(da->coordinates,&coors);CHKERRQ(ierr); 656636e97aSMatthew G Knepley ierr = VecGetLocalSize(da->coordinates,&last);CHKERRQ(ierr); 6647c6ae99SBarry Smith last = last - 3; 6757622a8eSBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Lower left corner %g %g %g : Upper right %g %g %g\n",(double)coors[0],(double)coors[1],(double)coors[2],(double)coors[last],(double)coors[last+1],(double)coors[last+2]);CHKERRQ(ierr); 686636e97aSMatthew G Knepley ierr = VecRestoreArrayRead(da->coordinates,&coors);CHKERRQ(ierr); 6947c6ae99SBarry Smith } 7047c6ae99SBarry Smith #endif 7147c6ae99SBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 721575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 738135c375SStefano Zampini } else if (format == PETSC_VIEWER_ASCII_GLVIS) { 748135c375SStefano Zampini ierr = DMView_DA_GLVis(da,viewer);CHKERRQ(ierr); 75616157d6SJed Brown } else { 76616157d6SJed Brown ierr = DMView_DA_VTK(da,viewer);CHKERRQ(ierr); 7747c6ae99SBarry Smith } 7847c6ae99SBarry Smith } else if (isdraw) { 7947c6ae99SBarry Smith PetscDraw draw; 8047c6ae99SBarry Smith PetscReal ymin = -1.0,ymax = (PetscReal)dd->N; 8147c6ae99SBarry Smith PetscReal xmin = -1.0,xmax = (PetscReal)((dd->M+2)*dd->P),x,y,ycoord,xcoord; 828ea3bf28SBarry Smith PetscInt k,plane,base; 838ea3bf28SBarry Smith const PetscInt *idx; 8447c6ae99SBarry Smith char node[10]; 8547c6ae99SBarry Smith PetscBool isnull; 8647c6ae99SBarry Smith 8747c6ae99SBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 885b399a63SLisandro Dalcin ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 895b399a63SLisandro Dalcin if (isnull) PetscFunctionReturn(0); 9047c6ae99SBarry Smith 915b399a63SLisandro Dalcin ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr); 925b399a63SLisandro Dalcin ierr = PetscDrawClear(draw);CHKERRQ(ierr); 935b399a63SLisandro Dalcin ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr); 945b399a63SLisandro Dalcin 955b399a63SLisandro Dalcin ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr); 9647c6ae99SBarry Smith /* first processor draw all node lines */ 9747c6ae99SBarry Smith if (!rank) { 9847c6ae99SBarry Smith for (k=0; k<dd->P; k++) { 9947c6ae99SBarry Smith ymin = 0.0; ymax = (PetscReal)(dd->N - 1); 10047c6ae99SBarry Smith for (xmin=(PetscReal)(k*(dd->M+1)); xmin<(PetscReal)(dd->M+(k*(dd->M+1))); xmin++) { 10147c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr); 10247c6ae99SBarry Smith } 10347c6ae99SBarry Smith xmin = (PetscReal)(k*(dd->M+1)); xmax = xmin + (PetscReal)(dd->M - 1); 10447c6ae99SBarry Smith for (ymin=0; ymin<(PetscReal)dd->N; ymin++) { 10547c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr); 10647c6ae99SBarry Smith } 10747c6ae99SBarry Smith } 10847c6ae99SBarry Smith } 1095b399a63SLisandro Dalcin ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr); 1105b399a63SLisandro Dalcin ierr = PetscDrawFlush(draw);CHKERRQ(ierr); 11147c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 11247c6ae99SBarry Smith 1135b399a63SLisandro Dalcin ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr); 1145b399a63SLisandro Dalcin /*Go through and draw for each plane*/ 1155b399a63SLisandro Dalcin for (k=0; k<dd->P; k++) { 11647c6ae99SBarry Smith if ((k >= dd->zs) && (k < dd->ze)) { 11747c6ae99SBarry Smith /* draw my box */ 11847c6ae99SBarry Smith ymin = dd->ys; 11947c6ae99SBarry Smith ymax = dd->ye - 1; 12047c6ae99SBarry Smith xmin = dd->xs/dd->w + (dd->M+1)*k; 12147c6ae99SBarry Smith xmax =(dd->xe-1)/dd->w + (dd->M+1)*k; 12247c6ae99SBarry Smith 12347c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr); 12447c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 12547c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 12647c6ae99SBarry Smith ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 12747c6ae99SBarry Smith 12847c6ae99SBarry Smith xmin = dd->xs/dd->w; 12947c6ae99SBarry Smith xmax =(dd->xe-1)/dd->w; 13047c6ae99SBarry Smith 131832b7cebSLisandro Dalcin /* identify which processor owns the box */ 132832b7cebSLisandro Dalcin ierr = PetscSNPrintf(node,sizeof(node),"%d",(int)rank);CHKERRQ(ierr); 133832b7cebSLisandro Dalcin ierr = PetscDrawString(draw,xmin+(dd->M+1)*k+.2,ymin+.3,PETSC_DRAW_RED,node);CHKERRQ(ierr); 13447c6ae99SBarry Smith /* put in numbers*/ 13547c6ae99SBarry Smith base = (dd->base+(dd->xe-dd->xs)*(dd->ye-dd->ys)*(k-dd->zs))/dd->w; 13647c6ae99SBarry Smith for (y=ymin; y<=ymax; y++) { 13747c6ae99SBarry Smith for (x=xmin+(dd->M+1)*k; x<=xmax+(dd->M+1)*k; x++) { 138832b7cebSLisandro Dalcin ierr = PetscSNPrintf(node,sizeof(node),"%d",(int)base++);CHKERRQ(ierr); 13947c6ae99SBarry Smith ierr = PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);CHKERRQ(ierr); 14047c6ae99SBarry Smith } 14147c6ae99SBarry Smith } 14247c6ae99SBarry Smith 14347c6ae99SBarry Smith } 14447c6ae99SBarry Smith } 1455b399a63SLisandro Dalcin ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr); 1465b399a63SLisandro Dalcin ierr = PetscDrawFlush(draw);CHKERRQ(ierr); 14747c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 14847c6ae99SBarry Smith 1495b399a63SLisandro Dalcin ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr); 15047c6ae99SBarry Smith for (k=0-dd->s; k<dd->P+dd->s; k++) { 15147c6ae99SBarry Smith /* Go through and draw for each plane */ 15247c6ae99SBarry Smith if ((k >= dd->Zs) && (k < dd->Ze)) { 15347c6ae99SBarry Smith /* overlay ghost numbers, useful for error checking */ 154565245c5SBarry Smith base = (dd->Xe-dd->Xs)*(dd->Ye-dd->Ys)*(k-dd->Zs)/dd->w; 15545b6f7e9SBarry Smith ierr = ISLocalToGlobalMappingGetBlockIndices(da->ltogmap,&idx);CHKERRQ(ierr); 15647c6ae99SBarry Smith plane=k; 157565245c5SBarry Smith /* Keep z wrap around points on the drawing */ 1588865f1eaSKarl Rupp if (k<0) plane=dd->P+k; 1598865f1eaSKarl Rupp if (k>=dd->P) plane=k-dd->P; 16047c6ae99SBarry Smith ymin = dd->Ys; ymax = dd->Ye; 16147c6ae99SBarry Smith xmin = (dd->M+1)*plane*dd->w; 16247c6ae99SBarry Smith xmax = (dd->M+1)*plane*dd->w+dd->M*dd->w; 16347c6ae99SBarry Smith for (y=ymin; y<ymax; y++) { 16447c6ae99SBarry Smith for (x=xmin+dd->Xs; x<xmin+dd->Xe; x+=dd->w) { 1658ea3bf28SBarry Smith sprintf(node,"%d",(int)(idx[base])); 16647c6ae99SBarry Smith ycoord = y; 16747c6ae99SBarry Smith /*Keep y wrap around points on drawing */ 1688865f1eaSKarl Rupp if (y<0) ycoord = dd->N+y; 1698865f1eaSKarl Rupp if (y>=dd->N) ycoord = y-dd->N; 17047c6ae99SBarry Smith xcoord = x; /* Keep x wrap points on drawing */ 1718865f1eaSKarl Rupp if (x<xmin) xcoord = xmax - (xmin-x); 1728865f1eaSKarl Rupp if (x>=xmax) xcoord = xmin + (x-xmax); 17347c6ae99SBarry Smith ierr = PetscDrawString(draw,xcoord/dd->w,ycoord,PETSC_DRAW_BLUE,node);CHKERRQ(ierr); 174565245c5SBarry Smith base++; 17547c6ae99SBarry Smith } 17647c6ae99SBarry Smith } 17745b6f7e9SBarry Smith ierr = ISLocalToGlobalMappingRestoreBlockIndices(da->ltogmap,&idx);CHKERRQ(ierr); 17847c6ae99SBarry Smith } 17947c6ae99SBarry Smith } 1805b399a63SLisandro Dalcin ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr); 1815b399a63SLisandro Dalcin ierr = PetscDrawFlush(draw);CHKERRQ(ierr); 18247c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 183832b7cebSLisandro Dalcin ierr = PetscDrawSave(draw);CHKERRQ(ierr); 184c9493c35SStefano Zampini } else if (isglvis) { 185c9493c35SStefano Zampini ierr = DMView_DA_GLVis(da,viewer);CHKERRQ(ierr); 1869a42bb27SBarry Smith } else if (isbinary) { 1879a42bb27SBarry Smith ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr); 1889a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 1899a42bb27SBarry Smith } else if (ismatlab) { 1909a42bb27SBarry Smith ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr); 1919a42bb27SBarry Smith #endif 19211aeaf0aSBarry Smith } 19347c6ae99SBarry Smith PetscFunctionReturn(0); 19447c6ae99SBarry Smith } 19547c6ae99SBarry Smith 1967087cfbeSBarry Smith PetscErrorCode DMSetUp_DA_3D(DM da) 19747c6ae99SBarry Smith { 19847c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 19947c6ae99SBarry Smith const PetscInt M = dd->M; 20047c6ae99SBarry Smith const PetscInt N = dd->N; 20147c6ae99SBarry Smith const PetscInt P = dd->P; 20247c6ae99SBarry Smith PetscInt m = dd->m; 20347c6ae99SBarry Smith PetscInt n = dd->n; 20447c6ae99SBarry Smith PetscInt p = dd->p; 20547c6ae99SBarry Smith const PetscInt dof = dd->w; 20647c6ae99SBarry Smith const PetscInt s = dd->s; 207bff4a2f0SMatthew G. Knepley DMBoundaryType bx = dd->bx; 208bff4a2f0SMatthew G. Knepley DMBoundaryType by = dd->by; 209bff4a2f0SMatthew G. Knepley DMBoundaryType bz = dd->bz; 21019fd82e9SBarry Smith DMDAStencilType stencil_type = dd->stencil_type; 21147c6ae99SBarry Smith PetscInt *lx = dd->lx; 21247c6ae99SBarry Smith PetscInt *ly = dd->ly; 21347c6ae99SBarry Smith PetscInt *lz = dd->lz; 21447c6ae99SBarry Smith MPI_Comm comm; 21547c6ae99SBarry Smith PetscMPIInt rank,size; 216ce00eea3SSatish Balay PetscInt xs = 0,xe,ys = 0,ye,zs = 0,ze,x = 0,y = 0,z = 0; 217bd1fc5aeSBarry Smith PetscInt Xs,Xe,Ys,Ye,Zs,Ze,IXs,IXe,IYs,IYe,IZs,IZe,pm; 2188ea3bf28SBarry Smith PetscInt left,right,up,down,bottom,top,i,j,k,*idx,nn; 21947c6ae99SBarry Smith PetscInt n0,n1,n2,n3,n4,n5,n6,n7,n8,n9,n10,n11,n12,n14; 22047c6ae99SBarry Smith PetscInt n15,n16,n17,n18,n19,n20,n21,n22,n23,n24,n25,n26; 221db87c5ecSEthan Coon PetscInt *bases,*ldims,base,x_t,y_t,z_t,s_t,count,s_x,s_y,s_z; 22247c6ae99SBarry Smith PetscInt sn0 = 0,sn1 = 0,sn2 = 0,sn3 = 0,sn5 = 0,sn6 = 0,sn7 = 0; 22347c6ae99SBarry Smith PetscInt sn8 = 0,sn9 = 0,sn11 = 0,sn15 = 0,sn24 = 0,sn25 = 0,sn26 = 0; 22447c6ae99SBarry Smith PetscInt sn17 = 0,sn18 = 0,sn19 = 0,sn20 = 0,sn21 = 0,sn23 = 0; 22547c6ae99SBarry Smith Vec local,global; 226bd1fc5aeSBarry Smith VecScatter gtol; 22745b6f7e9SBarry Smith IS to,from; 2286f951b95Secoon PetscBool twod; 22947c6ae99SBarry Smith PetscErrorCode ierr; 23047c6ae99SBarry Smith 2316f951b95Secoon 23247c6ae99SBarry Smith PetscFunctionBegin; 233bff4a2f0SMatthew G. Knepley if (stencil_type == DMDA_STENCIL_BOX && (bx == DM_BOUNDARY_MIRROR || by == DM_BOUNDARY_MIRROR || bz == DM_BOUNDARY_MIRROR)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Mirror boundary and box stencil"); 23440088605SBarry Smith if ((bx == DM_BOUNDARY_MIRROR) || (by == DM_BOUNDARY_MIRROR) || (bz == DM_BOUNDARY_MIRROR)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Mirror boundary not supported yet in 3d"); 23547c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject) da, &comm);CHKERRQ(ierr); 2363855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES) 237bafee8b4SSatish Balay if (((PetscInt64) M)*((PetscInt64) N)*((PetscInt64) P)*((PetscInt64) dof) > (PetscInt64) 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); 2383855c12bSBarry Smith #endif 2393855c12bSBarry Smith 24047c6ae99SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 24147c6ae99SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 24247c6ae99SBarry Smith 24347c6ae99SBarry Smith if (m != PETSC_DECIDE) { 24447c6ae99SBarry Smith if (m < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m); 24547c6ae99SBarry Smith else if (m > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size); 24647c6ae99SBarry Smith } 24747c6ae99SBarry Smith if (n != PETSC_DECIDE) { 24847c6ae99SBarry Smith if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n); 24947c6ae99SBarry Smith else if (n > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size); 25047c6ae99SBarry Smith } 25147c6ae99SBarry Smith if (p != PETSC_DECIDE) { 25247c6ae99SBarry Smith if (p < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Z direction: %D",p); 25347c6ae99SBarry Smith else if (p > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Z direction: %D %d",p,size); 25447c6ae99SBarry Smith } 2550716a85fSBarry 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); 25647c6ae99SBarry Smith 25747c6ae99SBarry Smith /* Partition the array among the processors */ 25847c6ae99SBarry Smith if (m == PETSC_DECIDE && n != PETSC_DECIDE && p != PETSC_DECIDE) { 25947c6ae99SBarry Smith m = size/(n*p); 26047c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) { 26147c6ae99SBarry Smith n = size/(m*p); 26247c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) { 26347c6ae99SBarry Smith p = size/(m*n); 26447c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) { 26547c6ae99SBarry Smith /* try for squarish distribution */ 266369cc0aeSBarry Smith m = (int)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)N*p))); 26747c6ae99SBarry Smith if (!m) m = 1; 26847c6ae99SBarry Smith while (m > 0) { 26947c6ae99SBarry Smith n = size/(m*p); 27047c6ae99SBarry Smith if (m*n*p == size) break; 27147c6ae99SBarry Smith m--; 27247c6ae99SBarry Smith } 27347c6ae99SBarry Smith if (!m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad p value: p = %D",p); 27447c6ae99SBarry Smith if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;} 27547c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) { 27647c6ae99SBarry Smith /* try for squarish distribution */ 277369cc0aeSBarry Smith m = (int)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n))); 27847c6ae99SBarry Smith if (!m) m = 1; 27947c6ae99SBarry Smith while (m > 0) { 28047c6ae99SBarry Smith p = size/(m*n); 28147c6ae99SBarry Smith if (m*n*p == size) break; 28247c6ae99SBarry Smith m--; 28347c6ae99SBarry Smith } 28447c6ae99SBarry Smith if (!m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad n value: n = %D",n); 28547c6ae99SBarry Smith if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;} 28647c6ae99SBarry Smith } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) { 28747c6ae99SBarry Smith /* try for squarish distribution */ 288369cc0aeSBarry Smith n = (int)(0.5 + PetscSqrtReal(((PetscReal)N)*((PetscReal)size)/((PetscReal)P*m))); 28947c6ae99SBarry Smith if (!n) n = 1; 29047c6ae99SBarry Smith while (n > 0) { 29147c6ae99SBarry Smith p = size/(m*n); 29247c6ae99SBarry Smith if (m*n*p == size) break; 29347c6ae99SBarry Smith n--; 29447c6ae99SBarry Smith } 29547c6ae99SBarry Smith if (!n) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad m value: m = %D",n); 29647c6ae99SBarry Smith if (N > P && n < p) {PetscInt _n = n; n = p; p = _n;} 29747c6ae99SBarry Smith } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) { 29847c6ae99SBarry Smith /* try for squarish distribution */ 299369cc0aeSBarry Smith n = (PetscInt)(0.5 + PetscPowReal(((PetscReal)N*N)*((PetscReal)size)/((PetscReal)P*M),(PetscReal)(1./3.))); 30047c6ae99SBarry Smith if (!n) n = 1; 30147c6ae99SBarry Smith while (n > 0) { 30247c6ae99SBarry Smith pm = size/n; 30347c6ae99SBarry Smith if (n*pm == size) break; 30447c6ae99SBarry Smith n--; 30547c6ae99SBarry Smith } 30647c6ae99SBarry Smith if (!n) n = 1; 307369cc0aeSBarry Smith m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n))); 30847c6ae99SBarry Smith if (!m) m = 1; 30947c6ae99SBarry Smith while (m > 0) { 31047c6ae99SBarry Smith p = size/(m*n); 31147c6ae99SBarry Smith if (m*n*p == size) break; 31247c6ae99SBarry Smith m--; 31347c6ae99SBarry Smith } 31447c6ae99SBarry Smith if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;} 315ce94432eSBarry Smith } else if (m*n*p != size) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition"); 31647c6ae99SBarry Smith 317ce94432eSBarry Smith if (m*n*p != size) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_PLIB,"Could not find good partition"); 318ce94432eSBarry Smith if (M < m) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m); 319ce94432eSBarry Smith if (N < n) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n); 320ce94432eSBarry Smith if (P < p) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Partition in z direction is too fine! %D %D",P,p); 32147c6ae99SBarry Smith 32247c6ae99SBarry Smith /* 32347c6ae99SBarry Smith Determine locally owned region 32447c6ae99SBarry Smith [x, y, or z]s is the first local node number, [x, y, z] is the number of local nodes 32547c6ae99SBarry Smith */ 32647c6ae99SBarry Smith 32747c6ae99SBarry Smith if (!lx) { 328785e854fSJed Brown ierr = PetscMalloc1(m, &dd->lx);CHKERRQ(ierr); 32947c6ae99SBarry Smith lx = dd->lx; 3308865f1eaSKarl Rupp for (i=0; i<m; i++) lx[i] = M/m + ((M % m) > (i % m)); 33147c6ae99SBarry Smith } 33247c6ae99SBarry Smith x = lx[rank % m]; 33347c6ae99SBarry Smith xs = 0; 3348865f1eaSKarl Rupp for (i=0; i<(rank%m); i++) xs += lx[i]; 335bff4a2f0SMatthew G. Knepley if ((x < s) && ((m > 1) || (bx == DM_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); 33647c6ae99SBarry Smith 33747c6ae99SBarry Smith if (!ly) { 338785e854fSJed Brown ierr = PetscMalloc1(n, &dd->ly);CHKERRQ(ierr); 33947c6ae99SBarry Smith ly = dd->ly; 3408865f1eaSKarl Rupp for (i=0; i<n; i++) ly[i] = N/n + ((N % n) > (i % n)); 34147c6ae99SBarry Smith } 34247c6ae99SBarry Smith y = ly[(rank % (m*n))/m]; 343bff4a2f0SMatthew G. Knepley if ((y < s) && ((n > 1) || (by == DM_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); 34430729d88SBarry Smith 34547c6ae99SBarry Smith ys = 0; 3468865f1eaSKarl Rupp for (i=0; i<(rank % (m*n))/m; i++) ys += ly[i]; 34747c6ae99SBarry Smith 34847c6ae99SBarry Smith if (!lz) { 349785e854fSJed Brown ierr = PetscMalloc1(p, &dd->lz);CHKERRQ(ierr); 35047c6ae99SBarry Smith lz = dd->lz; 3518865f1eaSKarl Rupp for (i=0; i<p; i++) lz[i] = P/p + ((P % p) > (i % p)); 35247c6ae99SBarry Smith } 35347c6ae99SBarry Smith z = lz[rank/(m*n)]; 354bcea557cSEthan Coon 355fdc81ce8SEthan Coon /* note this is different than x- and y-, as we will handle as an important special 356bff4a2f0SMatthew G. Knepley case when p=P=1 and DM_BOUNDARY_PERIODIC and s > z. This is to deal with 2D problems 357fdc81ce8SEthan Coon in a 3D code. Additional code for this case is noted with "2d case" comments */ 3586f951b95Secoon twod = PETSC_FALSE; 3598865f1eaSKarl Rupp if (P == 1) twod = PETSC_TRUE; 360bff4a2f0SMatthew G. Knepley else if ((z < s) && ((p > 1) || (bz == DM_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); 36147c6ae99SBarry Smith zs = 0; 3628865f1eaSKarl Rupp for (i=0; i<(rank/(m*n)); i++) zs += lz[i]; 36347c6ae99SBarry Smith ye = ys + y; 36447c6ae99SBarry Smith xe = xs + x; 36547c6ae99SBarry Smith ze = zs + z; 36647c6ae99SBarry Smith 36788661749SPeter Brune /* determine ghost region (Xs) and region scattered into (IXs) */ 368d9c9ebe5SPeter Brune if (xs-s > 0) { 369d9c9ebe5SPeter Brune Xs = xs - s; IXs = xs - s; 37088661749SPeter Brune } else { 3718865f1eaSKarl Rupp if (bx) Xs = xs - s; 3728865f1eaSKarl Rupp else Xs = 0; 37388661749SPeter Brune IXs = 0; 37488661749SPeter Brune } 375d9c9ebe5SPeter Brune if (xe+s <= M) { 376d9c9ebe5SPeter Brune Xe = xe + s; IXe = xe + s; 37788661749SPeter Brune } else { 37888661749SPeter Brune if (bx) { 379d9c9ebe5SPeter Brune Xs = xs - s; Xe = xe + s; 3808865f1eaSKarl Rupp } else Xe = M; 38188661749SPeter Brune IXe = M; 38288661749SPeter Brune } 38347c6ae99SBarry Smith 384bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) { 385d9c9ebe5SPeter Brune IXs = xs - s; 386d9c9ebe5SPeter Brune IXe = xe + s; 387d9c9ebe5SPeter Brune Xs = xs - s; 388d9c9ebe5SPeter Brune Xe = xe + s; 38988661749SPeter Brune } 39088661749SPeter Brune 391d9c9ebe5SPeter Brune if (ys-s > 0) { 392d9c9ebe5SPeter Brune Ys = ys - s; IYs = ys - s; 39388661749SPeter Brune } else { 3948865f1eaSKarl Rupp if (by) Ys = ys - s; 3958865f1eaSKarl Rupp else Ys = 0; 39688661749SPeter Brune IYs = 0; 39788661749SPeter Brune } 398d9c9ebe5SPeter Brune if (ye+s <= N) { 399d9c9ebe5SPeter Brune Ye = ye + s; IYe = ye + s; 40088661749SPeter Brune } else { 4018865f1eaSKarl Rupp if (by) Ye = ye + s; 4028865f1eaSKarl Rupp else Ye = N; 40388661749SPeter Brune IYe = N; 40488661749SPeter Brune } 40588661749SPeter Brune 406bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC || by == DM_BOUNDARY_MIRROR) { 407d9c9ebe5SPeter Brune IYs = ys - s; 408d9c9ebe5SPeter Brune IYe = ye + s; 409d9c9ebe5SPeter Brune Ys = ys - s; 410d9c9ebe5SPeter Brune Ye = ye + s; 41188661749SPeter Brune } 41288661749SPeter Brune 413d9c9ebe5SPeter Brune if (zs-s > 0) { 414d9c9ebe5SPeter Brune Zs = zs - s; IZs = zs - s; 41588661749SPeter Brune } else { 4168865f1eaSKarl Rupp if (bz) Zs = zs - s; 4178865f1eaSKarl Rupp else Zs = 0; 41888661749SPeter Brune IZs = 0; 41988661749SPeter Brune } 420d9c9ebe5SPeter Brune if (ze+s <= P) { 421d9c9ebe5SPeter Brune Ze = ze + s; IZe = ze + s; 42288661749SPeter Brune } else { 4238865f1eaSKarl Rupp if (bz) Ze = ze + s; 4248865f1eaSKarl Rupp else Ze = P; 42588661749SPeter Brune IZe = P; 42688661749SPeter Brune } 42788661749SPeter Brune 428bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_PERIODIC || bz == DM_BOUNDARY_MIRROR) { 429d9c9ebe5SPeter Brune IZs = zs - s; 430d9c9ebe5SPeter Brune IZe = ze + s; 431d9c9ebe5SPeter Brune Zs = zs - s; 432d9c9ebe5SPeter Brune Ze = ze + s; 43388661749SPeter Brune } 43447c6ae99SBarry Smith 43547c6ae99SBarry Smith /* Resize all X parameters to reflect w */ 436d9c9ebe5SPeter Brune s_x = s; 437d9c9ebe5SPeter Brune s_y = s; 438d9c9ebe5SPeter Brune s_z = s; 43947c6ae99SBarry Smith 44047c6ae99SBarry Smith /* determine starting point of each processor */ 44147c6ae99SBarry Smith nn = x*y*z; 442dcca6d9dSJed Brown ierr = PetscMalloc2(size+1,&bases,size,&ldims);CHKERRQ(ierr); 44347c6ae99SBarry Smith ierr = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr); 44447c6ae99SBarry Smith bases[0] = 0; 4458865f1eaSKarl Rupp for (i=1; i<=size; i++) bases[i] = ldims[i-1]; 4468865f1eaSKarl Rupp for (i=1; i<=size; i++) bases[i] += bases[i-1]; 447ce00eea3SSatish Balay base = bases[rank]*dof; 44847c6ae99SBarry Smith 44947c6ae99SBarry Smith /* allocate the base parallel and sequential vectors */ 450ce00eea3SSatish Balay dd->Nlocal = x*y*z*dof; 451b1fb7eb7SBarry Smith ierr = VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,NULL,&global);CHKERRQ(ierr); 452ce00eea3SSatish Balay dd->nlocal = (Xe-Xs)*(Ye-Ys)*(Ze-Zs)*dof; 453b1fb7eb7SBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,NULL,&local);CHKERRQ(ierr); 45447c6ae99SBarry Smith 45547c6ae99SBarry Smith /* generate appropriate vector scatters */ 45647c6ae99SBarry Smith /* local to global inserts non-ghost point region into global */ 45702fe608eSBarry Smith ierr = PetscMalloc1((IXe-IXs)*(IYe-IYs)*(IZe-IZs),&idx);CHKERRQ(ierr); 458ce00eea3SSatish Balay left = xs - Xs; right = left + x; 45947c6ae99SBarry Smith bottom = ys - Ys; top = bottom + y; 46047c6ae99SBarry Smith down = zs - Zs; up = down + z; 46147c6ae99SBarry Smith count = 0; 46247c6ae99SBarry Smith for (i=down; i<up; i++) { 46347c6ae99SBarry Smith for (j=bottom; j<top; j++) { 464ce00eea3SSatish Balay for (k=left; k<right; k++) { 465ce00eea3SSatish Balay idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 46647c6ae99SBarry Smith } 46747c6ae99SBarry Smith } 46847c6ae99SBarry Smith } 46947c6ae99SBarry Smith 470ce00eea3SSatish Balay /* global to local must include ghost points within the domain, 471ce00eea3SSatish Balay but not ghost points outside the domain that aren't periodic */ 472d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_BOX) { 473ce00eea3SSatish Balay left = IXs - Xs; right = left + (IXe-IXs); 474ce00eea3SSatish Balay bottom = IYs - Ys; top = bottom + (IYe-IYs); 475ce00eea3SSatish Balay down = IZs - Zs; up = down + (IZe-IZs); 476ce00eea3SSatish Balay count = 0; 477ce00eea3SSatish Balay for (i=down; i<up; i++) { 478ce00eea3SSatish Balay for (j=bottom; j<top; j++) { 479ce00eea3SSatish Balay for (k=left; k<right; k++) { 480ce00eea3SSatish Balay idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 481ce00eea3SSatish Balay } 482ce00eea3SSatish Balay } 483ce00eea3SSatish Balay } 484ce00eea3SSatish Balay ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr); 48547c6ae99SBarry Smith } else { 48647c6ae99SBarry Smith /* This is way ugly! We need to list the funny cross type region */ 487ce00eea3SSatish Balay left = xs - Xs; right = left + x; 48847c6ae99SBarry Smith bottom = ys - Ys; top = bottom + y; 48947c6ae99SBarry Smith down = zs - Zs; up = down + z; 49047c6ae99SBarry Smith count = 0; 491ce00eea3SSatish Balay /* the bottom chunck */ 492ce00eea3SSatish Balay for (i=(IZs-Zs); i<down; i++) { 49347c6ae99SBarry Smith for (j=bottom; j<top; j++) { 494ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 49547c6ae99SBarry Smith } 49647c6ae99SBarry Smith } 49747c6ae99SBarry Smith /* the middle piece */ 49847c6ae99SBarry Smith for (i=down; i<up; i++) { 49947c6ae99SBarry Smith /* front */ 500ce00eea3SSatish Balay for (j=(IYs-Ys); j<bottom; j++) { 501ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 50247c6ae99SBarry Smith } 50347c6ae99SBarry Smith /* middle */ 50447c6ae99SBarry Smith for (j=bottom; j<top; j++) { 505ce00eea3SSatish Balay for (k=IXs-Xs; k<IXe-Xs; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 50647c6ae99SBarry Smith } 50747c6ae99SBarry Smith /* back */ 508ce00eea3SSatish Balay for (j=top; j<top+IYe-ye; 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 top piece */ 513ce00eea3SSatish Balay for (i=up; i<up+IZe-ze; i++) { 51447c6ae99SBarry Smith for (j=bottom; j<top; j++) { 515ce00eea3SSatish Balay for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k; 51647c6ae99SBarry Smith } 51747c6ae99SBarry Smith } 51847c6ae99SBarry Smith ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr); 51947c6ae99SBarry Smith } 52047c6ae99SBarry Smith 52147c6ae99SBarry Smith /* determine who lies on each side of use stored in n24 n25 n26 52247c6ae99SBarry Smith n21 n22 n23 52347c6ae99SBarry Smith n18 n19 n20 52447c6ae99SBarry Smith 52547c6ae99SBarry Smith n15 n16 n17 52647c6ae99SBarry Smith n12 n14 52747c6ae99SBarry Smith n9 n10 n11 52847c6ae99SBarry Smith 52947c6ae99SBarry Smith n6 n7 n8 53047c6ae99SBarry Smith n3 n4 n5 53147c6ae99SBarry Smith n0 n1 n2 53247c6ae99SBarry Smith */ 53347c6ae99SBarry Smith 53447c6ae99SBarry Smith /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */ 53547c6ae99SBarry Smith /* Assume Nodes are Internal to the Cube */ 53647c6ae99SBarry Smith n0 = rank - m*n - m - 1; 53747c6ae99SBarry Smith n1 = rank - m*n - m; 53847c6ae99SBarry Smith n2 = rank - m*n - m + 1; 53947c6ae99SBarry Smith n3 = rank - m*n -1; 54047c6ae99SBarry Smith n4 = rank - m*n; 54147c6ae99SBarry Smith n5 = rank - m*n + 1; 54247c6ae99SBarry Smith n6 = rank - m*n + m - 1; 54347c6ae99SBarry Smith n7 = rank - m*n + m; 54447c6ae99SBarry Smith n8 = rank - m*n + m + 1; 54547c6ae99SBarry Smith 54647c6ae99SBarry Smith n9 = rank - m - 1; 54747c6ae99SBarry Smith n10 = rank - m; 54847c6ae99SBarry Smith n11 = rank - m + 1; 54947c6ae99SBarry Smith n12 = rank - 1; 55047c6ae99SBarry Smith n14 = rank + 1; 55147c6ae99SBarry Smith n15 = rank + m - 1; 55247c6ae99SBarry Smith n16 = rank + m; 55347c6ae99SBarry Smith n17 = rank + m + 1; 55447c6ae99SBarry Smith 55547c6ae99SBarry Smith n18 = rank + m*n - m - 1; 55647c6ae99SBarry Smith n19 = rank + m*n - m; 55747c6ae99SBarry Smith n20 = rank + m*n - m + 1; 55847c6ae99SBarry Smith n21 = rank + m*n - 1; 55947c6ae99SBarry Smith n22 = rank + m*n; 56047c6ae99SBarry Smith n23 = rank + m*n + 1; 56147c6ae99SBarry Smith n24 = rank + m*n + m - 1; 56247c6ae99SBarry Smith n25 = rank + m*n + m; 56347c6ae99SBarry Smith n26 = rank + m*n + m + 1; 56447c6ae99SBarry Smith 56547c6ae99SBarry Smith /* Assume Pieces are on Faces of Cube */ 56647c6ae99SBarry Smith 56747c6ae99SBarry Smith if (xs == 0) { /* First assume not corner or edge */ 56847c6ae99SBarry Smith n0 = rank -1 - (m*n); 56947c6ae99SBarry Smith n3 = rank + m -1 - (m*n); 57047c6ae99SBarry Smith n6 = rank + 2*m -1 - (m*n); 57147c6ae99SBarry Smith n9 = rank -1; 57247c6ae99SBarry Smith n12 = rank + m -1; 57347c6ae99SBarry Smith n15 = rank + 2*m -1; 57447c6ae99SBarry Smith n18 = rank -1 + (m*n); 57547c6ae99SBarry Smith n21 = rank + m -1 + (m*n); 57647c6ae99SBarry Smith n24 = rank + 2*m -1 + (m*n); 57747c6ae99SBarry Smith } 57847c6ae99SBarry Smith 579ce00eea3SSatish Balay if (xe == M) { /* First assume not corner or edge */ 58047c6ae99SBarry Smith n2 = rank -2*m +1 - (m*n); 58147c6ae99SBarry Smith n5 = rank - m +1 - (m*n); 58247c6ae99SBarry Smith n8 = rank +1 - (m*n); 58347c6ae99SBarry Smith n11 = rank -2*m +1; 58447c6ae99SBarry Smith n14 = rank - m +1; 58547c6ae99SBarry Smith n17 = rank +1; 58647c6ae99SBarry Smith n20 = rank -2*m +1 + (m*n); 58747c6ae99SBarry Smith n23 = rank - m +1 + (m*n); 58847c6ae99SBarry Smith n26 = rank +1 + (m*n); 58947c6ae99SBarry Smith } 59047c6ae99SBarry Smith 59147c6ae99SBarry Smith if (ys==0) { /* First assume not corner or edge */ 59247c6ae99SBarry Smith n0 = rank + m * (n-1) -1 - (m*n); 59347c6ae99SBarry Smith n1 = rank + m * (n-1) - (m*n); 59447c6ae99SBarry Smith n2 = rank + m * (n-1) +1 - (m*n); 59547c6ae99SBarry Smith n9 = rank + m * (n-1) -1; 59647c6ae99SBarry Smith n10 = rank + m * (n-1); 59747c6ae99SBarry Smith n11 = rank + m * (n-1) +1; 59847c6ae99SBarry Smith n18 = rank + m * (n-1) -1 + (m*n); 59947c6ae99SBarry Smith n19 = rank + m * (n-1) + (m*n); 60047c6ae99SBarry Smith n20 = rank + m * (n-1) +1 + (m*n); 60147c6ae99SBarry Smith } 60247c6ae99SBarry Smith 60347c6ae99SBarry Smith if (ye == N) { /* First assume not corner or edge */ 60447c6ae99SBarry Smith n6 = rank - m * (n-1) -1 - (m*n); 60547c6ae99SBarry Smith n7 = rank - m * (n-1) - (m*n); 60647c6ae99SBarry Smith n8 = rank - m * (n-1) +1 - (m*n); 60747c6ae99SBarry Smith n15 = rank - m * (n-1) -1; 60847c6ae99SBarry Smith n16 = rank - m * (n-1); 60947c6ae99SBarry Smith n17 = rank - m * (n-1) +1; 61047c6ae99SBarry Smith n24 = rank - m * (n-1) -1 + (m*n); 61147c6ae99SBarry Smith n25 = rank - m * (n-1) + (m*n); 61247c6ae99SBarry Smith n26 = rank - m * (n-1) +1 + (m*n); 61347c6ae99SBarry Smith } 61447c6ae99SBarry Smith 61547c6ae99SBarry Smith if (zs == 0) { /* First assume not corner or edge */ 61647c6ae99SBarry Smith n0 = size - (m*n) + rank - m - 1; 61747c6ae99SBarry Smith n1 = size - (m*n) + rank - m; 61847c6ae99SBarry Smith n2 = size - (m*n) + rank - m + 1; 61947c6ae99SBarry Smith n3 = size - (m*n) + rank - 1; 62047c6ae99SBarry Smith n4 = size - (m*n) + rank; 62147c6ae99SBarry Smith n5 = size - (m*n) + rank + 1; 62247c6ae99SBarry Smith n6 = size - (m*n) + rank + m - 1; 62347c6ae99SBarry Smith n7 = size - (m*n) + rank + m; 62447c6ae99SBarry Smith n8 = size - (m*n) + rank + m + 1; 62547c6ae99SBarry Smith } 62647c6ae99SBarry Smith 62747c6ae99SBarry Smith if (ze == P) { /* First assume not corner or edge */ 62847c6ae99SBarry Smith n18 = (m*n) - (size-rank) - m - 1; 62947c6ae99SBarry Smith n19 = (m*n) - (size-rank) - m; 63047c6ae99SBarry Smith n20 = (m*n) - (size-rank) - m + 1; 63147c6ae99SBarry Smith n21 = (m*n) - (size-rank) - 1; 63247c6ae99SBarry Smith n22 = (m*n) - (size-rank); 63347c6ae99SBarry Smith n23 = (m*n) - (size-rank) + 1; 63447c6ae99SBarry Smith n24 = (m*n) - (size-rank) + m - 1; 63547c6ae99SBarry Smith n25 = (m*n) - (size-rank) + m; 63647c6ae99SBarry Smith n26 = (m*n) - (size-rank) + m + 1; 63747c6ae99SBarry Smith } 63847c6ae99SBarry Smith 63947c6ae99SBarry Smith if ((xs==0) && (zs==0)) { /* Assume an edge, not corner */ 64047c6ae99SBarry Smith n0 = size - m*n + rank + m-1 - m; 64147c6ae99SBarry Smith n3 = size - m*n + rank + m-1; 64247c6ae99SBarry Smith n6 = size - m*n + rank + m-1 + m; 64347c6ae99SBarry Smith } 64447c6ae99SBarry Smith 64547c6ae99SBarry Smith if ((xs==0) && (ze==P)) { /* Assume an edge, not corner */ 64647c6ae99SBarry Smith n18 = m*n - (size - rank) + m-1 - m; 64747c6ae99SBarry Smith n21 = m*n - (size - rank) + m-1; 64847c6ae99SBarry Smith n24 = m*n - (size - rank) + m-1 + m; 64947c6ae99SBarry Smith } 65047c6ae99SBarry Smith 65147c6ae99SBarry Smith if ((xs==0) && (ys==0)) { /* Assume an edge, not corner */ 65247c6ae99SBarry Smith n0 = rank + m*n -1 - m*n; 65347c6ae99SBarry Smith n9 = rank + m*n -1; 65447c6ae99SBarry Smith n18 = rank + m*n -1 + m*n; 65547c6ae99SBarry Smith } 65647c6ae99SBarry Smith 65747c6ae99SBarry Smith if ((xs==0) && (ye==N)) { /* Assume an edge, not corner */ 65847c6ae99SBarry Smith n6 = rank - m*(n-1) + m-1 - m*n; 65947c6ae99SBarry Smith n15 = rank - m*(n-1) + m-1; 66047c6ae99SBarry Smith n24 = rank - m*(n-1) + m-1 + m*n; 66147c6ae99SBarry Smith } 66247c6ae99SBarry Smith 663ce00eea3SSatish Balay if ((xe==M) && (zs==0)) { /* Assume an edge, not corner */ 66447c6ae99SBarry Smith n2 = size - (m*n-rank) - (m-1) - m; 66547c6ae99SBarry Smith n5 = size - (m*n-rank) - (m-1); 66647c6ae99SBarry Smith n8 = size - (m*n-rank) - (m-1) + m; 66747c6ae99SBarry Smith } 66847c6ae99SBarry Smith 669ce00eea3SSatish Balay if ((xe==M) && (ze==P)) { /* Assume an edge, not corner */ 67047c6ae99SBarry Smith n20 = m*n - (size - rank) - (m-1) - m; 67147c6ae99SBarry Smith n23 = m*n - (size - rank) - (m-1); 67247c6ae99SBarry Smith n26 = m*n - (size - rank) - (m-1) + m; 67347c6ae99SBarry Smith } 67447c6ae99SBarry Smith 675ce00eea3SSatish Balay if ((xe==M) && (ys==0)) { /* Assume an edge, not corner */ 67647c6ae99SBarry Smith n2 = rank + m*(n-1) - (m-1) - m*n; 67747c6ae99SBarry Smith n11 = rank + m*(n-1) - (m-1); 67847c6ae99SBarry Smith n20 = rank + m*(n-1) - (m-1) + m*n; 67947c6ae99SBarry Smith } 68047c6ae99SBarry Smith 681ce00eea3SSatish Balay if ((xe==M) && (ye==N)) { /* Assume an edge, not corner */ 68247c6ae99SBarry Smith n8 = rank - m*n +1 - m*n; 68347c6ae99SBarry Smith n17 = rank - m*n +1; 68447c6ae99SBarry Smith n26 = rank - m*n +1 + m*n; 68547c6ae99SBarry Smith } 68647c6ae99SBarry Smith 68747c6ae99SBarry Smith if ((ys==0) && (zs==0)) { /* Assume an edge, not corner */ 68847c6ae99SBarry Smith n0 = size - m + rank -1; 68947c6ae99SBarry Smith n1 = size - m + rank; 69047c6ae99SBarry Smith n2 = size - m + rank +1; 69147c6ae99SBarry Smith } 69247c6ae99SBarry Smith 69347c6ae99SBarry Smith if ((ys==0) && (ze==P)) { /* Assume an edge, not corner */ 69447c6ae99SBarry Smith n18 = m*n - (size - rank) + m*(n-1) -1; 69547c6ae99SBarry Smith n19 = m*n - (size - rank) + m*(n-1); 69647c6ae99SBarry Smith n20 = m*n - (size - rank) + m*(n-1) +1; 69747c6ae99SBarry Smith } 69847c6ae99SBarry Smith 69947c6ae99SBarry Smith if ((ye==N) && (zs==0)) { /* Assume an edge, not corner */ 70047c6ae99SBarry Smith n6 = size - (m*n-rank) - m * (n-1) -1; 70147c6ae99SBarry Smith n7 = size - (m*n-rank) - m * (n-1); 70247c6ae99SBarry Smith n8 = size - (m*n-rank) - m * (n-1) +1; 70347c6ae99SBarry Smith } 70447c6ae99SBarry Smith 70547c6ae99SBarry Smith if ((ye==N) && (ze==P)) { /* Assume an edge, not corner */ 70647c6ae99SBarry Smith n24 = rank - (size-m) -1; 70747c6ae99SBarry Smith n25 = rank - (size-m); 70847c6ae99SBarry Smith n26 = rank - (size-m) +1; 70947c6ae99SBarry Smith } 71047c6ae99SBarry Smith 71147c6ae99SBarry Smith /* Check for Corners */ 7128865f1eaSKarl Rupp if ((xs==0) && (ys==0) && (zs==0)) n0 = size -1; 7138865f1eaSKarl Rupp if ((xs==0) && (ys==0) && (ze==P)) n18 = m*n-1; 7148865f1eaSKarl Rupp if ((xs==0) && (ye==N) && (zs==0)) n6 = (size-1)-m*(n-1); 7158865f1eaSKarl Rupp if ((xs==0) && (ye==N) && (ze==P)) n24 = m-1; 7168865f1eaSKarl Rupp if ((xe==M) && (ys==0) && (zs==0)) n2 = size-m; 7178865f1eaSKarl Rupp if ((xe==M) && (ys==0) && (ze==P)) n20 = m*n-m; 7188865f1eaSKarl Rupp if ((xe==M) && (ye==N) && (zs==0)) n8 = size-m*n; 7198865f1eaSKarl Rupp if ((xe==M) && (ye==N) && (ze==P)) n26 = 0; 72047c6ae99SBarry Smith 72147c6ae99SBarry Smith /* Check for when not X,Y, and Z Periodic */ 72247c6ae99SBarry Smith 72347c6ae99SBarry Smith /* If not X periodic */ 724bff4a2f0SMatthew G. Knepley if (bx != DM_BOUNDARY_PERIODIC) { 7258865f1eaSKarl Rupp if (xs==0) n0 = n3 = n6 = n9 = n12 = n15 = n18 = n21 = n24 = -2; 7268865f1eaSKarl Rupp if (xe==M) n2 = n5 = n8 = n11 = n14 = n17 = n20 = n23 = n26 = -2; 72747c6ae99SBarry Smith } 72847c6ae99SBarry Smith 72947c6ae99SBarry Smith /* If not Y periodic */ 730bff4a2f0SMatthew G. Knepley if (by != DM_BOUNDARY_PERIODIC) { 7318865f1eaSKarl Rupp if (ys==0) n0 = n1 = n2 = n9 = n10 = n11 = n18 = n19 = n20 = -2; 7328865f1eaSKarl Rupp if (ye==N) n6 = n7 = n8 = n15 = n16 = n17 = n24 = n25 = n26 = -2; 73347c6ae99SBarry Smith } 73447c6ae99SBarry Smith 73547c6ae99SBarry Smith /* If not Z periodic */ 736bff4a2f0SMatthew G. Knepley if (bz != DM_BOUNDARY_PERIODIC) { 7378865f1eaSKarl Rupp if (zs==0) n0 = n1 = n2 = n3 = n4 = n5 = n6 = n7 = n8 = -2; 7388865f1eaSKarl Rupp if (ze==P) n18 = n19 = n20 = n21 = n22 = n23 = n24 = n25 = n26 = -2; 73947c6ae99SBarry Smith } 74047c6ae99SBarry Smith 741785e854fSJed Brown ierr = PetscMalloc1(27,&dd->neighbors);CHKERRQ(ierr); 7428865f1eaSKarl Rupp 74347c6ae99SBarry Smith dd->neighbors[0] = n0; 74447c6ae99SBarry Smith dd->neighbors[1] = n1; 74547c6ae99SBarry Smith dd->neighbors[2] = n2; 74647c6ae99SBarry Smith dd->neighbors[3] = n3; 74747c6ae99SBarry Smith dd->neighbors[4] = n4; 74847c6ae99SBarry Smith dd->neighbors[5] = n5; 74947c6ae99SBarry Smith dd->neighbors[6] = n6; 75047c6ae99SBarry Smith dd->neighbors[7] = n7; 75147c6ae99SBarry Smith dd->neighbors[8] = n8; 75247c6ae99SBarry Smith dd->neighbors[9] = n9; 75347c6ae99SBarry Smith dd->neighbors[10] = n10; 75447c6ae99SBarry Smith dd->neighbors[11] = n11; 75547c6ae99SBarry Smith dd->neighbors[12] = n12; 75647c6ae99SBarry Smith dd->neighbors[13] = rank; 75747c6ae99SBarry Smith dd->neighbors[14] = n14; 75847c6ae99SBarry Smith dd->neighbors[15] = n15; 75947c6ae99SBarry Smith dd->neighbors[16] = n16; 76047c6ae99SBarry Smith dd->neighbors[17] = n17; 76147c6ae99SBarry Smith dd->neighbors[18] = n18; 76247c6ae99SBarry Smith dd->neighbors[19] = n19; 76347c6ae99SBarry Smith dd->neighbors[20] = n20; 76447c6ae99SBarry Smith dd->neighbors[21] = n21; 76547c6ae99SBarry Smith dd->neighbors[22] = n22; 76647c6ae99SBarry Smith dd->neighbors[23] = n23; 76747c6ae99SBarry Smith dd->neighbors[24] = n24; 76847c6ae99SBarry Smith dd->neighbors[25] = n25; 76947c6ae99SBarry Smith dd->neighbors[26] = n26; 77047c6ae99SBarry Smith 77147c6ae99SBarry Smith /* If star stencil then delete the corner neighbors */ 772d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 77347c6ae99SBarry Smith /* save information about corner neighbors */ 77447c6ae99SBarry Smith sn0 = n0; sn1 = n1; sn2 = n2; sn3 = n3; sn5 = n5; sn6 = n6; sn7 = n7; 77547c6ae99SBarry Smith sn8 = n8; sn9 = n9; sn11 = n11; sn15 = n15; sn17 = n17; sn18 = n18; 77647c6ae99SBarry Smith sn19 = n19; sn20 = n20; sn21 = n21; sn23 = n23; sn24 = n24; sn25 = n25; 77747c6ae99SBarry Smith sn26 = n26; 7788865f1eaSKarl Rupp n0 = n1 = n2 = n3 = n5 = n6 = n7 = n8 = n9 = n11 = n15 = n17 = n18 = n19 = n20 = n21 = n23 = n24 = n25 = n26 = -1; 77947c6ae99SBarry Smith } 78047c6ae99SBarry Smith 781785e854fSJed Brown ierr = PetscMalloc1((Xe-Xs)*(Ye-Ys)*(Ze-Zs),&idx);CHKERRQ(ierr); 78247c6ae99SBarry Smith 78347c6ae99SBarry Smith nn = 0; 78447c6ae99SBarry Smith /* Bottom Level */ 78547c6ae99SBarry Smith for (k=0; k<s_z; k++) { 78647c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 78747c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 788ce00eea3SSatish Balay x_t = lx[n0 % m]; 78947c6ae99SBarry Smith y_t = ly[(n0 % (m*n))/m]; 79047c6ae99SBarry Smith z_t = lz[n0 / (m*n)]; 79147c6ae99SBarry 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; 7928865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n0] + x_t*y_t*z_t - (s_y-i)*x_t - s_x; /* 2D case */ 7938865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 79447c6ae99SBarry Smith } 79547c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 79647c6ae99SBarry Smith x_t = x; 79747c6ae99SBarry Smith y_t = ly[(n1 % (m*n))/m]; 79847c6ae99SBarry Smith z_t = lz[n1 / (m*n)]; 79947c6ae99SBarry 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; 8008865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t; /* 2D case */ 8018865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 80247c6ae99SBarry Smith } 80347c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 804ce00eea3SSatish Balay x_t = lx[n2 % m]; 80547c6ae99SBarry Smith y_t = ly[(n2 % (m*n))/m]; 80647c6ae99SBarry Smith z_t = lz[n2 / (m*n)]; 80747c6ae99SBarry 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; 8088865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t; /* 2D case */ 8098865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 81047c6ae99SBarry Smith } 81147c6ae99SBarry Smith } 81247c6ae99SBarry Smith 81347c6ae99SBarry Smith for (i=0; i<y; i++) { 81447c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 815ce00eea3SSatish Balay x_t = lx[n3 % m]; 81647c6ae99SBarry Smith y_t = y; 81747c6ae99SBarry Smith z_t = lz[n3 / (m*n)]; 81847c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8198865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - x_t*y_t; /* 2D case */ 8208865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 82147c6ae99SBarry Smith } 82247c6ae99SBarry Smith 82347c6ae99SBarry Smith if (n4 >= 0) { /* middle */ 82447c6ae99SBarry Smith x_t = x; 82547c6ae99SBarry Smith y_t = y; 82647c6ae99SBarry Smith z_t = lz[n4 / (m*n)]; 82747c6ae99SBarry Smith s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8288865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n4] + i*x_t + x_t*y_t*z_t - x_t*y_t; /* 2D case */ 8298865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 830bff4a2f0SMatthew G. Knepley } else if (bz == DM_BOUNDARY_MIRROR) { 8318865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 83247c6ae99SBarry Smith } 83347c6ae99SBarry Smith 83447c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 835ce00eea3SSatish Balay x_t = lx[n5 % m]; 83647c6ae99SBarry Smith y_t = y; 83747c6ae99SBarry Smith z_t = lz[n5 / (m*n)]; 83847c6ae99SBarry Smith s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8398865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n5] + i*x_t + x_t*y_t*z_t - x_t*y_t; /* 2D case */ 8408865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 84147c6ae99SBarry Smith } 84247c6ae99SBarry Smith } 84347c6ae99SBarry Smith 84447c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 84547c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 846ce00eea3SSatish Balay x_t = lx[n6 % m]; 84747c6ae99SBarry Smith y_t = ly[(n6 % (m*n))/m]; 84847c6ae99SBarry Smith z_t = lz[n6 / (m*n)]; 84947c6ae99SBarry Smith s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8508865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - x_t*y_t; /* 2D case */ 8518865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 85247c6ae99SBarry Smith } 85347c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 85447c6ae99SBarry Smith x_t = x; 85547c6ae99SBarry Smith y_t = ly[(n7 % (m*n))/m]; 85647c6ae99SBarry Smith z_t = lz[n7 / (m*n)]; 85747c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8588865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - x_t*y_t; /* 2D case */ 8598865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 86047c6ae99SBarry Smith } 86147c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 862ce00eea3SSatish Balay x_t = lx[n8 % m]; 86347c6ae99SBarry Smith y_t = ly[(n8 % (m*n))/m]; 86447c6ae99SBarry Smith z_t = lz[n8 / (m*n)]; 86547c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 8668865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - x_t*y_t; /* 2D case */ 8678865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 86847c6ae99SBarry Smith } 86947c6ae99SBarry Smith } 87047c6ae99SBarry Smith } 87147c6ae99SBarry Smith 87247c6ae99SBarry Smith /* Middle Level */ 87347c6ae99SBarry Smith for (k=0; k<z; k++) { 87447c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 87547c6ae99SBarry Smith if (n9 >= 0) { /* left below */ 876ce00eea3SSatish Balay x_t = lx[n9 % m]; 87747c6ae99SBarry Smith y_t = ly[(n9 % (m*n))/m]; 87847c6ae99SBarry Smith /* z_t = z; */ 87947c6ae99SBarry Smith s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 8808865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 88147c6ae99SBarry Smith } 88247c6ae99SBarry Smith if (n10 >= 0) { /* directly below */ 88347c6ae99SBarry Smith x_t = x; 88447c6ae99SBarry Smith y_t = ly[(n10 % (m*n))/m]; 88547c6ae99SBarry Smith /* z_t = z; */ 88647c6ae99SBarry Smith s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 8878865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 888bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_MIRROR) { 8898865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 89047c6ae99SBarry Smith } 89147c6ae99SBarry Smith if (n11 >= 0) { /* right below */ 892ce00eea3SSatish Balay x_t = lx[n11 % m]; 89347c6ae99SBarry Smith y_t = ly[(n11 % (m*n))/m]; 89447c6ae99SBarry Smith /* z_t = z; */ 89547c6ae99SBarry Smith s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 8968865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 89747c6ae99SBarry Smith } 89847c6ae99SBarry Smith } 89947c6ae99SBarry Smith 90047c6ae99SBarry Smith for (i=0; i<y; i++) { 90147c6ae99SBarry Smith if (n12 >= 0) { /* directly left */ 902ce00eea3SSatish Balay x_t = lx[n12 % m]; 90347c6ae99SBarry Smith y_t = y; 90447c6ae99SBarry Smith /* z_t = z; */ 90547c6ae99SBarry Smith s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t; 9068865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 907bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { 9088865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = 0; 90947c6ae99SBarry Smith } 91047c6ae99SBarry Smith 91147c6ae99SBarry Smith /* Interior */ 91247c6ae99SBarry Smith s_t = bases[rank] + i*x + k*x*y; 9138865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = s_t++; 91447c6ae99SBarry Smith 91547c6ae99SBarry Smith if (n14 >= 0) { /* directly right */ 916ce00eea3SSatish Balay x_t = lx[n14 % m]; 91747c6ae99SBarry Smith y_t = y; 91847c6ae99SBarry Smith /* z_t = z; */ 91947c6ae99SBarry Smith s_t = bases[n14] + i*x_t + k*x_t*y_t; 9208865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 921bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { 9228865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = 0; 92347c6ae99SBarry Smith } 92447c6ae99SBarry Smith } 92547c6ae99SBarry Smith 92647c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 92747c6ae99SBarry Smith if (n15 >= 0) { /* left above */ 928ce00eea3SSatish Balay x_t = lx[n15 % m]; 92947c6ae99SBarry Smith y_t = ly[(n15 % (m*n))/m]; 93047c6ae99SBarry Smith /* z_t = z; */ 93147c6ae99SBarry Smith s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t; 9328865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 93347c6ae99SBarry Smith } 93447c6ae99SBarry Smith if (n16 >= 0) { /* directly above */ 93547c6ae99SBarry Smith x_t = x; 93647c6ae99SBarry Smith y_t = ly[(n16 % (m*n))/m]; 93747c6ae99SBarry Smith /* z_t = z; */ 93847c6ae99SBarry Smith s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t; 9398865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 940bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_MIRROR) { 9418865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 94247c6ae99SBarry Smith } 94347c6ae99SBarry Smith if (n17 >= 0) { /* right above */ 944ce00eea3SSatish Balay x_t = lx[n17 % m]; 94547c6ae99SBarry Smith y_t = ly[(n17 % (m*n))/m]; 94647c6ae99SBarry Smith /* z_t = z; */ 94747c6ae99SBarry Smith s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t; 9488865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 94947c6ae99SBarry Smith } 95047c6ae99SBarry Smith } 95147c6ae99SBarry Smith } 95247c6ae99SBarry Smith 95347c6ae99SBarry Smith /* Upper Level */ 95447c6ae99SBarry Smith for (k=0; k<s_z; k++) { 95547c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 95647c6ae99SBarry Smith if (n18 >= 0) { /* left below */ 957ce00eea3SSatish Balay x_t = lx[n18 % m]; 95847c6ae99SBarry Smith y_t = ly[(n18 % (m*n))/m]; 95947c6ae99SBarry Smith /* z_t = lz[n18 / (m*n)]; */ 96047c6ae99SBarry Smith s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 9618865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n18] - (s_y-i)*x_t -s_x + x_t*y_t; /* 2d case */ 9628865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 96347c6ae99SBarry Smith } 96447c6ae99SBarry Smith if (n19 >= 0) { /* directly below */ 96547c6ae99SBarry Smith x_t = x; 96647c6ae99SBarry Smith y_t = ly[(n19 % (m*n))/m]; 96747c6ae99SBarry Smith /* z_t = lz[n19 / (m*n)]; */ 96847c6ae99SBarry Smith s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 9698865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n19] - (s_y+1-i)*x_t + x_t*y_t; /* 2d case */ 9708865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 97147c6ae99SBarry Smith } 97247c6ae99SBarry Smith if (n20 >= 0) { /* right below */ 973ce00eea3SSatish Balay x_t = lx[n20 % m]; 97447c6ae99SBarry Smith y_t = ly[(n20 % (m*n))/m]; 97547c6ae99SBarry Smith /* z_t = lz[n20 / (m*n)]; */ 97647c6ae99SBarry Smith s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 9778865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n20] - (s_y+1-i)*x_t + x_t*y_t; /* 2d case */ 9788865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 97947c6ae99SBarry Smith } 98047c6ae99SBarry Smith } 98147c6ae99SBarry Smith 98247c6ae99SBarry Smith for (i=0; i<y; i++) { 98347c6ae99SBarry Smith if (n21 >= 0) { /* directly left */ 984ce00eea3SSatish Balay x_t = lx[n21 % m]; 98547c6ae99SBarry Smith y_t = y; 98647c6ae99SBarry Smith /* z_t = lz[n21 / (m*n)]; */ 98747c6ae99SBarry Smith s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t; 9888865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n21] + (i+1)*x_t - s_x; /* 2d case */ 9898865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 99047c6ae99SBarry Smith } 99147c6ae99SBarry Smith 99247c6ae99SBarry Smith if (n22 >= 0) { /* middle */ 99347c6ae99SBarry Smith x_t = x; 99447c6ae99SBarry Smith y_t = y; 99547c6ae99SBarry Smith /* z_t = lz[n22 / (m*n)]; */ 99647c6ae99SBarry Smith s_t = bases[n22] + i*x_t + k*x_t*y_t; 9978865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n22] + i*x_t; /* 2d case */ 9988865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 999bff4a2f0SMatthew G. Knepley } else if (bz == DM_BOUNDARY_MIRROR) { 10008865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 100147c6ae99SBarry Smith } 100247c6ae99SBarry Smith 100347c6ae99SBarry Smith if (n23 >= 0) { /* directly right */ 1004ce00eea3SSatish Balay x_t = lx[n23 % m]; 100547c6ae99SBarry Smith y_t = y; 100647c6ae99SBarry Smith /* z_t = lz[n23 / (m*n)]; */ 100747c6ae99SBarry Smith s_t = bases[n23] + i*x_t + k*x_t*y_t; 10088865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n23] + i*x_t; /* 2d case */ 10098865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 101047c6ae99SBarry Smith } 101147c6ae99SBarry Smith } 101247c6ae99SBarry Smith 101347c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 101447c6ae99SBarry Smith if (n24 >= 0) { /* left above */ 1015ce00eea3SSatish Balay x_t = lx[n24 % m]; 101647c6ae99SBarry Smith y_t = ly[(n24 % (m*n))/m]; 101747c6ae99SBarry Smith /* z_t = lz[n24 / (m*n)]; */ 101847c6ae99SBarry Smith s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t; 10198865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n24] + i*x_t - s_x; /* 2d case */ 10208865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 102147c6ae99SBarry Smith } 102247c6ae99SBarry Smith if (n25 >= 0) { /* directly above */ 102347c6ae99SBarry Smith x_t = x; 102447c6ae99SBarry Smith y_t = ly[(n25 % (m*n))/m]; 102547c6ae99SBarry Smith /* z_t = lz[n25 / (m*n)]; */ 102647c6ae99SBarry Smith s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t; 10278865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n25] + (i-1)*x_t; /* 2d case */ 10288865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 102947c6ae99SBarry Smith } 103047c6ae99SBarry Smith if (n26 >= 0) { /* right above */ 1031ce00eea3SSatish Balay x_t = lx[n26 % m]; 103247c6ae99SBarry Smith y_t = ly[(n26 % (m*n))/m]; 103347c6ae99SBarry Smith /* z_t = lz[n26 / (m*n)]; */ 103447c6ae99SBarry Smith s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t; 10358865f1eaSKarl Rupp if (twod && (s_t >= M*N*P)) s_t = bases[n26] + (i-1)*x_t; /* 2d case */ 10368865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 103747c6ae99SBarry Smith } 103847c6ae99SBarry Smith } 103947c6ae99SBarry Smith } 1040ce00eea3SSatish Balay 1041b1fb7eb7SBarry Smith ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_USE_POINTER,&from);CHKERRQ(ierr); 104247c6ae99SBarry Smith ierr = VecScatterCreate(global,from,local,to,>ol);CHKERRQ(ierr); 10433bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)gtol);CHKERRQ(ierr); 1044fcfd50ebSBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 1045fcfd50ebSBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 104647c6ae99SBarry Smith 1047d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 104847c6ae99SBarry Smith n0 = sn0; n1 = sn1; n2 = sn2; n3 = sn3; n5 = sn5; n6 = sn6; n7 = sn7; 104947c6ae99SBarry Smith n8 = sn8; n9 = sn9; n11 = sn11; n15 = sn15; n17 = sn17; n18 = sn18; 105047c6ae99SBarry Smith n19 = sn19; n20 = sn20; n21 = sn21; n23 = sn23; n24 = sn24; n25 = sn25; 105147c6ae99SBarry Smith n26 = sn26; 1052ce00eea3SSatish Balay } 105347c6ae99SBarry Smith 1054288e7d53SBarry Smith if (((stencil_type == DMDA_STENCIL_STAR) || (bx != DM_BOUNDARY_PERIODIC && bx) || (by != DM_BOUNDARY_PERIODIC && by) || (bz != DM_BOUNDARY_PERIODIC && bz))) { 1055ce00eea3SSatish Balay /* 1056ce00eea3SSatish Balay Recompute the local to global mappings, this time keeping the 1057ce00eea3SSatish Balay information about the cross corner processor numbers. 1058ce00eea3SSatish Balay */ 105947c6ae99SBarry Smith nn = 0; 106047c6ae99SBarry Smith /* Bottom Level */ 106147c6ae99SBarry Smith for (k=0; k<s_z; k++) { 106247c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 106347c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 1064ce00eea3SSatish Balay x_t = lx[n0 % m]; 106547c6ae99SBarry Smith y_t = ly[(n0 % (m*n))/m]; 106647c6ae99SBarry Smith z_t = lz[n0 / (m*n)]; 106747c6ae99SBarry 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; 10688865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1069ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0 && Zs-zs < 0) { 10708865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 107147c6ae99SBarry Smith } 107247c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 107347c6ae99SBarry Smith x_t = x; 107447c6ae99SBarry Smith y_t = ly[(n1 % (m*n))/m]; 107547c6ae99SBarry Smith z_t = lz[n1 / (m*n)]; 107647c6ae99SBarry 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; 10778865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1078ce00eea3SSatish Balay } else if (Ys-ys < 0 && Zs-zs < 0) { 10798865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 108047c6ae99SBarry Smith } 108147c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 1082ce00eea3SSatish Balay x_t = lx[n2 % m]; 108347c6ae99SBarry Smith y_t = ly[(n2 % (m*n))/m]; 108447c6ae99SBarry Smith z_t = lz[n2 / (m*n)]; 108547c6ae99SBarry 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; 10868865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1087ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0 && Zs-zs < 0) { 10888865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 108947c6ae99SBarry Smith } 109047c6ae99SBarry Smith } 109147c6ae99SBarry Smith 109247c6ae99SBarry Smith for (i=0; i<y; i++) { 109347c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 1094ce00eea3SSatish Balay x_t = lx[n3 % m]; 109547c6ae99SBarry Smith y_t = y; 109647c6ae99SBarry Smith z_t = lz[n3 / (m*n)]; 109747c6ae99SBarry Smith s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 10988865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1099ce00eea3SSatish Balay } else if (Xs-xs < 0 && Zs-zs < 0) { 11008865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 110147c6ae99SBarry Smith } 110247c6ae99SBarry Smith 110347c6ae99SBarry Smith if (n4 >= 0) { /* middle */ 110447c6ae99SBarry Smith x_t = x; 110547c6ae99SBarry Smith y_t = y; 110647c6ae99SBarry Smith z_t = lz[n4 / (m*n)]; 110747c6ae99SBarry Smith s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11088865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1109ce00eea3SSatish Balay } else if (Zs-zs < 0) { 1110bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_MIRROR) { 11118865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 1112c2859e5eSBarry Smith } else { 11138865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 111447c6ae99SBarry Smith } 1115c2859e5eSBarry Smith } 111647c6ae99SBarry Smith 111747c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 1118ce00eea3SSatish Balay x_t = lx[n5 % m]; 111947c6ae99SBarry Smith y_t = y; 112047c6ae99SBarry Smith z_t = lz[n5 / (m*n)]; 112147c6ae99SBarry Smith s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11228865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1123ce00eea3SSatish Balay } else if (xe-Xe < 0 && Zs-zs < 0) { 11248865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 112547c6ae99SBarry Smith } 112647c6ae99SBarry Smith } 112747c6ae99SBarry Smith 112847c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 112947c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 1130ce00eea3SSatish Balay x_t = lx[n6 % m]; 113147c6ae99SBarry Smith y_t = ly[(n6 % (m*n))/m]; 113247c6ae99SBarry Smith z_t = lz[n6 / (m*n)]; 113347c6ae99SBarry Smith s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11348865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1135ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0 && Zs-zs < 0) { 11368865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 113747c6ae99SBarry Smith } 113847c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 113947c6ae99SBarry Smith x_t = x; 114047c6ae99SBarry Smith y_t = ly[(n7 % (m*n))/m]; 114147c6ae99SBarry Smith z_t = lz[n7 / (m*n)]; 114247c6ae99SBarry Smith s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11438865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1144ce00eea3SSatish Balay } else if (ye-Ye < 0 && Zs-zs < 0) { 11458865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 114647c6ae99SBarry Smith } 114747c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 1148ce00eea3SSatish Balay x_t = lx[n8 % m]; 114947c6ae99SBarry Smith y_t = ly[(n8 % (m*n))/m]; 115047c6ae99SBarry Smith z_t = lz[n8 / (m*n)]; 115147c6ae99SBarry Smith s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t; 11528865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1153ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0 && Zs-zs < 0) { 11548865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 115547c6ae99SBarry Smith } 115647c6ae99SBarry Smith } 115747c6ae99SBarry Smith } 115847c6ae99SBarry Smith 115947c6ae99SBarry Smith /* Middle Level */ 116047c6ae99SBarry Smith for (k=0; k<z; k++) { 116147c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 116247c6ae99SBarry Smith if (n9 >= 0) { /* left below */ 1163ce00eea3SSatish Balay x_t = lx[n9 % m]; 116447c6ae99SBarry Smith y_t = ly[(n9 % (m*n))/m]; 116547c6ae99SBarry Smith /* z_t = z; */ 116647c6ae99SBarry Smith s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 11678865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1168ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0) { 11698865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 117047c6ae99SBarry Smith } 117147c6ae99SBarry Smith if (n10 >= 0) { /* directly below */ 117247c6ae99SBarry Smith x_t = x; 117347c6ae99SBarry Smith y_t = ly[(n10 % (m*n))/m]; 117447c6ae99SBarry Smith /* z_t = z; */ 117547c6ae99SBarry Smith s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 11768865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1177ce00eea3SSatish Balay } else if (Ys-ys < 0) { 1178bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_MIRROR) { 11798865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 1180c2859e5eSBarry Smith } else { 11818865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 1182c2859e5eSBarry Smith } 118347c6ae99SBarry Smith } 118447c6ae99SBarry Smith if (n11 >= 0) { /* right below */ 1185ce00eea3SSatish Balay x_t = lx[n11 % m]; 118647c6ae99SBarry Smith y_t = ly[(n11 % (m*n))/m]; 118747c6ae99SBarry Smith /* z_t = z; */ 118847c6ae99SBarry Smith s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 11898865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1190ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0) { 11918865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 119247c6ae99SBarry Smith } 119347c6ae99SBarry Smith } 119447c6ae99SBarry Smith 119547c6ae99SBarry Smith for (i=0; i<y; i++) { 119647c6ae99SBarry Smith if (n12 >= 0) { /* directly left */ 1197ce00eea3SSatish Balay x_t = lx[n12 % m]; 119847c6ae99SBarry Smith y_t = y; 119947c6ae99SBarry Smith /* z_t = z; */ 120047c6ae99SBarry Smith s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t; 12018865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1202ce00eea3SSatish Balay } else if (Xs-xs < 0) { 1203bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_MIRROR) { 12048865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = 0; 1205c2859e5eSBarry Smith } else { 12068865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 120747c6ae99SBarry Smith } 1208c2859e5eSBarry Smith } 120947c6ae99SBarry Smith 121047c6ae99SBarry Smith /* Interior */ 121147c6ae99SBarry Smith s_t = bases[rank] + i*x + k*x*y; 12128865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = s_t++; 121347c6ae99SBarry Smith 121447c6ae99SBarry Smith if (n14 >= 0) { /* directly right */ 1215ce00eea3SSatish Balay x_t = lx[n14 % m]; 121647c6ae99SBarry Smith y_t = y; 121747c6ae99SBarry Smith /* z_t = z; */ 121847c6ae99SBarry Smith s_t = bases[n14] + i*x_t + k*x_t*y_t; 12198865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1220ce00eea3SSatish Balay } else if (xe-Xe < 0) { 1221bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_MIRROR) { 12228865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = 0; 1223c2859e5eSBarry Smith } else { 12248865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 122547c6ae99SBarry Smith } 122647c6ae99SBarry Smith } 1227c2859e5eSBarry Smith } 122847c6ae99SBarry Smith 122947c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 123047c6ae99SBarry Smith if (n15 >= 0) { /* left above */ 1231ce00eea3SSatish Balay x_t = lx[n15 % m]; 123247c6ae99SBarry Smith y_t = ly[(n15 % (m*n))/m]; 123347c6ae99SBarry Smith /* z_t = z; */ 123447c6ae99SBarry Smith s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t; 12358865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1236ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0) { 12378865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 123847c6ae99SBarry Smith } 123947c6ae99SBarry Smith if (n16 >= 0) { /* directly above */ 124047c6ae99SBarry Smith x_t = x; 124147c6ae99SBarry Smith y_t = ly[(n16 % (m*n))/m]; 124247c6ae99SBarry Smith /* z_t = z; */ 124347c6ae99SBarry Smith s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t; 12448865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1245ce00eea3SSatish Balay } else if (ye-Ye < 0) { 1246bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_MIRROR) { 12478865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 1248c2859e5eSBarry Smith } else { 12498865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 125047c6ae99SBarry Smith } 1251c2859e5eSBarry Smith } 125247c6ae99SBarry Smith if (n17 >= 0) { /* right above */ 1253ce00eea3SSatish Balay x_t = lx[n17 % m]; 125447c6ae99SBarry Smith y_t = ly[(n17 % (m*n))/m]; 125547c6ae99SBarry Smith /* z_t = z; */ 125647c6ae99SBarry Smith s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t; 12578865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1258ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0) { 12598865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 126047c6ae99SBarry Smith } 126147c6ae99SBarry Smith } 126247c6ae99SBarry Smith } 126347c6ae99SBarry Smith 126447c6ae99SBarry Smith /* Upper Level */ 126547c6ae99SBarry Smith for (k=0; k<s_z; k++) { 126647c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 126747c6ae99SBarry Smith if (n18 >= 0) { /* left below */ 1268ce00eea3SSatish Balay x_t = lx[n18 % m]; 126947c6ae99SBarry Smith y_t = ly[(n18 % (m*n))/m]; 127047c6ae99SBarry Smith /* z_t = lz[n18 / (m*n)]; */ 127147c6ae99SBarry Smith s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t; 12728865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1273ce00eea3SSatish Balay } else if (Xs-xs < 0 && Ys-ys < 0 && ze-Ze < 0) { 12748865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 127547c6ae99SBarry Smith } 127647c6ae99SBarry Smith if (n19 >= 0) { /* directly below */ 127747c6ae99SBarry Smith x_t = x; 127847c6ae99SBarry Smith y_t = ly[(n19 % (m*n))/m]; 127947c6ae99SBarry Smith /* z_t = lz[n19 / (m*n)]; */ 128047c6ae99SBarry Smith s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 12818865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1282ce00eea3SSatish Balay } else if (Ys-ys < 0 && ze-Ze < 0) { 12838865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 128447c6ae99SBarry Smith } 128547c6ae99SBarry Smith if (n20 >= 0) { /* right below */ 1286ce00eea3SSatish Balay x_t = lx[n20 % m]; 128747c6ae99SBarry Smith y_t = ly[(n20 % (m*n))/m]; 128847c6ae99SBarry Smith /* z_t = lz[n20 / (m*n)]; */ 128947c6ae99SBarry Smith s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t; 12908865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1291ce00eea3SSatish Balay } else if (xe-Xe < 0 && Ys-ys < 0 && ze-Ze < 0) { 12928865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 129347c6ae99SBarry Smith } 129447c6ae99SBarry Smith } 129547c6ae99SBarry Smith 129647c6ae99SBarry Smith for (i=0; i<y; i++) { 129747c6ae99SBarry Smith if (n21 >= 0) { /* directly left */ 1298ce00eea3SSatish Balay x_t = lx[n21 % m]; 129947c6ae99SBarry Smith y_t = y; 130047c6ae99SBarry Smith /* z_t = lz[n21 / (m*n)]; */ 130147c6ae99SBarry Smith s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t; 13028865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1303ce00eea3SSatish Balay } else if (Xs-xs < 0 && ze-Ze < 0) { 13048865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 130547c6ae99SBarry Smith } 130647c6ae99SBarry Smith 130747c6ae99SBarry Smith if (n22 >= 0) { /* middle */ 130847c6ae99SBarry Smith x_t = x; 130947c6ae99SBarry Smith y_t = y; 131047c6ae99SBarry Smith /* z_t = lz[n22 / (m*n)]; */ 131147c6ae99SBarry Smith s_t = bases[n22] + i*x_t + k*x_t*y_t; 13128865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1313ce00eea3SSatish Balay } else if (ze-Ze < 0) { 1314bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_MIRROR) { 13158865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = 0; 1316c2859e5eSBarry Smith } else { 13178865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 131847c6ae99SBarry Smith } 1319c2859e5eSBarry Smith } 132047c6ae99SBarry Smith 132147c6ae99SBarry Smith if (n23 >= 0) { /* directly right */ 1322ce00eea3SSatish Balay x_t = lx[n23 % m]; 132347c6ae99SBarry Smith y_t = y; 132447c6ae99SBarry Smith /* z_t = lz[n23 / (m*n)]; */ 132547c6ae99SBarry Smith s_t = bases[n23] + i*x_t + k*x_t*y_t; 13268865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1327ce00eea3SSatish Balay } else if (xe-Xe < 0 && ze-Ze < 0) { 13288865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 132947c6ae99SBarry Smith } 133047c6ae99SBarry Smith } 133147c6ae99SBarry Smith 133247c6ae99SBarry Smith for (i=1; i<=s_y; i++) { 133347c6ae99SBarry Smith if (n24 >= 0) { /* left above */ 1334ce00eea3SSatish Balay x_t = lx[n24 % m]; 133547c6ae99SBarry Smith y_t = ly[(n24 % (m*n))/m]; 133647c6ae99SBarry Smith /* z_t = lz[n24 / (m*n)]; */ 133747c6ae99SBarry Smith s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t; 13388865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1339ce00eea3SSatish Balay } else if (Xs-xs < 0 && ye-Ye < 0 && ze-Ze < 0) { 13408865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 134147c6ae99SBarry Smith } 134247c6ae99SBarry Smith if (n25 >= 0) { /* directly above */ 134347c6ae99SBarry Smith x_t = x; 134447c6ae99SBarry Smith y_t = ly[(n25 % (m*n))/m]; 134547c6ae99SBarry Smith /* z_t = lz[n25 / (m*n)]; */ 134647c6ae99SBarry Smith s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t; 13478865f1eaSKarl Rupp for (j=0; j<x_t; j++) idx[nn++] = s_t++; 1348ce00eea3SSatish Balay } else if (ye-Ye < 0 && ze-Ze < 0) { 13498865f1eaSKarl Rupp for (j=0; j<x; j++) idx[nn++] = -1; 135047c6ae99SBarry Smith } 135147c6ae99SBarry Smith if (n26 >= 0) { /* right above */ 1352ce00eea3SSatish Balay x_t = lx[n26 % m]; 135347c6ae99SBarry Smith y_t = ly[(n26 % (m*n))/m]; 135447c6ae99SBarry Smith /* z_t = lz[n26 / (m*n)]; */ 135547c6ae99SBarry Smith s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t; 13568865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = s_t++; 1357ce00eea3SSatish Balay } else if (xe-Xe < 0 && ye-Ye < 0 && ze-Ze < 0) { 13588865f1eaSKarl Rupp for (j=0; j<s_x; j++) idx[nn++] = -1; 135947c6ae99SBarry Smith } 136047c6ae99SBarry Smith } 136147c6ae99SBarry Smith } 136247c6ae99SBarry Smith } 136347c6ae99SBarry Smith /* 136447c6ae99SBarry Smith Set the local to global ordering in the global vector, this allows use 136547c6ae99SBarry Smith of VecSetValuesLocal(). 136647c6ae99SBarry Smith */ 136745b6f7e9SBarry Smith ierr = ISLocalToGlobalMappingCreate(comm,dof,nn,idx,PETSC_OWN_POINTER,&da->ltogmap);CHKERRQ(ierr); 13683bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)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->base = base; 1381ce00eea3SSatish Balay da->ops->view = DMView_DA_3d; 13820298fd71SBarry Smith dd->ltol = NULL; 13830298fd71SBarry Smith dd->ao = NULL; 138447c6ae99SBarry Smith PetscFunctionReturn(0); 138547c6ae99SBarry Smith } 1386564755cdSBarry Smith 138747c6ae99SBarry Smith 138847c6ae99SBarry Smith /*@C 1389aa219208SBarry Smith DMDACreate3d - Creates an object that will manage the communication of three-dimensional 139047c6ae99SBarry Smith regular array data that is distributed across some processors. 139147c6ae99SBarry Smith 139247c6ae99SBarry Smith Collective on MPI_Comm 139347c6ae99SBarry Smith 139447c6ae99SBarry Smith Input Parameters: 139547c6ae99SBarry Smith + comm - MPI communicator 13961321219cSEthan Coon . bx,by,bz - type of ghost nodes the array have. 1397bff4a2f0SMatthew G. Knepley Use one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC. 1398aa219208SBarry Smith . stencil_type - Type of stencil (DMDA_STENCIL_STAR or DMDA_STENCIL_BOX) 1399897f7067SBarry Smith . M,N,P - global dimension in each direction of the array 140047c6ae99SBarry Smith . m,n,p - corresponding number of processors in each dimension 140147c6ae99SBarry Smith (or PETSC_DECIDE to have calculated) 140247c6ae99SBarry Smith . dof - number of degrees of freedom per node 140310d7c030SMatthew G Knepley . s - stencil width 140410d7c030SMatthew G Knepley - lx, ly, lz - arrays containing the number of nodes in each cell along 14050298fd71SBarry Smith the x, y, and z coordinates, or NULL. If non-null, these 140647c6ae99SBarry Smith must be of length as m,n,p and the corresponding 140747c6ae99SBarry Smith m,n, or p cannot be PETSC_DECIDE. Sum of the lx[] entries must be M, sum of 140847c6ae99SBarry Smith the ly[] must N, sum of the lz[] must be P 140947c6ae99SBarry Smith 141047c6ae99SBarry Smith Output Parameter: 141147c6ae99SBarry Smith . da - the resulting distributed array object 141247c6ae99SBarry Smith 141347c6ae99SBarry Smith Options Database Key: 1414706a11cbSBarry Smith + -dm_view - Calls DMView() at the conclusion of DMDACreate3d() 141547c6ae99SBarry Smith . -da_grid_x <nx> - number of grid points in x direction, if M < 0 141647c6ae99SBarry Smith . -da_grid_y <ny> - number of grid points in y direction, if N < 0 141747c6ae99SBarry Smith . -da_grid_z <nz> - number of grid points in z direction, if P < 0 141847c6ae99SBarry Smith . -da_processors_x <MX> - number of processors in x direction 141947c6ae99SBarry Smith . -da_processors_y <MY> - number of processors in y direction 142047c6ae99SBarry Smith . -da_processors_z <MZ> - number of processors in z direction 1421e0f5d30fSBarry Smith . -da_refine_x <rx> - refinement ratio in x direction 1422e0f5d30fSBarry Smith . -da_refine_y <ry> - refinement ratio in y direction 1423e0f5d30fSBarry Smith . -da_refine_z <rz>- refinement ratio in z directio 1424e0f5d30fSBarry Smith - -da_refine <n> - refine the DMDA n times before creating it, , if M, N, or P < 0 142547c6ae99SBarry Smith 142647c6ae99SBarry Smith Level: beginner 142747c6ae99SBarry Smith 142847c6ae99SBarry Smith Notes: 1429aa219208SBarry Smith The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the 1430aa219208SBarry Smith standard 7-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes 143147c6ae99SBarry Smith the standard 27-pt stencil. 143247c6ae99SBarry Smith 1433aa219208SBarry Smith The array data itself is NOT stored in the DMDA, it is stored in Vec objects; 1434564755cdSBarry Smith The appropriate vector objects can be obtained with calls to DMCreateGlobalVector() 1435564755cdSBarry Smith and DMCreateLocalVector() and calls to VecDuplicate() if more are needed. 143647c6ae99SBarry Smith 1437897f7067SBarry Smith You must call DMSetUp() after this call before using this DM. 1438897f7067SBarry Smith 1439897f7067SBarry Smith If you wish to use the options database to change values in the DMDA call DMSetFromOptions() after this call 1440897f7067SBarry Smith but before DMSetUp(). 1441897f7067SBarry Smith 144247c6ae99SBarry Smith .keywords: distributed array, create, three-dimensional 144347c6ae99SBarry Smith 1444aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate2d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(), 144599f0b487SRichard Tran Mills DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDASetRefinementFactor(), 1446d461ba97SBarry Smith DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges() 144747c6ae99SBarry Smith 144847c6ae99SBarry Smith @*/ 1449bff4a2f0SMatthew G. Knepley PetscErrorCode DMDACreate3d(MPI_Comm comm,DMBoundaryType bx,DMBoundaryType by,DMBoundaryType bz,DMDAStencilType stencil_type,PetscInt M, 14509a42bb27SBarry 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) 145147c6ae99SBarry Smith { 145247c6ae99SBarry Smith PetscErrorCode ierr; 145347c6ae99SBarry Smith 145447c6ae99SBarry Smith PetscFunctionBegin; 1455aa219208SBarry Smith ierr = DMDACreate(comm, da);CHKERRQ(ierr); 1456c73cfb54SMatthew G. Knepley ierr = DMSetDimension(*da, 3);CHKERRQ(ierr); 1457aa219208SBarry Smith ierr = DMDASetSizes(*da, M, N, P);CHKERRQ(ierr); 1458aa219208SBarry Smith ierr = DMDASetNumProcs(*da, m, n, p);CHKERRQ(ierr); 14591321219cSEthan Coon ierr = DMDASetBoundaryType(*da, bx, by, bz);CHKERRQ(ierr); 1460aa219208SBarry Smith ierr = DMDASetDof(*da, dof);CHKERRQ(ierr); 1461aa219208SBarry Smith ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr); 1462aa219208SBarry Smith ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr); 1463aa219208SBarry Smith ierr = DMDASetOwnershipRanges(*da, lx, ly, lz);CHKERRQ(ierr); 146447c6ae99SBarry Smith PetscFunctionReturn(0); 146547c6ae99SBarry Smith } 1466